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

#include <G4PSPopulation.hh>

+ Inheritance diagram for G4PSPopulation:

Public Member Functions

 G4PSPopulation (G4String name, G4int depth=0)
 
virtual ~G4PSPopulation ()
 
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 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)
 

Additional Inherited Members

- 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 49 of file G4PSPopulation.hh.

Constructor & Destructor Documentation

◆ G4PSPopulation()

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

Definition at line 42 of file G4PSPopulation.cc.

43 :G4VPrimitiveScorer(name,depth),HCID(-1),weighted(false)
44{
45 SetUnit("");
46}
virtual void SetUnit(const G4String &unit)

◆ ~G4PSPopulation()

G4PSPopulation::~G4PSPopulation ( )
virtual

Definition at line 48 of file G4PSPopulation.cc.

49{;}

Member Function Documentation

◆ clear()

void G4PSPopulation::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 77 of file G4PSPopulation.cc.

77 {
78 EvtMap->clear();
79 fCellTrackLogger.clear();
80}
void clear()
Definition: G4THitsMap.hh:209

◆ DrawAll()

void G4PSPopulation::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 82 of file G4PSPopulation.cc.

83{;}

◆ EndOfEvent()

void G4PSPopulation::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 72 of file G4PSPopulation.cc.

73{
74 fCellTrackLogger.clear();
75}

◆ Initialize()

void G4PSPopulation::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 65 of file G4PSPopulation.cc.

66{
68 if(HCID < 0) {HCID = GetCollectionID(0);}
69 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
70}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ PrintAll()

void G4PSPopulation::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 85 of file G4PSPopulation.cc.

86{
87 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
88 G4cout << " PrimitiveScorer " << GetName() << G4endl;
89 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
90 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
91 for(; itr != EvtMap->GetMap()->end(); itr++) {
92 G4cout << " copy no.: " << itr->first
93 << " population: " << *(itr->second)/GetUnitValue()
94 << " [tracks]"
95 << G4endl;
96 }
97}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
G4int entries() const
Definition: G4THitsMap.hh:79
G4double GetUnitValue() const

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 51 of file G4PSPopulation.cc.

52{
53
54 G4int index = GetIndex(aStep);
55 G4TrackLogger& tlog = fCellTrackLogger[index];
56 if (tlog.FirstEnterance(aStep->GetTrack()->GetTrackID())) {
57 G4double val = 1.0;
58 if(weighted) val *= aStep->GetPreStepPoint()->GetWeight();
59 EvtMap->add(index,val);
60 }
61
62 return TRUE;
63}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:138
G4bool FirstEnterance(G4int trid)
G4int GetTrackID() const
virtual G4int GetIndex(G4Step *)
#define TRUE
Definition: globals.hh:55

◆ SetUnit()

void G4PSPopulation::SetUnit ( const G4String unit)
virtual

Definition at line 99 of file G4PSPopulation.cc.

100{
101 if (unit == "" ){
102 unitName = unit;
103 unitValue = 1.0;
104 }else{
105 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
106 G4Exception("G4PSPopulation::SetUnit","DetPS0014",JustWarning,msg);
107 }
108
109}
@ JustWarning
const G4String & GetUnit() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4PSPopulation().

◆ Weighted()

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

Definition at line 61 of file G4PSPopulation.hh.

61{ weighted = flg; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


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