Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PSCylinderSurfaceFlux.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// // G4PSCylinderSurfaceFlux
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
41// ////////////////////////////////////////////////////////////////////////////////
42// (Description)
43// This is a primitive scorer class for scoring Surface Flux.
44// Current version assumes only for G4Tubs shape, and the surface
45// is fixed on inner plane of the tube.
46//
47// Surface is defined at the innner surface of the tube.
48// Direction R R+dR
49// 0 IN || OUT ->|<- |
50// 1 IN ->| |
51// 2 OUT |<- |
52//
53// Created: 2007-03-29 Tsukasa ASO
54// 2010-07-22 Introduce Unit specification.
55// 2010-07-22 Add weighted and divideByArea options
56// 2011-02-21 Get correct momentum direction in Flux_Out.
57// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
58// vs. surface flux * track weight (y) (Makoto Asai)
59//
61
63 G4int depth)
65 , HCID(-1)
66 , fDirection(direction)
67 , EvtMap(0)
68 , weighted(true)
69 , divideByArea(true)
70{
72 SetUnit("percm2");
73}
74
76 const G4String& unit,
77 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 G4VSolid* solid = ComputeCurrentSolid(aStep);
95 assert(dynamic_cast<G4Tubs*>(solid));
96
97 G4Tubs* tubsSolid = static_cast<G4Tubs*>(solid);
98
99 G4int dirFlag = IsSelectedSurface(aStep, tubsSolid);
100
101 if(dirFlag > 0)
102 {
103 if(fDirection == fFlux_InOut || dirFlag == fDirection)
104 {
105 G4StepPoint* thisStep = 0;
106 if(dirFlag == fFlux_In)
107 {
108 thisStep = preStep;
109 }
110 else if(dirFlag == fFlux_Out)
111 {
112 thisStep = aStep->GetPostStepPoint();
113 }
114 else
115 {
116 return FALSE;
117 }
118
119 G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
120 G4ThreeVector pdirection = thisStep->GetMomentumDirection();
121 G4ThreeVector localdir =
122 theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
123 G4ThreeVector position = thisStep->GetPosition();
124 G4ThreeVector localpos =
126 G4double angleFactor =
127 (localdir.x() * localpos.x() + localdir.y() * localpos.y()) /
128 std::sqrt(localdir.x() * localdir.x() + localdir.y() * localdir.y() +
129 localdir.z() * localdir.z()) /
130 std::sqrt(localpos.x() * localpos.x() + localpos.y() * localpos.y());
131
132 if(angleFactor < 0)
133 angleFactor *= -1.;
134 G4double square = 2. * tubsSolid->GetZHalfLength() *
135 tubsSolid->GetInnerRadius() *
136 tubsSolid->GetDeltaPhiAngle() / radian;
137
138 G4double flux = 1.0;
139 if(weighted)
140 flux *= preStep->GetWeight();
141 // Current (Particle Weight)
142
143 flux = flux / angleFactor;
144 if(divideByArea)
145 flux /= square;
146 // Flux with angle.
147 G4int index = GetIndex(aStep);
148 EvtMap->add(index, flux);
149
150 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
151 {
152 auto filler = G4VScoreHistFiller::Instance();
153 if(!filler)
154 {
155 G4Exception("G4PSCylinderSurfaceFlux::ProcessHits", "SCORER0123",
157 "G4TScoreHistFiller is not instantiated!! Histogram is "
158 "not filled.");
159 }
160 else
161 {
162 filler->FillH1(hitIDMap[index], thisStep->GetKineticEnergy(), flux);
163 }
164 }
165
166 return TRUE;
167 }
168 else
169 {
170 return FALSE;
171 }
172 }
173 else
174 {
175 return FALSE;
176 }
177}
178
180 G4Tubs* tubsSolid)
181{
182 G4TouchableHandle theTouchable =
186
188 {
189 // Entering Geometry
190 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition();
191 G4ThreeVector localpos1 =
192 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
193 if(std::fabs(localpos1.z()) > tubsSolid->GetZHalfLength())
194 return -1;
195 // if(std::fabs( localpos1.x()*localpos1.x()+localpos1.y()*localpos1.y()
196 // - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius()))
197 // <kCarTolerance ){
198 G4double localR2 =
199 localpos1.x() * localpos1.x() + localpos1.y() * localpos1.y();
200 G4double InsideRadius = tubsSolid->GetInnerRadius();
201 if(localR2 >
202 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
203 localR2 <
204 (InsideRadius + kCarTolerance) * (InsideRadius + kCarTolerance))
205 {
206 return fFlux_In;
207 }
208 }
209
211 {
212 // Exiting Geometry
213 G4ThreeVector stppos2 = aStep->GetPostStepPoint()->GetPosition();
214 G4ThreeVector localpos2 =
215 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
216 if(std::fabs(localpos2.z()) > tubsSolid->GetZHalfLength())
217 return -1;
218 // if(std::fabs( localpos2.x()*localpos2.x()+localpos2.y()*localpos2.y()
219 // - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius()))
220 // <kCarTolerance ){
221 G4double localR2 =
222 localpos2.x() * localpos2.x() + localpos2.y() * localpos2.y();
223 G4double InsideRadius = tubsSolid->GetInnerRadius();
224 if(localR2 >
225 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
226 localR2 <
227 (InsideRadius + kCarTolerance) * (InsideRadius + kCarTolerance))
228 {
229 return fFlux_Out;
230 }
231 }
232
233 return -1;
234}
235
237{
239 GetName());
240 if(HCID < 0)
243}
244
246
248
250
252{
253 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
254 G4cout << " PrimitiveScorer" << GetName() << G4endl;
255 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
256 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
257 for(; itr != EvtMap->GetMap()->end(); itr++)
258 {
259 G4cout << " copy no.: " << itr->first
260 << " flux : " << *(itr->second) / GetUnitValue() << " ["
261 << GetUnit() << "]" << G4endl;
262 }
263}
264
266{
267 if(divideByArea)
268 {
269 CheckAndSetUnit(unit, "Per Unit Surface");
270 }
271 else
272 {
273 if(unit == "")
274 {
275 unitName = unit;
276 unitValue = 1.0;
277 }
278 else
279 {
280 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
281 GetUnit() + "] ) for " + GetName();
282 G4Exception("G4PSCylinderSurfaceFlux::SetUnit", "DetPS0003", JustWarning,
283 msg);
284 }
285 }
286}
287
289{
290 // Per Unit Surface
291 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
292 (1. / cm2));
293 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
294 (1. / mm2));
295 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
296}
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
static constexpr double radian
Definition: G4SIunits.hh:122
@ 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
double x() const
double y() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
const G4AffineTransform & GetTopTransform() const
virtual void EndOfEvent(G4HCofThisEvent *)
virtual void Initialize(G4HCofThisEvent *)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
G4THitsMap< G4double > * EvtMap
G4PSCylinderSurfaceFlux(G4String name, G4int direction, G4int depth=0)
virtual void SetUnit(const G4String &unit)
G4int IsSelectedSurface(G4Step *, G4Tubs *)
G4StepStatus GetStepStatus() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
Definition: G4Tubs.hh:75
G4double GetZHalfLength() const
G4double GetInnerRadius() const
G4double GetDeltaPhiAngle() const
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
G4VSolid * ComputeCurrentSolid(G4Step *aStep)
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
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:82
const char * name(G4int ptype)