Geant4-11
G4PSMinKinEAtGeneration.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28// G4PSMinKinEAtGeneration
30#include "G4UnitsTable.hh"
31#include "G4VScoreHistFiller.hh"
32
33// (Description)
34// This is a primitive scorer class for scoring minimum energy of
35// newly produced particles in the geometry.
36//
37// Created: 2005-11-17 Tsukasa ASO, Akinori Kimura.
38// 2010-07-22 Introduce Unit specification.
39// 2011-09-09 Modify comment in PrintAll(). T.Aso
40// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
41// vs. track weight (y) (Makoto Asai)
42//
43
46 , HCID(-1)
47 , EvtMap(0)
48{
49 SetUnit("MeV");
50}
51
53 const G4String& unit,
54 G4int depth)
56 , HCID(-1)
57 , EvtMap(0)
58{
59 SetUnit(unit);
60}
61
63
65{
66 // ===============
67 // First Step,
68 // Confirm this is a newly produced secondary.
69 //
70 //- check for newly produced particle. e.g. Step number is 1.
71 if(aStep->GetTrack()->GetCurrentStepNumber() != 1)
72 return FALSE;
73 //- check for this is not a primary particle. e.g. ParentID != 0 .
74 if(aStep->GetTrack()->GetParentID() == 0)
75 return FALSE;
76
77 //===============================================
78 //- This is a newly produced secondary particle.
79 //===============================================
80
81 // ===============
82 // Second Step,
83 // Confirm this track has lower energy than previous one.
84 //
85
86 // -Kinetic energy of this particle at the starting point.
87 G4int index = GetIndex(aStep);
88 G4double kinetic = aStep->GetPreStepPoint()->GetKineticEnergy();
89
90 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
91 {
92 auto filler = G4VScoreHistFiller::Instance();
93 if(!filler)
94 {
96 "G4PSMinKinEAtGeneration::ProcessHits", "SCORER0123", JustWarning,
97 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
98 }
99 else
100 {
101 filler->FillH1(hitIDMap[index], kinetic,
102 aStep->GetPreStepPoint()->GetWeight());
103 }
104 }
105
106 // -Stored value in the current HitsMap.
107 G4double* mapValue = ((*EvtMap)[index]);
108
109 // -
110 // If mapValue exits (e.g not NULL ), compare it with
111 // current track's kinetic energy.
112 if(mapValue && (kinetic > *mapValue))
113 return FALSE;
114
115 // -
116 // Current Track is a newly produced secondary and has lower
117 // kinetic energy than previous one in this geometry.
118 //
119 EvtMap->set(index, kinetic);
120 return TRUE;
121}
122
124{
126 if(HCID < 0)
127 {
129 }
131}
132
134
136
138
140{
141 G4cout << " PrimitiveScorer " << GetName() << G4endl;
142 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
143 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
144 for(; itr != EvtMap->GetMap()->end(); itr++)
145 {
146 G4cout << " copy no.: " << itr->first
147 << " energy: " << *(itr->second) / GetUnitValue() << " ["
148 << GetUnit() << "]" << G4endl;
149 }
150}
151
153{
154 CheckAndSetUnit(unit, "Energy");
155}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void SetUnit(const G4String &unit)
G4PSMinKinEAtGeneration(G4String name, G4int depth=0)
G4THitsMap< G4double > * EvtMap
virtual void Initialize(G4HCofThisEvent *)
virtual void EndOfEvent(G4HCofThisEvent *)
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4int GetCurrentStepNumber() const
G4int GetParentID() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
void CheckAndSetUnit(const G4String &unit, const G4String &category)
G4double GetUnitValue() const
static G4VScoreHistFiller * Instance()
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t set(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:280
size_t entries() const
Definition: G4THitsMap.hh:158
void clear()
Definition: G4THitsMap.hh:524
const char * name(G4int ptype)