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

#include <G4PSNofSecondary.hh>

+ Inheritance diagram for G4PSNofSecondary:

Public Member Functions

 G4PSNofSecondary (G4String name, G4int depth=0)
 
 ~G4PSNofSecondary () override=default
 
void SetParticle (const G4String &particleName)
 
void Weighted (G4bool flg=true)
 
void Initialize (G4HCofThisEvent *) override
 
void clear () override
 
void PrintAll () override
 
virtual void SetUnit (const G4String &unit)
 
- Public Member Functions inherited from G4VPrimitivePlotter
 ~G4VPrimitivePlotter () override=default
 
void Plot (G4int copyNo, G4int histID)
 
G4int GetNumberOfHist () const
 
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
- Public Member Functions inherited from G4VPrimitiveScorer
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
virtual ~G4VPrimitiveScorer ()=default
 
G4int GetCollectionID (G4int)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
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 G4VPrimitivePlotter
std::map< G4int, G4inthitIDMap
 
- 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 53 of file G4PSNofSecondary.hh.

Constructor & Destructor Documentation

◆ G4PSNofSecondary()

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

Definition at line 43 of file G4PSNofSecondary.cc.

◆ ~G4PSNofSecondary()

G4PSNofSecondary::~G4PSNofSecondary ( )
overridedefault

Member Function Documentation

◆ clear()

void G4PSNofSecondary::clear ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 110 of file G4PSNofSecondary.cc.

110{ EvtMap->clear(); }

◆ Initialize()

void G4PSNofSecondary::Initialize ( G4HCofThisEvent * HCE)
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 100 of file G4PSNofSecondary.cc.

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

◆ PrintAll()

void G4PSNofSecondary::PrintAll ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 112 of file G4PSNofSecondary.cc.

113{
114 G4cout << " PrimitiveScorer " << GetName() << G4endl;
115 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
116 for(const auto& [copy, secondaries] : *(EvtMap->GetMap()))
117 {
118 G4cout << " copy no.: " << copy
119 << " num of secondaries: " << *(secondaries) / GetUnitValue()
120 << G4endl;
121 }
122}
#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 G4PSNofSecondary::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )
overrideprotectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 47 of file G4PSNofSecondary.cc.

48{
49 //- check for newly produced particle. e.g. first step.
50 if(aStep->GetTrack()->GetCurrentStepNumber() != 1)
51 return false;
52 //- check for this is not a primary particle. e.g. ParentID > 0 .
53 if(aStep->GetTrack()->GetParentID() == 0)
54 return false;
55 //- check the particle if the partifle definition is given.
56 if((particleDef != nullptr) && particleDef != aStep->GetTrack()->GetDefinition())
57 return false;
58 //
59 //- This is a newly produced secondary particle.
60 G4int index = GetIndex(aStep);
61 G4double weight = 1.0;
62 if(weighted)
63 weight *= aStep->GetPreStepPoint()->GetWeight();
64 EvtMap->add(index, weight);
65
66 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend())
67 {
68 auto filler = G4VScoreHistFiller::Instance();
69 if(filler == nullptr)
70 {
72 "G4PSVolumeFlux::ProcessHits", "SCORER0123", JustWarning,
73 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
74 }
75 else
76 {
77 filler->FillH1(hitIDMap[index],
78 aStep->GetPreStepPoint()->GetKineticEnergy(), weight);
79 }
80 }
81
82 return true;
83}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double GetKineticEnergy() const
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4int GetCurrentStepNumber() const
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
size_t add(const G4int &key, U *&aHit) const

◆ SetParticle()

void G4PSNofSecondary::SetParticle ( const G4String & particleName)

Definition at line 85 of file G4PSNofSecondary.cc.

86{
89 if(pd == nullptr)
90 {
91 G4String msg = "Particle <";
92 msg += particleName;
93 msg += "> not found.";
94 G4Exception("G4PSNofSecondary::SetParticle", "DetPS0101", FatalException,
95 msg);
96 }
97 particleDef = pd;
98}
@ FatalException
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

◆ SetUnit()

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

Definition at line 124 of file G4PSNofSecondary.cc.

125{
126 if(unit.empty())
127 {
128 unitName = unit;
129 unitValue = 1.0;
130 }
131 else
132 {
133 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
134 GetUnit() + "] ) for " + GetName();
135 G4Exception("G4PSNofSecondary::SetUnit", "DetPS0010", JustWarning, msg);
136 }
137}
const G4String & GetUnit() const

◆ Weighted()

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

Definition at line 62 of file G4PSNofSecondary.hh.

62{ weighted = flg; }

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