Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PSPassageTrackLength.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28// G4PSPassageTrackLength
30#include "G4StepStatus.hh"
31#include "G4Track.hh"
32#include "G4UnitsTable.hh"
33#include "G4VScoreHistFiller.hh"
34
35////////////////////////////////////////////////////////////////////////////////
36// (Description)
37// This is a primitive scorer class for scoring only track length.
38// The tracks which passed a geometry is taken into account.
39//
40//
41// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
42// 2010-07-22 Introduce Unit specification.
43// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of track length
44// vs. number of tracks (not weighted) (Makoto Asai)
45//
46///////////////////////////////////////////////////////////////////////////////
47
49 :G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
50 EvtMap(0),weighted(false)
51{
52 SetUnit("mm");
53}
54
56 const G4String& unit,
57 G4int depth)
58 :G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
59 EvtMap(0),weighted(false)
60{
61 SetUnit(unit);
62}
63
64
66{;}
67
69{
70
71 if ( IsPassed(aStep) ) {
72 G4int index = GetIndex(aStep);
73 EvtMap->add(index,fTrackLength);
74
75 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
76 {
77 auto filler = G4VScoreHistFiller::Instance();
78 if(!filler)
79 {
80 G4Exception("G4PSPassageTrackLength::ProcessHits","SCORER0123",JustWarning,
81 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
82 }
83 else
84 {
85 filler->FillH1(hitIDMap[index],fTrackLength,1.);
86 }
87 }
88
89 }
90
91 return TRUE;
92}
93
95 G4bool Passed = FALSE;
96
97 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
98 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
99
100 G4int trkid = aStep->GetTrack()->GetTrackID();
101 G4double trklength = aStep->GetStepLength();
102 if(weighted) trklength *= aStep->GetPreStepPoint()->GetWeight();
103
104 if ( IsEnter &&IsExit ){ // Passed at one step
105 fTrackLength = trklength; // Track length is absolutely given.
106 Passed = TRUE;
107 }else if ( IsEnter ){ // Enter a new geometry
108 fCurrentTrkID = trkid; // Resetting the current track.
109 fTrackLength = trklength;
110 }else if ( IsExit ){ // Exit a current geometry
111 if ( fCurrentTrkID == trkid ) {
112 fTrackLength += trklength; // Adding the track length to current one,
113 Passed = TRUE; // if the track is same as entered.
114 }
115 }else{ // Inside geometry
116 if ( fCurrentTrkID == trkid ){ // Adding the track length to current one ,
117 fTrackLength += trklength; // if the track is same as entered.
118 }
119 }
120
121 return Passed;
122}
123
125{
126 fCurrentTrkID = -1;
127
128 EvtMap = new G4THitsMap<G4double>(detector->GetName(),GetName());
129 if ( HCID < 0 ) HCID = GetCollectionID(0);
130 HCE->AddHitsCollection(HCID,EvtMap);
131}
132
134{;}
135
137 EvtMap->clear();
138}
139
141{;}
142
144{
145 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
146 G4cout << " PrimitiveSenstivity " << GetName() << G4endl;
147 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
148 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
149 for(; itr != EvtMap->GetMap()->end(); itr++) {
150 G4cout << " copy no.: " << itr->first
151 << " track length : "
152 << *(itr->second)/GetUnitValue()
153 << " [" << GetUnit()<< "]"
154 << G4endl;
155 }
156}
157
159{
160 CheckAndSetUnit(unit,"Length");
161}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
virtual void SetUnit(const G4String &unit)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void EndOfEvent(G4HCofThisEvent *)
virtual void Initialize(G4HCofThisEvent *)
G4PSPassageTrackLength(G4String name, G4int depth=0)
G4StepStatus GetStepStatus() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
void CheckAndSetUnit(const G4String &unit, const G4String &category)
G4double GetUnitValue() const
static G4VScoreHistFiller * Instance()
Map_t * GetMap() const
Definition: G4THitsMap.hh:151
size_t entries() const
Definition: G4THitsMap.hh:155
void clear()
Definition: G4THitsMap.hh:537
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177