Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PSPassageTrackLength Class Reference

#include <G4PSPassageTrackLength.hh>

+ Inheritance diagram for G4PSPassageTrackLength:

Public Member Functions

 G4PSPassageTrackLength (G4String name, G4int depth=0)
 
 G4PSPassageTrackLength (G4String name, const G4String &unit, G4int depth=0)
 
virtual ~G4PSPassageTrackLength ()
 
void Weighted (G4bool flg=true)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
virtual void SetUnit (const G4String &unit)
 
- Public Member Functions inherited from G4VPrimitivePlotter
 G4VPrimitivePlotter (G4String name, G4int depth=0)
 
virtual ~G4VPrimitivePlotter ()
 
void Plot (G4int copyNo, G4int histID)
 
G4int GetNumberOfHist () const
 
- Public Member Functions inherited from G4VPrimitiveScorer
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
virtual ~G4VPrimitiveScorer ()
 
G4int GetCollectionID (G4int)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
void SetUnit (const G4String &unit)
 
const G4StringGetUnit () const
 
G4double GetUnitValue () const
 
void SetMultiFunctionalDetector (G4MultiFunctionalDetector *d)
 
G4MultiFunctionalDetectorGetMultiFunctionalDetector () const
 
G4String GetName () const
 
void SetFilter (G4VSDFilter *f)
 
G4VSDFilterGetFilter () const
 
void SetVerboseLevel (G4int vl)
 
G4int GetVerboseLevel () const
 
void SetNijk (G4int i, G4int j, G4int k)
 

Protected Member Functions

virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
 
G4bool IsPassed (G4Step *)
 
- Protected Member Functions inherited from G4VPrimitiveScorer
virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)=0
 
virtual G4int GetIndex (G4Step *)
 
void CheckAndSetUnit (const G4String &unit, const G4String &category)
 
G4VSolidComputeSolid (G4Step *aStep, G4int replicaIdx)
 
G4VSolidComputeCurrentSolid (G4Step *aStep)
 

Additional Inherited Members

- Protected Attributes inherited from G4VPrimitivePlotter
std::map< G4int, G4inthitIDMap
 
- Protected Attributes inherited from G4VPrimitiveScorer
G4String primitiveName
 
G4MultiFunctionalDetectordetector
 
G4VSDFilterfilter
 
G4int verboseLevel
 
G4int indexDepth
 
G4String unitName
 
G4double unitValue
 
G4int fNi
 
G4int fNj
 
G4int fNk
 

Detailed Description

Definition at line 48 of file G4PSPassageTrackLength.hh.

Constructor & Destructor Documentation

◆ G4PSPassageTrackLength() [1/2]

G4PSPassageTrackLength::G4PSPassageTrackLength ( G4String  name,
G4int  depth = 0 
)

Definition at line 48 of file G4PSPassageTrackLength.cc.

49 :G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
50 EvtMap(0),weighted(false)
51{
52 SetUnit("mm");
53}
virtual void SetUnit(const G4String &unit)

◆ G4PSPassageTrackLength() [2/2]

G4PSPassageTrackLength::G4PSPassageTrackLength ( G4String  name,
const G4String unit,
G4int  depth = 0 
)

Definition at line 55 of file G4PSPassageTrackLength.cc.

58 :G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
59 EvtMap(0),weighted(false)
60{
61 SetUnit(unit);
62}

◆ ~G4PSPassageTrackLength()

G4PSPassageTrackLength::~G4PSPassageTrackLength ( )
virtual

Definition at line 65 of file G4PSPassageTrackLength.cc.

66{;}

Member Function Documentation

◆ clear()

void G4PSPassageTrackLength::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 136 of file G4PSPassageTrackLength.cc.

136 {
137 EvtMap->clear();
138}
void clear()
Definition: G4THitsMap.hh:537

◆ DrawAll()

void G4PSPassageTrackLength::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 140 of file G4PSPassageTrackLength.cc.

141{;}

◆ EndOfEvent()

void G4PSPassageTrackLength::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 133 of file G4PSPassageTrackLength.cc.

134{;}

◆ Initialize()

void G4PSPassageTrackLength::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 124 of file G4PSPassageTrackLength.cc.

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}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ IsPassed()

G4bool G4PSPassageTrackLength::IsPassed ( G4Step aStep)
protected

Definition at line 94 of file G4PSPassageTrackLength.cc.

94 {
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}
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
G4StepStatus GetStepStatus() const
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const

Referenced by ProcessHits().

◆ PrintAll()

void G4PSPassageTrackLength::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 143 of file G4PSPassageTrackLength.cc.

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}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetUnit() const
G4double GetUnitValue() const
Map_t * GetMap() const
Definition: G4THitsMap.hh:151
size_t entries() const
Definition: G4THitsMap.hh:155

◆ ProcessHits()

G4bool G4PSPassageTrackLength::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 68 of file G4PSPassageTrackLength.cc.

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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

◆ SetUnit()

void G4PSPassageTrackLength::SetUnit ( const G4String unit)
virtual

Definition at line 158 of file G4PSPassageTrackLength.cc.

159{
160 CheckAndSetUnit(unit,"Length");
161}
void CheckAndSetUnit(const G4String &unit, const G4String &category)

Referenced by G4PSPassageTrackLength(), G4PSPassageTrackLength3D::G4PSPassageTrackLength3D(), and G4ScoreQuantityMessenger::SetNewValue().

◆ Weighted()

void G4PSPassageTrackLength::Weighted ( G4bool  flg = true)
inline

Definition at line 57 of file G4PSPassageTrackLength.hh.

57{ weighted = flg; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


The documentation for this class was generated from the following files: