Geant4 9.6.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 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)
 

Additional Inherited Members

- 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 42 of file G4PSMinKinEAtGeneration.cc.

43 :G4VPrimitiveScorer(name,depth),HCID(-1)
44{
45 SetUnit("MeV");
46}
virtual void SetUnit(const G4String &unit)

◆ G4PSMinKinEAtGeneration() [2/2]

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

Definition at line 48 of file G4PSMinKinEAtGeneration.cc.

51 :G4VPrimitiveScorer(name,depth),HCID(-1)
52{
53 SetUnit(unit);
54}

◆ ~G4PSMinKinEAtGeneration()

G4PSMinKinEAtGeneration::~G4PSMinKinEAtGeneration ( )
virtual

Definition at line 56 of file G4PSMinKinEAtGeneration.cc.

57{;}

Member Function Documentation

◆ clear()

void G4PSMinKinEAtGeneration::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 109 of file G4PSMinKinEAtGeneration.cc.

109 {
110 EvtMap->clear();
111}
void clear()
Definition: G4THitsMap.hh:209

◆ DrawAll()

void G4PSMinKinEAtGeneration::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 113 of file G4PSMinKinEAtGeneration.cc.

114{;}

◆ EndOfEvent()

void G4PSMinKinEAtGeneration::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 106 of file G4PSMinKinEAtGeneration.cc.

107{;}

◆ Initialize()

void G4PSMinKinEAtGeneration::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 99 of file G4PSMinKinEAtGeneration.cc.

100{
101 EvtMap = new G4THitsMap<G4double>(detector->GetName(),GetName());
102 if(HCID < 0) {HCID = GetCollectionID(0);}
103 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
104}
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 116 of file G4PSMinKinEAtGeneration.cc.

117{
118 G4cout << " PrimitiveScorer " << GetName() << G4endl;
119 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
120 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
121 for(; itr != EvtMap->GetMap()->end(); itr++) {
122 G4cout << " copy no.: " << itr->first
123 << " energy: " << *(itr->second)/GetUnitValue()
124 << " ["<<GetUnit()<<"]"
125 << G4endl;
126 }
127}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
G4int entries() const
Definition: G4THitsMap.hh:79
const G4String & GetUnit() const
G4double GetUnitValue() const

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 59 of file G4PSMinKinEAtGeneration.cc.

60{
61 // ===============
62 // First Step,
63 // Confirm this is a newly produced secondary.
64 //
65 //- check for newly produced particle. e.g. Step number is 1.
66 if ( aStep->GetTrack()->GetCurrentStepNumber() != 1) return FALSE;
67 //- check for this is not a primary particle. e.g. ParentID != 0 .
68 if ( aStep->GetTrack()->GetParentID() == 0 ) return FALSE;
69
70 //===============================================
71 //- This is a newly produced secondary particle.
72 //===============================================
73
74 // ===============
75 // Second Step,
76 // Confirm this track has lower energy than previous one.
77 //
78
79 // -Kinetic energy of this particle at the starting point.
80 G4double kinetic = aStep->GetPreStepPoint()->GetKineticEnergy();
81
82 // -Stored value in the current HitsMap.
83 G4int index = GetIndex(aStep);
84 G4double* mapValue= ((*EvtMap)[index]);
85
86 // -
87 // If mapValue exits (e.g not NULL ), compare it with
88 // current track's kinetic energy.
89 if ( mapValue && ( kinetic > *mapValue ) ) return FALSE;
90
91 // -
92 // Current Track is a newly produced secondary and has lower
93 // kinetic energy than previous one in this geometry.
94 //
95 EvtMap->set(index, kinetic);
96 return TRUE;
97}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4int set(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:165
G4int GetCurrentStepNumber() const
G4int GetParentID() const
virtual G4int GetIndex(G4Step *)
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52

◆ SetUnit()

void G4PSMinKinEAtGeneration::SetUnit ( const G4String unit)
virtual

Definition at line 129 of file G4PSMinKinEAtGeneration.cc.

130{
131 CheckAndSetUnit(unit,"Energy");
132}
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: