Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PSSphereSurfaceFlux.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: G4PSSphereSurfaceFlux.cc 81087 2014-05-20 15:44:27Z gcosmo $
28 //
29 // G4PSSphereSurfaceFlux
30 #include "G4PSSphereSurfaceFlux.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 only Surface Flux.
43 // Flux version assumes only for G4Sphere shape.
44 //
45 // Surface is defined at the inside of sphere.
46 // Direction -Rmin +Rmax
47 // 0 IN || OUT ->|<- |
48 // 1 IN ->| |
49 // 2 OUT |<- |
50 //
51 // Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
52 // 29-Mar-2007 T.Aso, Bug fix for momentum direction at outgoing flux.
53 // 2010-07-22 Introduce Unit specification.
54 // 2010-07-22 Add weighted and divideByAre options
55 // 2011-02-21 Get correct momentum direction in Flux_Out.
56 // 2011-09-09 Modify comment in PrintAll().
57 // 2014-03-03 T.Aso, To use always positive value for anglefactor.
59 
61  G4int direction, G4int depth)
62  : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
63  EvtMap(0),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),
74  EvtMap(0),weighted(true),divideByArea(true)
75 {
77  SetUnit(unit);
78 }
79 
81 {;}
82 
84 {
85  G4StepPoint* preStep = aStep->GetPreStepPoint();
86 
87  G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
88  G4VPVParameterisation* physParam = physVol->GetParameterisation();
89  G4VSolid * solid = 0;
90  if(physParam)
91  { // for parameterized volume
93  ->GetReplicaNumber(indexDepth);
94  solid = physParam->ComputeSolid(idx, physVol);
95  solid->ComputeDimensions(physParam,idx,physVol);
96  }
97  else
98  { // for ordinary volume
99  solid = physVol->GetLogicalVolume()->GetSolid();
100  }
101 
102  G4Sphere* sphereSolid = (G4Sphere*)(solid);
103 
104  G4int dirFlag =IsSelectedSurface(aStep,sphereSolid);
105  if ( dirFlag > 0 ) {
106  if ( fDirection == fFlux_InOut || fDirection == dirFlag ){
107 
108  G4StepPoint* thisStep=0;
109  if ( dirFlag == fFlux_In ){
110  thisStep = preStep;
111  }else if ( dirFlag == fFlux_Out ){
112  thisStep = aStep->GetPostStepPoint();
113  }else{
114  return FALSE;
115  }
116 
117  G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
118  G4ThreeVector pdirection = thisStep->GetMomentumDirection();
119  G4ThreeVector localdir =
120  theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
121  G4double localdirL2 = localdir.x()*localdir.x()
122  +localdir.y()*localdir.y()
123  +localdir.z()*localdir.z();
124  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
125  G4ThreeVector localpos1 =
126  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
127  G4double localR2 = localpos1.x()*localpos1.x()
128  +localpos1.y()*localpos1.y()
129  +localpos1.z()*localpos1.z();
130  G4double anglefactor = (localdir.x()*localpos1.x()
131  +localdir.y()*localpos1.y()
132  +localdir.z()*localpos1.z())
133  /std::sqrt(localdirL2)/std::sqrt(localR2);
134  if ( anglefactor < 0.0 ) anglefactor *= -1.0;
135 
136  G4double radi = sphereSolid->GetInnerRadius();
137  G4double dph = sphereSolid->GetDeltaPhiAngle()/radian;
138  G4double stth = sphereSolid->GetStartThetaAngle()/radian;
139  G4double enth = stth+sphereSolid->GetDeltaThetaAngle()/radian;
140  G4double square = radi*radi*dph*( -std::cos(enth) + std::cos(stth) );
141 
142  G4double current = 1.0;
143  if ( weighted ) thisStep->GetWeight(); // Flux (Particle Weight)
144  if ( divideByArea ) current = current/square; // Flux with angle.
145 
146  current /= anglefactor;
147 
148  G4int index = GetIndex(aStep);
149  EvtMap->add(index,current);
150  }
151  }
152 
153  return TRUE;
154 }
155 
157 
158  G4TouchableHandle theTouchable =
161 
162  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
163  // Entering Geometry
164  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
165  G4ThreeVector localpos1 =
166  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
167  G4double localR2 = localpos1.x()*localpos1.x()
168  +localpos1.y()*localpos1.y()
169  +localpos1.z()*localpos1.z();
170  //G4double InsideRadius2 =
171  // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
172  //if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
173  G4double InsideRadius = sphereSolid->GetInnerRadius();
174  if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
175  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
176  return fFlux_In;
177  }
178  }
179 
180  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
181  // Exiting Geometry
182  G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
183  G4ThreeVector localpos2 =
184  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
185  G4double localR2 = localpos2.x()*localpos2.x()
186  +localpos2.y()*localpos2.y()
187  +localpos2.z()*localpos2.z();
188  //G4double InsideRadius2 =
189  // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
190  //if(std::facb(localR2 - InsideRadius2) ) < kCarTolerance ){
191  G4double InsideRadius = sphereSolid->GetInnerRadius();
192  if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
193  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
194  return fFlux_Out;
195  }
196  }
197 
198  return -1;
199 }
200 
202 {
203  EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
204  if ( HCID < 0 ) HCID = GetCollectionID(0);
205  HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
206 }
207 
209 {;}
210 
212  EvtMap->clear();
213 }
214 
216 {;}
217 
219 {
220  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
221  G4cout << " PrimitiveScorer " << GetName() <<G4endl;
222  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
223  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
224  for(; itr != EvtMap->GetMap()->end(); itr++) {
225  G4cout << " copy no.: " << itr->first
226  << " Flux : " << *(itr->second)/GetUnitValue()
227  << " ["<<GetUnit()<<"]"
228  << G4endl;
229  }
230 }
231 
233 {
234  if ( divideByArea ) {
235  CheckAndSetUnit(unit,"Per Unit Surface");
236  } else {
237  if (unit == "" ){
238  unitName = unit;
239  unitValue = 1.0;
240  }else{
241  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
242  G4Exception("G4PSSphereSurfaceFlux::SetUnit","DetPS0016",JustWarning,msg);
243  }
244  }
245 }
246 
248  // Per Unit Surface
249  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
250  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
251  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
252 }
253 
const XML_Char * name
Definition: expat.h:151
void clear()
Definition: G4THitsMap.hh:267
G4double GetWeight() const
G4String GetName() const
static constexpr double cm2
Definition: G4SIunits.hh:120
double x() const
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:138
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4StepStatus GetStepStatus() const
G4double GetSurfaceTolerance() 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)
G4double GetDeltaPhiAngle() const
int G4int
Definition: G4Types.hh:78
double z() const
static constexpr double mm2
Definition: G4SIunits.hh:116
virtual void Initialize(G4HCofThisEvent *)
G4StepPoint * GetPreStepPoint() const
virtual void DefineUnitAndCategory()
const G4ThreeVector & GetMomentumDirection() const
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetStartThetaAngle() const
const G4ThreeVector & GetPosition() const
G4double GetUnitValue() const
virtual void SetUnit(const G4String &unit)
bool G4bool
Definition: G4Types.hh:79
#define FALSE
Definition: globals.hh:52
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4int GetIndex(G4Step *)
G4int IsSelectedSurface(G4Step *, G4Sphere *)
#define TRUE
Definition: globals.hh:55
const G4double kCarTolerance
Definition: G4Step.hh:76
G4int GetCollectionID(G4int)
G4PSSphereSurfaceFlux(G4String name, G4int direction, G4int depth=0)
G4double GetInnerRadius() const
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
static constexpr double radian
Definition: G4SIunits.hh:142
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
virtual void EndOfEvent(G4HCofThisEvent *)
double y() 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
G4double GetDeltaThetaAngle() const
static constexpr double m2
Definition: G4SIunits.hh:130
const G4TouchableHandle & GetTouchableHandle() const
static G4GeometryTolerance * GetInstance()
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)