Geant4 11.1.1
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 (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

G4bool ProcessHits (G4Step *, G4TouchableHistory *) override
 
- 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 G4PSMinKinEAtGeneration.hh.

Constructor & Destructor Documentation

◆ G4PSMinKinEAtGeneration() [1/2]

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

Definition at line 44 of file G4PSMinKinEAtGeneration.cc.

◆ 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(); }
void clear()
Definition: G4THitsMap.hh:524

◆ 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
G4int GetCollectionID(G4int)

◆ 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:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetUnit() const
G4double GetUnitValue() const
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t entries() const
Definition: G4THitsMap.hh:158
G4double energy(const ThreeVector &p, const G4double m)
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

◆ 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)
Definition: G4Exception.cc:59
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
Definition: G4THitsMap.hh:280

◆ 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: