Geant4-11
G4PSPassageCellFlux.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// G4PSPassageCellFlux
30
31#include "G4SystemOfUnits.hh"
32#include "G4StepStatus.hh"
33#include "G4Track.hh"
34#include "G4VSolid.hh"
35#include "G4VPhysicalVolume.hh"
37#include "G4UnitsTable.hh"
38#include "G4VScoreHistFiller.hh"
39
41// (Description)
42// This is a primitive scorer class for scoring cell flux.
43// The Cell Flux is defined by a track length divided by a geometry
44// volume, where only tracks passing through the geometry are taken
45// into account. e.g. the unit of Cell Flux is mm/mm3.
46//
47// If you want to score all tracks in the geometry volume,
48// please use G4PSCellFlux.
49//
50// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
51// 2010-07-22 Introduce Unit specification.
52// 2010-07-22 Add weighted option
53// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
54// vs. cell flux * track weight (y) (Makoto Asai)
55//
57
60 , HCID(-1)
61 , fCurrentTrkID(-1)
62 , fCellFlux(0)
63 , EvtMap(0)
64 , weighted(true)
65{
67 SetUnit("percm2");
68}
69
71 G4int depth)
73 , HCID(-1)
74 , fCurrentTrkID(-1)
75 , fCellFlux(0)
76 , EvtMap(0)
77 , weighted(true)
78{
80 SetUnit(unit);
81}
82
84
86{
87 if(IsPassed(aStep))
88 {
89 G4int idx =
91 ->GetReplicaNumber(indexDepth);
92 G4double cubicVolume = ComputeVolume(aStep, idx);
93
94 fCellFlux /= cubicVolume;
95 G4int index = GetIndex(aStep);
96 EvtMap->add(index, fCellFlux);
97
98 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
99 {
100 auto filler = G4VScoreHistFiller::Instance();
101 if(!filler)
102 {
104 "G4PSPassageCellFlux::ProcessHits", "SCORER0123", JustWarning,
105 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
106 }
107 else
108 {
109 filler->FillH1(hitIDMap[index],
111 }
112 }
113 }
114
115 return TRUE;
116}
117
119{
120 G4bool Passed = FALSE;
121
122 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
123 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
124
125 G4int trkid = aStep->GetTrack()->GetTrackID();
126 G4double trklength = aStep->GetStepLength();
127 if(weighted)
128 trklength *= aStep->GetPreStepPoint()->GetWeight();
129
130 if(IsEnter && IsExit)
131 { // Passed at one step
132 fCellFlux = trklength; // Track length is absolutely given.
133 Passed = TRUE;
134 }
135 else if(IsEnter)
136 { // Enter a new geometry
137 fCurrentTrkID = trkid; // Resetting the current track.
138 fCellFlux = trklength;
139 }
140 else if(IsExit)
141 { // Exit a current geometry
142 if(fCurrentTrkID == trkid)
143 { // if the track is same as entered,
144 fCellFlux += trklength; // add the track length to current one.
145 Passed = TRUE;
146 }
147 }
148 else
149 { // Inside geometry
150 if(fCurrentTrkID == trkid)
151 { // if the track is same as entered,
152 fCellFlux += trklength; // adding the track length to current one.
153 }
154 }
155
156 return Passed;
157}
158
160{
161 fCurrentTrkID = -1;
162
164 if(HCID < 0)
167}
168
170
172
174
176{
177 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
178 G4cout << " PrimitiveScorer " << GetName() << G4endl;
179 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
180 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
181 for(; itr != EvtMap->GetMap()->end(); itr++)
182 {
183 G4cout << " copy no.: " << itr->first
184 << " cell flux : " << *(itr->second) / GetUnitValue() << " ["
185 << GetUnit() << G4endl;
186 }
187}
188
190{
191 CheckAndSetUnit(unit, "Per Unit Surface");
192}
193
195{
196 // Per Unit Surface
197 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
198 (1. / cm2));
199 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
200 (1. / mm2));
201 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
202}
203
205{
206 return ComputeSolid(aStep, idx)->GetCubicVolume();
207}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double mm2
Definition: G4SIunits.hh:96
static constexpr double cm2
Definition: G4SIunits.hh:100
static constexpr double m2
Definition: G4SIunits.hh:110
@ fGeomBoundary
Definition: G4StepStatus.hh:43
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 G4double ComputeVolume(G4Step *, G4int idx)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void Initialize(G4HCofThisEvent *)
virtual void SetUnit(const G4String &unit)
G4PSPassageCellFlux(G4String name, G4int depth=0)
G4THitsMap< G4double > * EvtMap
virtual void EndOfEvent(G4HCofThisEvent *)
virtual void DefineUnitAndCategory()
virtual G4bool IsPassed(G4Step *)
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
G4VSolid * ComputeSolid(G4Step *aStep, G4int replicaIdx)
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
void CheckAndSetUnit(const G4String &unit, const G4String &category)
G4double GetUnitValue() const
static G4VScoreHistFiller * Instance()
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
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
const char * name(G4int ptype)