Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
28//
29// G4PSFlatSurfaceCurrent
31
32#include "G4SystemOfUnits.hh"
33#include "G4StepStatus.hh"
34#include "G4Track.hh"
35#include "G4VSolid.hh"
36#include "G4VPhysicalVolume.hh"
38#include "G4UnitsTable.hh"
40////////////////////////////////////////////////////////////////////////////////
41// (Description)
42// This is a primitive scorer class for scoring only Surface Flux.
43// Current version assumes only for G4Box shape.
44//
45// Surface is defined at the -Z surface.
46// Direction -Z +Z
47// 0 IN || OUT ->|<- |
48// 1 IN ->| |
49// 2 OUT |<- |
50//
51// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
52// 17-Nov-2005 T.Aso, Bug fix for area definition.
53// 31-Mar-2007 T.Aso, Add option for normalizing by the area.
54// 2010-07-22 Introduce Unit specification.
55//
56///////////////////////////////////////////////////////////////////////////////
57
58
60 G4int direction, G4int depth)
61 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
62 weighted(true),divideByArea(true)
63{
65 SetUnit("percm2");
66}
67
69 G4int direction,
70 const G4String& unit,
71 G4int depth)
72 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
73 weighted(true),divideByArea(true)
74{
76 SetUnit(unit);
77}
78
80{;}
81
83{
84 G4StepPoint* preStep = aStep->GetPreStepPoint();
85 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
86 G4VPVParameterisation* physParam = physVol->GetParameterisation();
87 G4VSolid * solid = 0;
88 if(physParam)
89 { // for parameterized volume
91 ->GetReplicaNumber(indexDepth);
92 solid = physParam->ComputeSolid(idx, physVol);
93 solid->ComputeDimensions(physParam,idx,physVol);
94 }
95 else
96 { // for ordinary volume
97 solid = physVol->GetLogicalVolume()->GetSolid();
98 }
99
100 G4Box* boxSolid = (G4Box*)(solid);
101
102 G4int dirFlag =IsSelectedSurface(aStep,boxSolid);
103 if ( dirFlag > 0 ) {
104 if ( fDirection == fCurrent_InOut || fDirection == dirFlag ){
105 G4int index = GetIndex(aStep);
106 G4TouchableHandle theTouchable = preStep->GetTouchableHandle();
107 G4double current = 1.0;
108 if ( weighted ) current=preStep->GetWeight(); // Current (Particle Weight)
109 if ( divideByArea ){
110 G4double square = 4.*boxSolid->GetXHalfLength()*boxSolid->GetYHalfLength();
111 current = current/square; // Normalized by Area
112 }
113 EvtMap->add(index,current);
114 }
115 }
116
117 return TRUE;
118}
119
121
122 G4TouchableHandle theTouchable =
125
126 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
127 // Entering Geometry
128 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
129 G4ThreeVector localpos1 =
130 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
131 if(std::fabs( localpos1.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
132 return fCurrent_In;
133 }
134 }
135
136 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
137 // Exiting Geometry
138 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
139 G4ThreeVector localpos2 =
140 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
141 if(std::fabs( localpos2.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
142 return fCurrent_Out;
143 }
144 }
145
146 return -1;
147}
148
150{
151 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
152 if ( HCID < 0 ) HCID = GetCollectionID(0);
153 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
154}
155
157{;}
158
160 EvtMap->clear();
161}
162
164{;}
165
167{
168 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
169 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
170 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
171 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
172 for(; itr != EvtMap->GetMap()->end(); itr++) {
173 G4cout << " copy no.: " << itr->first << " current : " ;
174 if ( divideByArea ) {
175 G4cout << *(itr->second)/GetUnitValue()
176 << " ["<<GetUnit()<<"]";
177 }else {
178 G4cout << *(itr->second)/GetUnitValue() << " [tracks]";
179 }
180 G4cout << G4endl;
181 }
182}
183
185{
186 if ( divideByArea ) {
187 CheckAndSetUnit(unit,"Per Unit Surface");
188 } else {
189 if (unit == "" ){
190 unitName = unit;
191 unitValue = 1.0;
192 }else{
193 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
194 G4Exception("G4PSFlatSurfaceCurrent::SetUnit","DetPS0007",JustWarning,msg);
195 }
196 }
197}
198
200 // Per Unit Surface
201 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
202 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
203 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
204}
205
@ JustWarning
@ fCurrent_Out
@ fCurrent_In
@ fCurrent_InOut
@ fGeomBoundary
Definition: G4StepStatus.hh:54
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double z() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
Definition: G4Box.hh:55
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 GetWeight() const
Definition: G4Step.hh:78
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:138
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
G4int entries() const
Definition: G4THitsMap.hh:79
void clear()
Definition: G4THitsMap.hh:209
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0
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
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:86
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define TRUE
Definition: globals.hh:55