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

#include <G4PSMinKinEAtGeneration.hh>

+ Inheritance diagram for G4PSMinKinEAtGeneration:

Public Member Functions

 G4PSMinKinEAtGeneration (G4String name, G4int depth=0)
 
 G4PSMinKinEAtGeneration (G4String name, const G4String &unit, G4int depth=0)
 
 ~G4PSMinKinEAtGeneration () override=default
 
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 48 of file G4PSMinKinEAtGeneration.hh.

Constructor & Destructor Documentation

◆ G4PSMinKinEAtGeneration() [1/2]

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

Definition at line 44 of file G4PSMinKinEAtGeneration.cc.

45 : G4PSMinKinEAtGeneration(name, "MeV", depth)
46{}
G4PSMinKinEAtGeneration(G4String name, G4int depth=0)

◆ G4PSMinKinEAtGeneration() [2/2]

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

Definition at line 48 of file G4PSMinKinEAtGeneration.cc.

51 : G4VPrimitivePlotter(name, depth)
52 , HCID(-1)
53 , EvtMap(nullptr)
54{
55 SetUnit(unit);
56}
virtual void SetUnit(const G4String &unit)

◆ ~G4PSMinKinEAtGeneration()

G4PSMinKinEAtGeneration::~G4PSMinKinEAtGeneration ( )
overridedefault

Member Function Documentation

◆ clear()

void G4PSMinKinEAtGeneration::clear ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 127 of file G4PSMinKinEAtGeneration.cc.

127{ EvtMap->clear(); }

◆ Initialize()

void G4PSMinKinEAtGeneration::Initialize ( G4HCofThisEvent * HCE)
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 117 of file G4PSMinKinEAtGeneration.cc.

118{
119 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
120 if(HCID < 0)
121 {
122 HCID = GetCollectionID(0);
123 }
124 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap);
125}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector

◆ PrintAll()

void G4PSMinKinEAtGeneration::PrintAll ( )
overridevirtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 129 of file G4PSMinKinEAtGeneration.cc.

130{
131 G4cout << " PrimitiveScorer " << GetName() << G4endl;
132 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
133 for(const auto& [copy, energy] : *(EvtMap->GetMap()))
134 {
135 G4cout << " copy no.: " << copy
136 << " energy: " << *(energy) / GetUnitValue() << " ["
137 << GetUnit() << "]" << G4endl;
138 }
139}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetUnit() const
G4double GetUnitValue() const
Map_t * GetMap() const
size_t entries() const
G4double energy(const ThreeVector &p, const G4double m)
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 58 of file G4PSMinKinEAtGeneration.cc.

59{
60 // ===============
61 // First Step,
62 // Confirm this is a newly produced secondary.
63 //
64 //- check for newly produced particle. e.g. Step number is 1.
65 if(aStep->GetTrack()->GetCurrentStepNumber() != 1)
66 return false;
67 //- check for this is not a primary particle. e.g. ParentID != 0 .
68 if(aStep->GetTrack()->GetParentID() == 0)
69 return false;
70
71 //===============================================
72 //- This is a newly produced secondary particle.
73 //===============================================
74
75 // ===============
76 // Second Step,
77 // Confirm this track has lower energy than previous one.
78 //
79
80 // -Kinetic energy of this particle at the starting point.
81 G4int index = GetIndex(aStep);
82 G4double kinetic = aStep->GetPreStepPoint()->GetKineticEnergy();
83
84 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend())
85 {
86 auto filler = G4VScoreHistFiller::Instance();
87 if(filler == nullptr)
88 {
90 "G4PSMinKinEAtGeneration::ProcessHits", "SCORER0123", JustWarning,
91 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
92 }
93 else
94 {
95 filler->FillH1(hitIDMap[index], kinetic,
96 aStep->GetPreStepPoint()->GetWeight());
97 }
98 }
99
100 // -Stored value in the current HitsMap.
101 G4double* mapValue = ((*EvtMap)[index]);
102
103 // -
104 // If mapValue exits (e.g not NULL ), compare it with
105 // current track's kinetic energy.
106 if((mapValue != nullptr) && (kinetic > *mapValue))
107 return false;
108
109 // -
110 // Current Track is a newly produced secondary and has lower
111 // kinetic energy than previous one in this geometry.
112 //
113 EvtMap->set(index, kinetic);
114 return true;
115}
@ 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
G4int GetParentID() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
size_t set(const G4int &key, U *&aHit) const

◆ SetUnit()

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

Definition at line 141 of file G4PSMinKinEAtGeneration.cc.

142{
143 CheckAndSetUnit(unit, "Energy");
144}
void CheckAndSetUnit(const G4String &unit, const G4String &category)

Referenced by G4PSMinKinEAtGeneration(), G4PSMinKinEAtGeneration3D::G4PSMinKinEAtGeneration3D(), and G4ScoreQuantityMessenger::SetNewValue().


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