Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
41////////////////////////////////////////////////////////////////////////////////
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)
61///////////////////////////////////////////////////////////////////////////////
62
64 G4int direction, G4int depth)
65 : G4VPrimitivePlotter(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
66 weighted(true),divideByArea(true)
67{
69 SetUnit("percm2");
70}
71
73 G4int direction,
74 const G4String& unit,
75 G4int depth)
76 : G4VPrimitivePlotter(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
77 weighted(true),divideByArea(true)
78{
80 SetUnit(unit);
81}
82
84{;}
85
87{
88 G4StepPoint* preStep = aStep->GetPreStepPoint();
89 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
90 G4VPVParameterisation* physParam = physVol->GetParameterisation();
91 G4VSolid * solid = 0;
92 if(physParam)
93 { // for parameterized volume
95 ->GetReplicaNumber(indexDepth);
96 solid = physParam->ComputeSolid(idx, physVol);
97 solid->ComputeDimensions(physParam,idx,physVol);
98 }
99 else
100 { // for ordinary volume
101 solid = physVol->GetLogicalVolume()->GetSolid();
102 }
103
104 G4Box* boxSolid = (G4Box*)(solid);
105
106 G4int dirFlag =IsSelectedSurface(aStep,boxSolid);
107 if ( dirFlag > 0 ) {
108 if ( fDirection == fFlux_InOut || fDirection == dirFlag ){
109
110 G4StepPoint* thisStep=0;
111 if ( dirFlag == fFlux_In ){
112 thisStep = preStep;
113 }else if ( dirFlag == fFlux_Out ){
114 thisStep = aStep->GetPostStepPoint();
115 }else{
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 //
124 G4double angleFactor = localdir.z();
125 if ( angleFactor < 0 ) angleFactor *= -1.;
126 G4double flux = 1.0;
127 if ( weighted ) flux *=preStep->GetWeight(); // Current (Particle Weight)
128 //
129 G4double square = 4.*boxSolid->GetXHalfLength()*boxSolid->GetYHalfLength();
130 //
131 flux = flux/angleFactor; // Flux with angle.
132 if ( divideByArea ) flux /= square;
133 //
134 G4int index = GetIndex(aStep);
135 EvtMap->add(index,flux);
136
137 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
138 {
139 auto filler = G4VScoreHistFiller::Instance();
140 if(!filler)
141 {
142 G4Exception("G4PSFlatSurfaceFlux::ProcessHits","SCORER0123",JustWarning,
143 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
144 }
145 else
146 {
147 filler->FillH1(hitIDMap[index],preStep->GetKineticEnergy(),flux);
148 }
149 }
150
151 }
152 }
153#ifdef debug
154 G4cout << " PASSED vol "
155 << index << " trk "<<trkid<<" len " << fFlatSurfaceFlux<<G4endl;
156#endif
157
158 return TRUE;
159}
160
162
163 G4TouchableHandle theTouchable =
166
167 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
168 // Entering Geometry
169 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
170 G4ThreeVector localpos1 =
171 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
172 if(std::fabs( localpos1.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
173 return fFlux_In;
174 }
175 }
176
177 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
178 // Exiting Geometry
179 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
180 G4ThreeVector localpos2 =
181 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
182 if(std::fabs( localpos2.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
183 return fFlux_Out;
184 }
185 }
186
187 return -1;
188}
189
191{
193 GetName());
194 if ( HCID < 0 ) HCID = GetCollectionID(0);
195 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
196}
197
199{;}
200
202 EvtMap->clear();
203}
204
206{;}
207
209{
210 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
211 G4cout << " PrimitiveScorer" << GetName() <<G4endl;
212 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
213 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
214 for(; itr != EvtMap->GetMap()->end(); itr++) {
215 G4cout << " copy no.: " << itr->first
216 << " flux : " << *(itr->second)/GetUnitValue()
217 << " [" << GetUnit() <<"]"
218 << G4endl;
219 }
220}
221
223{
224 if ( divideByArea ) {
225 CheckAndSetUnit(unit,"Per Unit Surface");
226 } else {
227 if (unit == "" ){
228 unitName = unit;
229 unitValue = 1.0;
230 }else{
231 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
232 G4Exception("G4PSFlatSurfaceFlux::SetUnit","DetPS0008",JustWarning,msg);
233 }
234 }
235}
236
238 // Per Unit Surface
239 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
240 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
241 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
242}
243
244
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
@ 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
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:125
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
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:83