Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
45 : G4PSMinKinEAtGeneration(name, "MeV", depth)
46{}
47
49 const G4String& unit,
50 G4int depth)
51 : G4VPrimitivePlotter(name, depth)
52 , HCID(-1)
53 , EvtMap(nullptr)
54{
55 SetUnit(unit);
56}
57
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}
116
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}
126
128
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}
140
142{
143 CheckAndSetUnit(unit, "Energy");
144}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
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
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
virtual void SetUnit(const G4String &unit)
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
G4PSMinKinEAtGeneration(G4String name, G4int depth=0)
void Initialize(G4HCofThisEvent *) override
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