Geant4 11.2.2
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)
 
 ~G4PSPopulation () override=default
 
void Weighted (G4bool flg=true)
 
void Initialize (G4HCofThisEvent *) override
 
void EndOfEvent (G4HCofThisEvent *) override
 
void clear () override
 
void PrintAll () override
 
virtual void SetUnit (const G4String &unit)
 
- Public Member Functions inherited from G4VPrimitiveScorer
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
virtual ~G4VPrimitiveScorer ()=default
 
G4int GetCollectionID (G4int)
 
virtual void DrawAll ()
 
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

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

Additional Inherited Members

- Protected Attributes inherited from G4VPrimitiveScorer
G4String primitiveName
 
G4MultiFunctionalDetectordetector {nullptr}
 
G4VSDFilterfilter {nullptr}
 
G4int verboseLevel {0}
 
G4int indexDepth
 
G4String unitName {"NoUnit"}
 
G4double unitValue {1.0}
 
G4int fNi {0}
 
G4int fNj {0}
 
G4int fNk {0}
 

Detailed Description

Definition at line 47 of file G4PSPopulation.hh.

Constructor & Destructor Documentation

◆ G4PSPopulation()

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

Definition at line 40 of file G4PSPopulation.cc.

41 : G4VPrimitiveScorer(name, depth)
42{
43 SetUnit("");
44}
virtual void SetUnit(const G4String &unit)
G4VPrimitiveScorer(G4String name, G4int depth=0)

◆ ~G4PSPopulation()

G4PSPopulation::~G4PSPopulation ( )
overridedefault

Member Function Documentation

◆ clear()

void G4PSPopulation::clear ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 73 of file G4PSPopulation.cc.

74{
75 EvtMap->clear();
76 fCellTrackLogger.clear();
77}

◆ EndOfEvent()

void G4PSPopulation::EndOfEvent ( G4HCofThisEvent * )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 71 of file G4PSPopulation.cc.

71{ fCellTrackLogger.clear(); }

◆ Initialize()

void G4PSPopulation::Initialize ( G4HCofThisEvent * HCE)
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 61 of file G4PSPopulation.cc.

62{
63 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
64 if(HCID < 0)
65 {
66 HCID = GetCollectionID(0);
67 }
68 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap);
69}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector

◆ PrintAll()

void G4PSPopulation::PrintAll ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 79 of file G4PSPopulation.cc.

80{
81 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
82 G4cout << " PrimitiveScorer " << GetName() << G4endl;
83 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
84 for(const auto& [copy, population] : *(EvtMap->GetMap()))
85 {
86 G4cout << " copy no.: " << copy
87 << " population: " << *(population) / GetUnitValue() << " [tracks]"
88 << G4endl;
89 }
90}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetUnitValue() const
Map_t * GetMap() const
size_t entries() const
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)

◆ ProcessHits()

G4bool G4PSPopulation::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )
overrideprotectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 46 of file G4PSPopulation.cc.

47{
48 G4int index = GetIndex(aStep);
49 G4TrackLogger& tlog = fCellTrackLogger[index];
50 if(tlog.FirstEnterance(aStep->GetTrack()->GetTrackID()))
51 {
52 G4double val = 1.0;
53 if(weighted)
54 val *= aStep->GetPreStepPoint()->GetWeight();
55 EvtMap->add(index, val);
56 }
57
58 return true;
59}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4bool FirstEnterance(G4int trid)
G4int GetTrackID() const
virtual G4int GetIndex(G4Step *)
size_t add(const G4int &key, U *&aHit) const

◆ SetUnit()

void G4PSPopulation::SetUnit ( const G4String & unit)
virtual

Definition at line 92 of file G4PSPopulation.cc.

93{
94 if(unit.empty())
95 {
96 unitName = unit;
97 unitValue = 1.0;
98 }
99 else
100 {
101 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
102 GetUnit() + "] ) for " + GetName();
103 G4Exception("G4PSPopulation::SetUnit", "DetPS0014", JustWarning, msg);
104 }
105}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
const G4String & GetUnit() const

Referenced by G4PSPopulation().

◆ Weighted()

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

Definition at line 53 of file G4PSPopulation.hh.

53{ weighted = flg; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


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