Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
40////////////////////////////////////////////////////////////////////////////////
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//
56///////////////////////////////////////////////////////////////////////////////
57
59 : G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
60 EvtMap(0),weighted(true)
61{
63 SetUnit("percm2");
64}
65
67 G4int depth)
68 : G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
69 EvtMap(0),weighted(true)
70{
72 SetUnit(unit);
73}
74
76{;}
77
79{
80
81 if ( IsPassed(aStep) ) {
83 (aStep->GetPreStepPoint()->GetTouchable()))
84 ->GetReplicaNumber(indexDepth);
85 G4double cubicVolume = ComputeVolume(aStep, idx);
86
87 fCellFlux /= cubicVolume;
88 G4int index = GetIndex(aStep);
89 EvtMap->add(index,fCellFlux);
90
91 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
92 {
93 auto filler = G4VScoreHistFiller::Instance();
94 if(!filler)
95 {
96 G4Exception("G4PSPassageCellFlux::ProcessHits","SCORER0123",JustWarning,
97 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
98 }
99 else
100 {
101 filler->FillH1(hitIDMap[index],aStep->GetPreStepPoint()->GetKineticEnergy(),fCellFlux);
102 }
103 }
104
105 }
106
107 return TRUE;
108}
109
111 G4bool Passed = FALSE;
112
113 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
114 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
115
116 G4int trkid = aStep->GetTrack()->GetTrackID();
117 G4double trklength = aStep->GetStepLength();
118 if ( weighted ) trklength *= aStep->GetPreStepPoint()->GetWeight();
119
120 if ( IsEnter &&IsExit ){ // Passed at one step
121 fCellFlux = trklength; // Track length is absolutely given.
122 Passed = TRUE;
123 }else if ( IsEnter ){ // Enter a new geometry
124 fCurrentTrkID = trkid; // Resetting the current track.
125 fCellFlux = trklength;
126 }else if ( IsExit ){ // Exit a current geometry
127 if ( fCurrentTrkID == trkid ) {// if the track is same as entered,
128 fCellFlux += trklength; // add the track length to current one.
129 Passed = TRUE;
130 }
131 }else{ // Inside geometry
132 if ( fCurrentTrkID == trkid ){ // if the track is same as entered,
133 fCellFlux += trklength; // adding the track length to current one.
134 }
135 }
136
137 return Passed;
138}
139
141{
142 fCurrentTrkID = -1;
143
144 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
145 GetName());
146 if ( HCID < 0 ) HCID = GetCollectionID(0);
147 HCE->AddHitsCollection(HCID,EvtMap);
148
149}
150
152{;}
153
155 EvtMap->clear();
156}
157
159{;}
160
162{
163 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
164 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
165 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
166 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
167 for(; itr != EvtMap->GetMap()->end(); itr++) {
168 G4cout << " copy no.: " << itr->first
169 << " cell flux : " << *(itr->second)/GetUnitValue()
170 << " [" << GetUnit()
171 << G4endl;
172 }
173}
174
176{
177 CheckAndSetUnit(unit,"Per Unit Surface");
178}
179
181 // Per Unit Surface
182 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
183 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
184 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
185}
186
188{
189 return ComputeSolid(aStep, idx)->GetCubicVolume();
190}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ 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)
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:176
Map_t * GetMap() const
Definition: G4THitsMap.hh:151
size_t entries() const
Definition: G4THitsMap.hh:155
void clear()
Definition: G4THitsMap.hh:537
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177