Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // $Id: G4PSFlatSurfaceFlux.cc 81087 2014-05-20 15:44:27Z gcosmo $
28 //
29 // G4PSFlatSurfaceFlux
30 #include "G4PSFlatSurfaceFlux.hh"
31 
32 #include "G4SystemOfUnits.hh"
33 #include "G4StepStatus.hh"
34 #include "G4Track.hh"
35 #include "G4VSolid.hh"
36 #include "G4VPhysicalVolume.hh"
37 #include "G4VPVParameterisation.hh"
38 #include "G4UnitsTable.hh"
39 #include "G4GeometryTolerance.hh"
41 // (Description)
42 // This is a primitive scorer class for scoring Surface Flux.
43 // Current version assumes only for G4Box shape, and the surface
44 // is fixed on -Z plane.
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 //
54 // 18-Nov-2005 T.Aso, To use always positive value for anglefactor.
55 // 29-Mar-2007 T.Aso, Bug fix for momentum direction at outgoing flux.
56 // 2010-07-22 Introduce Unit specification.
57 // 2010-07-22 Add weighted and divideByAre options
59 
61  G4int direction, G4int depth)
62  : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
63  weighted(true),divideByArea(true)
64 {
66  SetUnit("percm2");
67 }
68 
70  G4int direction,
71  const G4String& unit,
72  G4int depth)
73  : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
74  weighted(true),divideByArea(true)
75 {
77  SetUnit(unit);
78 }
79 
81 {;}
82 
84 {
85  G4StepPoint* preStep = aStep->GetPreStepPoint();
86  G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
87  G4VPVParameterisation* physParam = physVol->GetParameterisation();
88  G4VSolid * solid = 0;
89  if(physParam)
90  { // for parameterized volume
92  ->GetReplicaNumber(indexDepth);
93  solid = physParam->ComputeSolid(idx, physVol);
94  solid->ComputeDimensions(physParam,idx,physVol);
95  }
96  else
97  { // for ordinary volume
98  solid = physVol->GetLogicalVolume()->GetSolid();
99  }
100 
101  G4Box* boxSolid = (G4Box*)(solid);
102 
103  G4int dirFlag =IsSelectedSurface(aStep,boxSolid);
104  if ( dirFlag > 0 ) {
105  if ( fDirection == fFlux_InOut || fDirection == dirFlag ){
106 
107  G4StepPoint* thisStep=0;
108  if ( dirFlag == fFlux_In ){
109  thisStep = preStep;
110  }else if ( dirFlag == fFlux_Out ){
111  thisStep = aStep->GetPostStepPoint();
112  }else{
113  return FALSE;
114  }
115 
116  G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
117  G4ThreeVector pdirection = thisStep->GetMomentumDirection();
118  G4ThreeVector localdir =
119  theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
120  //
121  G4double angleFactor = localdir.z();
122  if ( angleFactor < 0 ) angleFactor *= -1.;
123  G4double flux = 1.0;
124  if ( weighted ) flux *=preStep->GetWeight(); // Current (Particle Weight)
125  //
126  G4double square = 4.*boxSolid->GetXHalfLength()*boxSolid->GetYHalfLength();
127  //
128  flux = flux/angleFactor; // Flux with angle.
129  if ( divideByArea ) flux /= square;
130  //
131  G4int index = GetIndex(aStep);
132  EvtMap->add(index,flux);
133  }
134  }
135 #ifdef debug
136  G4cout << " PASSED vol "
137  << index << " trk "<<trkid<<" len " << fFlatSurfaceFlux<<G4endl;
138 #endif
139 
140  return TRUE;
141 }
142 
144 
145  G4TouchableHandle theTouchable =
148 
149  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
150  // Entering Geometry
151  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
152  G4ThreeVector localpos1 =
153  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
154  if(std::fabs( localpos1.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
155  return fFlux_In;
156  }
157  }
158 
159  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
160  // Exiting Geometry
161  G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
162  G4ThreeVector localpos2 =
163  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
164  if(std::fabs( localpos2.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
165  return fFlux_Out;
166  }
167  }
168 
169  return -1;
170 }
171 
173 {
175  GetName());
176  if ( HCID < 0 ) HCID = GetCollectionID(0);
177  HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
178 }
179 
181 {;}
182 
184  EvtMap->clear();
185 }
186 
188 {;}
189 
191 {
192  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
193  G4cout << " PrimitiveScorer" << GetName() <<G4endl;
194  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
195  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
196  for(; itr != EvtMap->GetMap()->end(); itr++) {
197  G4cout << " copy no.: " << itr->first
198  << " flux : " << *(itr->second)/GetUnitValue()
199  << " [" << GetUnit() <<"]"
200  << G4endl;
201  }
202 }
203 
205 {
206  if ( divideByArea ) {
207  CheckAndSetUnit(unit,"Per Unit Surface");
208  } else {
209  if (unit == "" ){
210  unitName = unit;
211  unitValue = 1.0;
212  }else{
213  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
214  G4Exception("G4PSFlatSurfaceFlux::SetUnit","DetPS0008",JustWarning,msg);
215  }
216  }
217 }
218 
220  // Per Unit Surface
221  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
222  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
223  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
224 }
225 
226 
const XML_Char * name
Definition: expat.h:151
G4double GetXHalfLength() const
void clear()
Definition: G4THitsMap.hh:267
G4double GetWeight() const
G4String GetName() const
virtual void Initialize(G4HCofThisEvent *)
G4int IsSelectedSurface(G4Step *, G4Box *)
static constexpr double cm2
Definition: G4SIunits.hh:120
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:138
Definition: G4Box.hh:64
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4StepStatus GetStepStatus() const
virtual void DefineUnitAndCategory()
G4double GetSurfaceTolerance() const
G4MultiFunctionalDetector * GetMultiFunctionalDetector() const
G4VSolid * GetSolid() const
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:86
const G4VTouchable * GetTouchable() const
void CheckAndSetUnit(const G4String &unit, const G4String &category)
int G4int
Definition: G4Types.hh:78
G4double GetZHalfLength() const
virtual void SetUnit(const G4String &unit)
double z() const
static constexpr double mm2
Definition: G4SIunits.hh:116
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
const G4ThreeVector & GetPosition() const
G4double GetUnitValue() const
bool G4bool
Definition: G4Types.hh:79
virtual void EndOfEvent(G4HCofThisEvent *)
#define FALSE
Definition: globals.hh:52
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4int GetIndex(G4Step *)
G4double GetYHalfLength() const
#define TRUE
Definition: globals.hh:55
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
const G4double kCarTolerance
Definition: G4Step.hh:76
G4int GetCollectionID(G4int)
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int entries() const
Definition: G4THitsMap.hh:200
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:106
G4LogicalVolume * GetLogicalVolume() const
G4StepPoint * GetPostStepPoint() const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:99
#define G4endl
Definition: G4ios.hh:61
G4MultiFunctionalDetector * detector
const G4AffineTransform & GetTopTransform() const
double G4double
Definition: G4Types.hh:76
const G4String & GetUnit() const
G4PSFlatSurfaceFlux(G4String name, G4int direction, G4int depth=0)
static constexpr double m2
Definition: G4SIunits.hh:130
const G4TouchableHandle & GetTouchableHandle() const
static G4GeometryTolerance * GetInstance()