Geant4-11
G4PSCylinderSurfaceCurrent.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// G4PSCylinderSurfaceCurrent
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 Current.
44// Current version assumes only for G4Tubs shape.
45//
46// Surface is defined at the inner surface of the tube.
47// Direction R R+dR
48// 0 IN || OUT ->|<- |
49// 1 IN ->| |
50// 2 OUT |<- |
51//
52// Created: 2007-03-21 Tsukasa ASO
53// 2010-07-22 Introduce Unit specification.
54// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
55// vs. surface current * track weight (y) (Makoto Asai)
56//
58
60 G4int direction,
61 G4int depth)
63 , HCID(-1)
64 , fDirection(direction)
65 , EvtMap(0)
66 , weighted(true)
67 , divideByArea(true)
68{
70 SetUnit("percm2");
71}
72
74 G4int direction,
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
92{
93 G4StepPoint* preStep = aStep->GetPreStepPoint();
94 G4VSolid* solid = ComputeCurrentSolid(aStep);
95 G4Tubs* tubsSolid = static_cast<G4Tubs*>(solid);
96
97 G4int dirFlag = IsSelectedSurface(aStep, tubsSolid);
98 // G4cout << " pos " << preStep->GetPosition() <<" dirFlag " << G4endl;
99 if(dirFlag > 0)
100 {
101 if(fDirection == fCurrent_InOut || fDirection == dirFlag)
102 {
103 G4TouchableHandle theTouchable = preStep->GetTouchableHandle();
104 //
105 G4double current = 1.0;
106 if(weighted)
107 current = preStep->GetWeight(); // Current (Particle Weight)
108 //
109 if(divideByArea)
110 {
111 G4double square = 2. * tubsSolid->GetZHalfLength() *
112 tubsSolid->GetInnerRadius() *
113 tubsSolid->GetDeltaPhiAngle() / radian;
114 current = current / square; // Current normalized by Area
115 }
116
117 G4int index = GetIndex(aStep);
118 EvtMap->add(index, current);
119
120 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
121 {
122 auto filler = G4VScoreHistFiller::Instance();
123 if(!filler)
124 {
125 G4Exception("G4PSCylinderSurfaceCurrent::ProcessHits", "SCORER0123",
127 "G4TScoreHistFiller is not instantiated!! Histogram is "
128 "not filled.");
129 }
130 else
131 {
132 filler->FillH1(hitIDMap[index], preStep->GetKineticEnergy(), current);
133 }
134 }
135 }
136 }
137
138 return TRUE;
139}
140
142 G4Tubs* tubsSolid)
143{
144 G4TouchableHandle theTouchable =
148
150 {
151 // Entering Geometry
152 G4ThreeVector stppos1 = aStep->GetPreStepPoint()->GetPosition();
153 G4ThreeVector localpos1 =
154 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
155 if(std::fabs(localpos1.z()) > tubsSolid->GetZHalfLength())
156 return -1;
157 G4double localR2 =
158 localpos1.x() * localpos1.x() + localpos1.y() * localpos1.y();
159 G4double InsideRadius = tubsSolid->GetInnerRadius();
160 if(localR2 >
161 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
162 localR2 <
163 (InsideRadius + kCarTolerance) * (InsideRadius + 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()) > tubsSolid->GetZHalfLength())
176 return -1;
177 G4double localR2 =
178 localpos2.x() * localpos2.x() + localpos2.y() * localpos2.y();
179 G4double InsideRadius = tubsSolid->GetInnerRadius();
180 if(localR2 >
181 (InsideRadius - kCarTolerance) * (InsideRadius - kCarTolerance) &&
182 localR2 <
183 (InsideRadius + kCarTolerance) * (InsideRadius + kCarTolerance))
184 {
185 return fCurrent_Out;
186 }
187 }
188
189 return -1;
190}
191
193{
195 if(HCID < 0)
198}
199
201
203
205
207{
208 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
209 G4cout << " PrimitiveScorer " << GetName() << G4endl;
210 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
211 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
212 for(; itr != EvtMap->GetMap()->end(); itr++)
213 {
214 G4cout << " copy no.: " << itr->first << " current : ";
215 if(divideByArea)
216 {
217 G4cout << *(itr->second) / GetUnitValue() << " [" << GetUnit() << "]";
218 }
219 else
220 {
221 G4cout << *(itr->second) << " [tracks]";
222 }
223 G4cout << G4endl;
224 }
225}
226
228{
229 if(divideByArea)
230 {
231 CheckAndSetUnit(unit, "Per Unit Surface");
232 }
233 else
234 {
235 if(unit == "")
236 {
237 unitName = unit;
238 unitValue = 1.0;
239 }
240 else
241 {
242 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
243 GetUnit() + "] ) for " + GetName();
244 G4Exception("G4PSCylinderSurfaceCurrent::SetUnit", "DetPS0002",
245 JustWarning, msg);
246 }
247 }
248}
249
251{
252 // Per Unit Surface
253 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface",
254 (1. / cm2));
255 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface",
256 (1. / mm2));
257 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2));
258}
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
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
double z() const
double x() const
double y() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
const G4AffineTransform & GetTopTransform() const
G4PSCylinderSurfaceCurrent(G4String name, G4int direction, G4int depth=0)
virtual void SetUnit(const G4String &unit)
G4int IsSelectedSurface(G4Step *, G4Tubs *)
virtual void Initialize(G4HCofThisEvent *)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void EndOfEvent(G4HCofThisEvent *)
G4StepStatus GetStepStatus() const
const G4ThreeVector & GetPosition() 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
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)