Geant4-11
G4PSFlatSurfaceFlux.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// G4PSFlatSurfaceFlux
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"
39#include "G4VScoreHistFiller.hh"
40
42// (Description)
43// This is a primitive scorer class for scoring Surface Flux.
44// Current version assumes only for G4Box shape, and the surface
45// is fixed on -Z plane.
46//
47// Surface is defined at the -Z surface.
48// Direction -Z +Z
49// 0 IN || OUT ->|<- |
50// 1 IN ->| |
51// 2 OUT |<- |
52//
53// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
54//
55// 18-Nov-2005 T.Aso, To use always positive value for anglefactor.
56// 29-Mar-2007 T.Aso, Bug fix for momentum direction at outgoing flux.
57// 2010-07-22 Introduce Unit specification.
58// 2010-07-22 Add weighted and divideByAre options
59// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
60// vs. cell flux * track weight (Makoto Asai)
62
64 G4int depth)
66 , HCID(-1)
67 , fDirection(direction)
68 , EvtMap(0)
69 , weighted(true)
70 , divideByArea(true)
71{
73 SetUnit("percm2");
74}
75
77 const G4String& unit, G4int depth)
79 , HCID(-1)
80 , fDirection(direction)
81 , EvtMap(0)
82 , weighted(true)
83 , divideByArea(true)
84{
86 SetUnit(unit);
87}
88
90
92{
93 G4StepPoint* preStep = aStep->GetPreStepPoint();
94 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
95 G4VPVParameterisation* physParam = physVol->GetParameterisation();
96 G4VSolid* solid = 0;
97 if(physParam)
98 { // for parameterized volume
99 G4int idx =
101 ->GetReplicaNumber(indexDepth);
102 solid = physParam->ComputeSolid(idx, physVol);
103 solid->ComputeDimensions(physParam, idx, physVol);
104 }
105 else
106 { // for ordinary volume
107 solid = physVol->GetLogicalVolume()->GetSolid();
108 }
109
110 G4Box* boxSolid = (G4Box*) (solid);
111
112 G4int dirFlag = IsSelectedSurface(aStep, boxSolid);
113 if(dirFlag > 0)
114 {
115 if(fDirection == fFlux_InOut || fDirection == dirFlag)
116 {
117 G4StepPoint* thisStep = 0;
118 if(dirFlag == fFlux_In)
119 {
120 thisStep = preStep;
121 }
122 else if(dirFlag == fFlux_Out)
123 {
124 thisStep = aStep->GetPostStepPoint();
125 }
126 else
127 {
128 return FALSE;
129 }
130
131 G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
132 G4ThreeVector pdirection = thisStep->GetMomentumDirection();
133 G4ThreeVector localdir =
134 theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
135 //
136 G4double angleFactor = localdir.z();
137 if(angleFactor < 0)
138 angleFactor *= -1.;
139 G4double flux = 1.0;
140 if(weighted)
141 flux *= preStep->GetWeight(); // Current (Particle Weight)
142 //
143 G4double square =
144 4. * boxSolid->GetXHalfLength() * boxSolid->GetYHalfLength();
145 //
146 flux = flux / angleFactor; // Flux with angle.
147 if(divideByArea)
148 flux /= square;
149 //
150 G4int index = GetIndex(aStep);
151 EvtMap->add(index, flux);
152
153 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
154 {
155 auto filler = G4VScoreHistFiller::Instance();
156 if(!filler)
157 {
158 G4Exception("G4PSFlatSurfaceFlux::ProcessHits", "SCORER0123",
160 "G4TScoreHistFiller is not instantiated!! Histogram is "
161 "not filled.");
162 }
163 else
164 {
165 filler->FillH1(hitIDMap[index], preStep->GetKineticEnergy(), flux);
166 }
167 }
168 }
169 }
170#ifdef debug
171 G4cout << " PASSED vol " << index << " trk " << trkid << " len "
172 << fFlatSurfaceFlux << G4endl;
173#endif
174
175 return TRUE;
176}
177
179{
180 G4TouchableHandle theTouchable =
184
186 {
187 // Entering Geometry
188 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition();
189 G4ThreeVector localpos1 =
190 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
191 if(std::fabs(localpos1.z() + boxSolid->GetZHalfLength()) < kCarTolerance)
192 {
193 return fFlux_In;
194 }
195 }
196
198 {
199 // Exiting Geometry
200 G4ThreeVector stppos2 = aStep->GetPostStepPoint()->GetPosition();
201 G4ThreeVector localpos2 =
202 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
203 if(std::fabs(localpos2.z() + boxSolid->GetZHalfLength()) < kCarTolerance)
204 {
205 return fFlux_Out;
206 }
207 }
208
209 return -1;
210}
211
213{
215 GetName());
216 if(HCID < 0)
219}
220
222
224
226
228{
229 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
230 G4cout << " PrimitiveScorer" << GetName() << G4endl;
231 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
232 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
233 for(; itr != EvtMap->GetMap()->end(); itr++)
234 {
235 G4cout << " copy no.: " << itr->first
236 << " flux : " << *(itr->second) / GetUnitValue() << " ["
237 << GetUnit() << "]" << G4endl;
238 }
239}
240
242{
243 if(divideByArea)
244 {
245 CheckAndSetUnit(unit, "Per Unit Surface");
246 }
247 else
248 {
249 if(unit == "")
250 {
251 unitName = unit;
252 unitValue = 1.0;
253 }
254 else
255 {
256 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
257 GetUnit() + "] ) for " + GetName();
258 G4Exception("G4PSFlatSurfaceFlux::SetUnit", "DetPS0008", JustWarning,
259 msg);
260 }
261 }
262}
263
265{
266 // Per Unit Surface
267 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
268 (1. / cm2));
269 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
270 (1. / mm2));
271 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
272}
const G4double kCarTolerance
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fFlux_InOut
@ fFlux_Out
@ fFlux_In
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
double z() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
Definition: G4Box.hh:56
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4VSolid * GetSolid() const
const G4AffineTransform & GetTopTransform() const
G4THitsMap< G4double > * EvtMap
G4int IsSelectedSurface(G4Step *, G4Box *)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
G4PSFlatSurfaceFlux(G4String name, G4int direction, G4int depth=0)
virtual void EndOfEvent(G4HCofThisEvent *)
virtual void SetUnit(const G4String &unit)
virtual void DefineUnitAndCategory()
virtual void Initialize(G4HCofThisEvent *)
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
const G4TouchableHandle & GetTouchableHandle() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
G4MultiFunctionalDetector * GetMultiFunctionalDetector() const
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
void CheckAndSetUnit(const G4String &unit, const G4String &category)
G4double GetUnitValue() const
static G4VScoreHistFiller * Instance()
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
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
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:82
const char * name(G4int ptype)