Geant4 10.7.0
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)
 
virtual ~G4PSMinKinEAtGeneration ()
 
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 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

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)
 
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 49 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 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0)
46{
47 SetUnit("MeV");
48}
virtual void SetUnit(const G4String &unit)

◆ G4PSMinKinEAtGeneration() [2/2]

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

Definition at line 50 of file G4PSMinKinEAtGeneration.cc.

53 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0)
54{
55 SetUnit(unit);
56}

◆ ~G4PSMinKinEAtGeneration()

G4PSMinKinEAtGeneration::~G4PSMinKinEAtGeneration ( )
virtual

Definition at line 58 of file G4PSMinKinEAtGeneration.cc.

59{;}

Member Function Documentation

◆ clear()

void G4PSMinKinEAtGeneration::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 125 of file G4PSMinKinEAtGeneration.cc.

125 {
126 EvtMap->clear();
127}
void clear()
Definition: G4THitsMap.hh:537

◆ DrawAll()

void G4PSMinKinEAtGeneration::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 129 of file G4PSMinKinEAtGeneration.cc.

130{;}

◆ EndOfEvent()

void G4PSMinKinEAtGeneration::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 122 of file G4PSMinKinEAtGeneration.cc.

123{;}

◆ Initialize()

void G4PSMinKinEAtGeneration::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 115 of file G4PSMinKinEAtGeneration.cc.

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

◆ PrintAll()

void G4PSMinKinEAtGeneration::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 132 of file G4PSMinKinEAtGeneration.cc.

133{
134 G4cout << " PrimitiveScorer " << GetName() << G4endl;
135 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
136 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
137 for(; itr != EvtMap->GetMap()->end(); itr++) {
138 G4cout << " copy no.: " << itr->first
139 << " energy: " << *(itr->second)/GetUnitValue()
140 << " ["<<GetUnit()<<"]"
141 << G4endl;
142 }
143}
#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:151
size_t entries() const
Definition: G4THitsMap.hh:155

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 61 of file G4PSMinKinEAtGeneration.cc.

62{
63 // ===============
64 // First Step,
65 // Confirm this is a newly produced secondary.
66 //
67 //- check for newly produced particle. e.g. Step number is 1.
68 if ( aStep->GetTrack()->GetCurrentStepNumber() != 1) return FALSE;
69 //- check for this is not a primary particle. e.g. ParentID != 0 .
70 if ( aStep->GetTrack()->GetParentID() == 0 ) return FALSE;
71
72 //===============================================
73 //- This is a newly produced secondary particle.
74 //===============================================
75
76 // ===============
77 // Second Step,
78 // Confirm this track has lower energy than previous one.
79 //
80
81 // -Kinetic energy of this particle at the starting point.
82 G4int index = GetIndex(aStep);
83 G4double kinetic = aStep->GetPreStepPoint()->GetKineticEnergy();
84
85 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
86 {
87 auto filler = G4VScoreHistFiller::Instance();
88 if(!filler)
89 {
90 G4Exception("G4PSMinKinEAtGeneration::ProcessHits","SCORER0123",JustWarning,
91 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
92 }
93 else
94 {
95 filler->FillH1(hitIDMap[index],kinetic,aStep->GetPreStepPoint()->GetWeight());
96 }
97 }
98
99 // -Stored value in the current HitsMap.
100 G4double* mapValue= ((*EvtMap)[index]);
101
102 // -
103 // If mapValue exits (e.g not NULL ), compare it with
104 // current track's kinetic energy.
105 if ( mapValue && ( kinetic > *mapValue ) ) return FALSE;
106
107 // -
108 // Current Track is a newly produced secondary and has lower
109 // kinetic energy than previous one in this geometry.
110 //
111 EvtMap->set(index, kinetic);
112 return TRUE;
113}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
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:287

◆ SetUnit()

void G4PSMinKinEAtGeneration::SetUnit ( const G4String unit)
virtual

Definition at line 145 of file G4PSMinKinEAtGeneration.cc.

146{
147 CheckAndSetUnit(unit,"Energy");
148}
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: