Geant4-11
G4PSFlatSurfaceCurrent.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// G4PSFlatSurfaceCurrent
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 only Surface Flux.
44// Current version assumes only for G4Box shape.
45//
46// Surface is defined at the -Z surface.
47// Direction -Z +Z
48// 0 IN || OUT ->|<- |
49// 1 IN ->| |
50// 2 OUT |<- |
51//
52// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
53// 17-Nov-2005 T.Aso, Bug fix for area definition.
54// 31-Mar-2007 T.Aso, Add option for normalizing by the area.
55// 2010-07-22 Introduce Unit specification.
56// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
57// vs. Surface Current * track weight (y) (Makoto Asai)
58//
60
62 G4int depth)
64 , HCID(-1)
65 , fDirection(direction)
66 , EvtMap(0)
67 , weighted(true)
68 , divideByArea(true)
69{
71 SetUnit("percm2");
72}
73
75 const G4String& unit,
76 G4int depth)
78 , HCID(-1)
79 , fDirection(direction)
80 , EvtMap(0)
81 , weighted(true)
82 , divideByArea(true)
83{
85 SetUnit(unit);
86}
87
89
91{
92 G4StepPoint* preStep = aStep->GetPreStepPoint();
93 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
94 G4VPVParameterisation* physParam = physVol->GetParameterisation();
95 G4VSolid* solid = 0;
96 if(physParam)
97 { // for parameterized volume
98 G4int idx =
100 ->GetReplicaNumber(indexDepth);
101 solid = physParam->ComputeSolid(idx, physVol);
102 solid->ComputeDimensions(physParam, idx, physVol);
103 }
104 else
105 { // for ordinary volume
106 solid = physVol->GetLogicalVolume()->GetSolid();
107 }
108
109 G4Box* boxSolid = (G4Box*) (solid);
110
111 G4int dirFlag = IsSelectedSurface(aStep, boxSolid);
112 if(dirFlag > 0)
113 {
114 if(fDirection == fCurrent_InOut || fDirection == dirFlag)
115 {
116 G4int index = GetIndex(aStep);
117 G4TouchableHandle theTouchable = preStep->GetTouchableHandle();
118 G4double current = 1.0;
119 if(weighted)
120 current = preStep->GetWeight(); // Current (Particle Weight)
121 if(divideByArea)
122 {
123 G4double square =
124 4. * boxSolid->GetXHalfLength() * boxSolid->GetYHalfLength();
125 current = current / square; // Normalized by Area
126 }
127 EvtMap->add(index, current);
128
129 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
130 {
131 auto filler = G4VScoreHistFiller::Instance();
132 if(!filler)
133 {
134 G4Exception("G4PSFlatSurfaceCurrent::ProcessHits", "SCORER0123",
136 "G4TScoreHistFiller is not instantiated!! Histogram is "
137 "not filled.");
138 }
139 else
140 {
141 filler->FillH1(hitIDMap[index], preStep->GetKineticEnergy(), current);
142 }
143 }
144 }
145 }
146
147 return TRUE;
148}
149
151{
152 G4TouchableHandle theTouchable =
156
158 {
159 // Entering Geometry
160 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition();
161 G4ThreeVector localpos1 =
162 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
163 if(std::fabs(localpos1.z() + boxSolid->GetZHalfLength()) < kCarTolerance)
164 {
165 return fCurrent_In;
166 }
167 }
168
170 {
171 // Exiting Geometry
172 G4ThreeVector stppos2 = aStep->GetPostStepPoint()->GetPosition();
173 G4ThreeVector localpos2 =
174 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
175 if(std::fabs(localpos2.z() + boxSolid->GetZHalfLength()) < kCarTolerance)
176 {
177 return fCurrent_Out;
178 }
179 }
180
181 return -1;
182}
183
185{
187 if(HCID < 0)
190}
191
193
195
197
199{
200 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
201 G4cout << " PrimitiveScorer " << GetName() << G4endl;
202 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
203 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
204 for(; itr != EvtMap->GetMap()->end(); itr++)
205 {
206 G4cout << " copy no.: " << itr->first << " current : ";
207 if(divideByArea)
208 {
209 G4cout << *(itr->second) / GetUnitValue() << " [" << GetUnit() << "]";
210 }
211 else
212 {
213 G4cout << *(itr->second) / GetUnitValue() << " [tracks]";
214 }
215 G4cout << G4endl;
216 }
217}
218
220{
221 if(divideByArea)
222 {
223 CheckAndSetUnit(unit, "Per Unit Surface");
224 }
225 else
226 {
227 if(unit == "")
228 {
229 unitName = unit;
230 unitValue = 1.0;
231 }
232 else
233 {
234 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
235 GetUnit() + "] ) for " + GetName();
236 G4Exception("G4PSFlatSurfaceCurrent::SetUnit", "DetPS0007", JustWarning,
237 msg);
238 }
239 }
240}
241
243{
244 // Per Unit Surface
245 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
246 (1. / cm2));
247 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
248 (1. / mm2));
249 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
250}
const G4double kCarTolerance
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fCurrent_Out
@ fCurrent_In
@ fCurrent_InOut
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
double z() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) 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 *)
G4THitsMap< G4double > * EvtMap
virtual void EndOfEvent(G4HCofThisEvent *)
virtual void SetUnit(const G4String &unit)
G4PSFlatSurfaceCurrent(G4String name, G4int direction, G4int depth=0)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void Initialize(G4HCofThisEvent *)
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetPosition() 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
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)