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

#include <G4PSTrackCounter.hh>

+ Inheritance diagram for G4PSTrackCounter:

Public Member Functions

 G4PSTrackCounter (G4String name, G4int direction, G4int depth=0)
 
virtual ~G4PSTrackCounter ()
 
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 *)
 
- 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 G4PSTrackCounter.hh.

Constructor & Destructor Documentation

◆ G4PSTrackCounter()

G4PSTrackCounter::G4PSTrackCounter ( G4String  name,
G4int  direction,
G4int  depth = 0 
)

Definition at line 42 of file G4PSTrackCounter.cc.

43 :G4VPrimitivePlotter(name,depth),HCID(-1),fDirection(direction),
44 EvtMap(0),weighted(false)
45{
46 SetUnit("");
47}
virtual void SetUnit(const G4String &unit)

◆ ~G4PSTrackCounter()

G4PSTrackCounter::~G4PSTrackCounter ( )
virtual

Definition at line 49 of file G4PSTrackCounter.cc.

50{;}

Member Function Documentation

◆ clear()

void G4PSTrackCounter::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 112 of file G4PSTrackCounter.cc.

112 {
113 EvtMap->clear();
114}
void clear()
Definition: G4THitsMap.hh:537

◆ DrawAll()

void G4PSTrackCounter::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 116 of file G4PSTrackCounter.cc.

117{;}

◆ EndOfEvent()

void G4PSTrackCounter::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 109 of file G4PSTrackCounter.cc.

110{;}

◆ Initialize()

void G4PSTrackCounter::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 101 of file G4PSTrackCounter.cc.

102{
103 EvtMap = new G4THitsMap<G4double>(detector->GetName(),GetName());
104 if(HCID < 0) {HCID = GetCollectionID(0);}
105 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
106 //G4cout << "----- PS Initialize " << G4endl;
107}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ PrintAll()

void G4PSTrackCounter::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 119 of file G4PSTrackCounter.cc.

120{
121 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
122 G4cout << " PrimitiveScorer " << GetName() << G4endl;
123 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
124 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
125 for(; itr != EvtMap->GetMap()->end(); itr++) {
126 G4cout << " copy no.: " << itr->first
127 << " track count: " << *(itr->second)
128 << " [tracks] "
129 << G4endl;
130 }
131}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Map_t * GetMap() const
Definition: G4THitsMap.hh:151
size_t entries() const
Definition: G4THitsMap.hh:155

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 52 of file G4PSTrackCounter.cc.

53{
54
55 G4StepPoint* preStep = aStep->GetPreStepPoint();
56 G4StepPoint* posStep = aStep->GetPostStepPoint();
57 G4bool IsEnter = preStep->GetStepStatus()==fGeomBoundary;
58 G4bool IsExit = posStep->GetStepStatus()==fGeomBoundary;
59
60 // Regular : count in prestep volume.
61 G4int index = GetIndex(aStep);
62
63 // G4cout << " trk " << aStep->GetTrack()->GetTrackID()
64 // << " index " << index << " In " << IsEnter << " Out " <<IsExit
65 // << " PreStatus " << preStep->GetStepStatus()
66 // << " PostStatus " << posStep->GetStepStatus()
67 // << " w " << aStep->GetPreStepPoint()->GetWeight()
68 // << G4endl;
69
70 G4bool flag = FALSE;
71
72 if ( IsEnter && fDirection == fCurrent_In ) flag = TRUE;
73 else if ( IsExit && fDirection == fCurrent_Out ) flag = TRUE;
74 else if ( (IsExit||IsEnter) && fDirection == fCurrent_InOut ) flag = TRUE;
75
76 if ( flag ){
77 //G4cout << " Count " << G4endl;
78 G4double val = 1.0;
79 if(weighted) val *= aStep->GetPreStepPoint()->GetWeight();
80 EvtMap->add(index,val);
81
82 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
83 {
84 auto filler = G4VScoreHistFiller::Instance();
85 if(!filler)
86 {
87 G4Exception("G4PSTrackCounter::ProcessHits","SCORER0123",JustWarning,
88 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
89 }
90 else
91 {
92 filler->FillH1(hitIDMap[index],aStep->GetPreStepPoint()->GetKineticEnergy(),val);
93 }
94 }
95
96 }
97
98 return TRUE;
99}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fCurrent_Out
@ fCurrent_In
@ fCurrent_InOut
@ 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 GetKineticEnergy() const
G4double GetWeight() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
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 G4PSTrackCounter::SetUnit ( const G4String unit)
virtual

Definition at line 133 of file G4PSTrackCounter.cc.

134{
135 if (unit == "" ){
136 unitName = unit;
137 unitValue = 1.0;
138 }else{
139 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
140 G4Exception("G4PSTrackCounter::SetUnit","DetPS0018",JustWarning,msg);
141 }
142
143}
const G4String & GetUnit() const

Referenced by G4PSTrackCounter().

◆ Weighted()

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

Definition at line 60 of file G4PSTrackCounter.hh.

60{ weighted = flg; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


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