Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
41////////////////////////////////////////////////////////////////////////////////
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//
59///////////////////////////////////////////////////////////////////////////////
60
61
63 G4int direction, G4int depth)
64 :G4VPrimitivePlotter(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
65 weighted(true),divideByArea(true)
66{
68 SetUnit("percm2");
69}
70
72 G4int direction,
73 const G4String& unit,
74 G4int depth)
75 :G4VPrimitivePlotter(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
76 weighted(true),divideByArea(true)
77{
79 SetUnit(unit);
80}
81
83{;}
84
86{
87 G4StepPoint* preStep = aStep->GetPreStepPoint();
88 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
89 G4VPVParameterisation* physParam = physVol->GetParameterisation();
90 G4VSolid * solid = 0;
91 if(physParam)
92 { // for parameterized volume
94 ->GetReplicaNumber(indexDepth);
95 solid = physParam->ComputeSolid(idx, physVol);
96 solid->ComputeDimensions(physParam,idx,physVol);
97 }
98 else
99 { // for ordinary volume
100 solid = physVol->GetLogicalVolume()->GetSolid();
101 }
102
103 G4Box* boxSolid = (G4Box*)(solid);
104
105 G4int dirFlag =IsSelectedSurface(aStep,boxSolid);
106 if ( dirFlag > 0 ) {
107 if ( fDirection == fCurrent_InOut || fDirection == dirFlag ){
108 G4int index = GetIndex(aStep);
109 G4TouchableHandle theTouchable = preStep->GetTouchableHandle();
110 G4double current = 1.0;
111 if ( weighted ) current=preStep->GetWeight(); // Current (Particle Weight)
112 if ( divideByArea ){
113 G4double square = 4.*boxSolid->GetXHalfLength()*boxSolid->GetYHalfLength();
114 current = current/square; // Normalized by Area
115 }
116 EvtMap->add(index,current);
117
118 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
119 {
120 auto filler = G4VScoreHistFiller::Instance();
121 if(!filler)
122 {
123 G4Exception("G4PSFlatSurfaceCurrent::ProcessHits","SCORER0123",JustWarning,
124 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
125 }
126 else
127 {
128 filler->FillH1(hitIDMap[index],preStep->GetKineticEnergy(),current);
129 }
130 }
131
132 }
133 }
134
135 return TRUE;
136}
137
139
140 G4TouchableHandle theTouchable =
143
144 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
145 // Entering Geometry
146 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
147 G4ThreeVector localpos1 =
148 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
149 if(std::fabs( localpos1.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
150 return fCurrent_In;
151 }
152 }
153
154 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
155 // Exiting Geometry
156 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
157 G4ThreeVector localpos2 =
158 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
159 if(std::fabs( localpos2.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
160 return fCurrent_Out;
161 }
162 }
163
164 return -1;
165}
166
168{
169 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
170 if ( HCID < 0 ) HCID = GetCollectionID(0);
171 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
172}
173
175{;}
176
178 EvtMap->clear();
179}
180
182{;}
183
185{
186 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
187 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
188 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
189 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
190 for(; itr != EvtMap->GetMap()->end(); itr++) {
191 G4cout << " copy no.: " << itr->first << " current : " ;
192 if ( divideByArea ) {
193 G4cout << *(itr->second)/GetUnitValue()
194 << " ["<<GetUnit()<<"]";
195 }else {
196 G4cout << *(itr->second)/GetUnitValue() << " [tracks]";
197 }
198 G4cout << G4endl;
199 }
200}
201
203{
204 if ( divideByArea ) {
205 CheckAndSetUnit(unit,"Per Unit Surface");
206 } else {
207 if (unit == "" ){
208 unitName = unit;
209 unitValue = 1.0;
210 }else{
211 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
212 G4Exception("G4PSFlatSurfaceCurrent::SetUnit","DetPS0007",JustWarning,msg);
213 }
214 }
215}
216
218 // Per Unit Surface
219 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
220 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
221 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
222}
223
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
@ 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 *)
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: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