Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PSPassageCellCurrent.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// G4PSPassageCellCurrent
30#include "G4StepStatus.hh"
31#include "G4Track.hh"
32#include "G4VSolid.hh"
33#include "G4UnitsTable.hh"
34#include "G4VScoreHistFiller.hh"
35
36////////////////////////////////////////////////////////////////////////////////
37// (Description)
38// This is a primitive scorer class for scoring number of tracks,
39// where only tracks passing through the geometry are taken
40// into account.
41//
42// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
43// 2010-07-22 Introduce Unit specification.
44// 2010-07-22 Add weighted option
45// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
46// vs. current * track weight (y) (Makoto Asai)
47//
48///////////////////////////////////////////////////////////////////////////////
49
51 : G4VPrimitivePlotter(name, depth)
52{
53 SetUnit("");
54}
55
57{
58 if(IsPassed(aStep))
59 {
60 fCurrent = 1.;
61 if(weighted)
62 fCurrent = aStep->GetPreStepPoint()->GetWeight();
63 G4int index = GetIndex(aStep);
64 EvtMap->add(index, fCurrent);
65
66 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend())
67 {
68 auto filler = G4VScoreHistFiller::Instance();
69 if(filler == nullptr)
70 {
72 "G4PSVolumeFlux::ProcessHits", "SCORER0123", JustWarning,
73 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
74 }
75 else
76 {
77 filler->FillH1(hitIDMap[index],
78 aStep->GetPreStepPoint()->GetKineticEnergy(), fCurrent);
79 }
80 }
81 }
82
83 return true;
84}
85
87{
88 G4bool Passed = false;
89
90 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
91 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
92
93 G4int trkid = aStep->GetTrack()->GetTrackID();
94
95 if(IsEnter && IsExit)
96 { // Passed at one step
97 Passed = true;
98 }
99 else if(IsEnter)
100 { // Enter a new geometry
101 fCurrentTrkID = trkid; // Resetting the current track.
102 }
103 else if(IsExit)
104 { // Exit a current geometry
105 if(fCurrentTrkID == trkid)
106 {
107 Passed = true; // if the track is same as entered.
108 }
109 }
110 else
111 { // Inside geometry
112 if(fCurrentTrkID == trkid)
113 { // Adding the track length to current one ,
114 }
115 }
116 return Passed;
117}
118
120{
121 fCurrentTrkID = -1;
122
123 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
124 if(HCID < 0)
125 HCID = GetCollectionID(0);
126 HCE->AddHitsCollection(HCID, EvtMap);
127}
128
130
132{
133 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
134 G4cout << " PrimitiveScorer " << GetName() << G4endl;
135 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
136 for(const auto& [copy, current] : *(EvtMap->GetMap()))
137 {
138 G4cout << " copy no.: " << copy
139 << " cell current : " << *(current) << " [tracks] " << G4endl;
140 }
141}
142
144{
145 if(unit.empty())
146 {
147 unitName = unit;
148 unitValue = 1.0;
149 }
150 else
151 {
152 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
153 GetUnit() + "] ) for " + GetName();
154 G4Exception("G4PSPassageCellCurrent::SetUnit", "DetPS0012", JustWarning,
155 msg);
156 }
157}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
@ fGeomBoundary
Definition: G4StepStatus.hh:43
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)
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
virtual void SetUnit(const G4String &unit)
virtual G4bool IsPassed(G4Step *)
G4PSPassageCellCurrent(G4String name, G4int depth=0)
void Initialize(G4HCofThisEvent *) override
G4StepStatus GetStepStatus() const
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
static G4VScoreHistFiller * Instance()
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t entries() const
Definition: G4THitsMap.hh:158
void clear()
Definition: G4THitsMap.hh:524
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177