Geant4  9.6.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$
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().
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 
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  G4Sphere* sphereSolid = (G4Sphere*)(solid);
102 
103  G4int dirFlag =IsSelectedSurface(aStep,sphereSolid);
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  G4double localdirL2 = localdir.x()*localdir.x()
121  +localdir.y()*localdir.y()
122  +localdir.z()*localdir.z();
123  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
124  G4ThreeVector localpos1 =
125  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
126  G4double localR2 = localpos1.x()*localpos1.x()
127  +localpos1.y()*localpos1.y()
128  +localpos1.z()*localpos1.z();
129  G4double anglefactor = (localdir.x()*localpos1.x()
130  +localdir.y()*localpos1.y()
131  +localdir.z()*localpos1.z())
132  /std::sqrt(localdirL2)/std::sqrt(localR2);
133 
134  G4double radi = sphereSolid->GetInsideRadius();
135  G4double dph = sphereSolid->GetDeltaPhiAngle()/radian;
136  G4double stth = sphereSolid->GetStartThetaAngle()/radian;
137  G4double enth = stth+sphereSolid->GetDeltaThetaAngle()/radian;
138  G4double square = radi*radi*dph*( -std::cos(enth) + std::cos(stth) );
139 
140  G4double current = 1.0;
141  if ( weighted ) thisStep->GetWeight(); // Flux (Particle Weight)
142  if ( divideByArea ) current = current/square; // Flux with angle.
143 
144  current /= anglefactor;
145 
146  G4int index = GetIndex(aStep);
147  EvtMap->add(index,current);
148  }
149  }
150 
151  return TRUE;
152 }
153 
155 
156  G4TouchableHandle theTouchable =
159 
160  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
161  // Entering Geometry
162  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
163  G4ThreeVector localpos1 =
164  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
165  G4double localR2 = localpos1.x()*localpos1.x()
166  +localpos1.y()*localpos1.y()
167  +localpos1.z()*localpos1.z();
168  //G4double InsideRadius2 =
169  // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
170  //if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
171  G4double InsideRadius = sphereSolid->GetInsideRadius();
172  if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
173  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
174  return fFlux_In;
175  }
176  }
177 
178  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
179  // Exiting Geometry
180  G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
181  G4ThreeVector localpos2 =
182  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
183  G4double localR2 = localpos2.x()*localpos2.x()
184  +localpos2.y()*localpos2.y()
185  +localpos2.z()*localpos2.z();
186  //G4double InsideRadius2 =
187  // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
188  //if(std::facb(localR2 - InsideRadius2) ) < kCarTolerance ){
189  G4double InsideRadius = sphereSolid->GetInsideRadius();
190  if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
191  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
192  return fFlux_Out;
193  }
194  }
195 
196  return -1;
197 }
198 
200 {
201  EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
202  if ( HCID < 0 ) HCID = GetCollectionID(0);
203  HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
204 }
205 
207 {;}
208 
210  EvtMap->clear();
211 }
212 
214 {;}
215 
217 {
218  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
219  G4cout << " PrimitiveScorer " << GetName() <<G4endl;
220  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
221  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
222  for(; itr != EvtMap->GetMap()->end(); itr++) {
223  G4cout << " copy no.: " << itr->first
224  << " Flux : " << *(itr->second)/GetUnitValue()
225  << " ["<<GetUnit()<<"]"
226  << G4endl;
227  }
228 }
229 
231 {
232  if ( divideByArea ) {
233  CheckAndSetUnit(unit,"Per Unit Surface");
234  } else {
235  if (unit == "" ){
236  unitName = unit;
237  unitValue = 1.0;
238  }else{
239  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
240  G4Exception("G4PSSphereSurfaceFlux::SetUnit","DetPS0016",JustWarning,msg);
241  }
242  }
243 }
244 
246  // Per Unit Surface
247  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
248  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
249  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
250 }
251