Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PSCylinderSurfaceFlux.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: G4PSCylinderSurfaceFlux.cc 81087 2014-05-20 15:44:27Z gcosmo $
28 //
29 // // G4PSCylinderSurfaceFlux
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"
40 // ////////////////////////////////////////////////////////////////////////////////
41 // (Description)
42 // This is a primitive scorer class for scoring Surface Flux.
43 // Current version assumes only for G4Tubs shape, and the surface
44 // is fixed on inner plane of the tube.
45 //
46 // Surface is defined at the innner surface of the tube.
47 // Direction R R+dR
48 // 0 IN || OUT ->|<- |
49 // 1 IN ->| |
50 // 2 OUT |<- |
51 //
52 // Created: 2007-03-29 Tsukasa ASO
53 // 2010-07-22 Introduce Unit specification.
54 // 2010-07-22 Add weighted and divideByArea options
55 // 2011-02-21 Get correct momentum direction in Flux_Out.
57 
59  G4int direction, G4int depth)
60  : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
61  weighted(true),divideByArea(true)
62 {
64  SetUnit("percm2");
65 }
66 
68  G4int direction,
69  const G4String& unit,
70  G4int depth)
71  : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
72  weighted(true),divideByArea(true)
73 {
75  SetUnit(unit);
76 }
77 
79 {;}
80 
82 {
83  G4StepPoint* preStep = aStep->GetPreStepPoint();
84 
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  G4Tubs* tubsSolid = (G4Tubs*)(solid);
101 
102  G4int dirFlag =IsSelectedSurface(aStep,tubsSolid);
103 
104  if ( dirFlag > 0 ){
105  if (fDirection == fFlux_InOut || dirFlag == fDirection ){
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  G4ThreeVector position = thisStep->GetPosition();
121  G4ThreeVector localpos =
122  theTouchable->GetHistory()->GetTopTransform().TransformAxis(position);
123  G4double angleFactor = (localdir.x()*localpos.x()+localdir.y()*localpos.y())
124  /std::sqrt(localdir.x()*localdir.x()
125  +localdir.y()*localdir.y()+localdir.z()*localdir.z())
126  /std::sqrt(localpos.x()*localpos.x()+localpos.y()*localpos.y());
127 
128  if ( angleFactor < 0 ) angleFactor *= -1.;
129  G4double square = 2.*tubsSolid->GetZHalfLength()
130  *tubsSolid->GetInnerRadius()* tubsSolid->GetDeltaPhiAngle()/radian;
131 
132  G4double flux = 1.0;
133  if ( weighted ) flux *=preStep->GetWeight();
134  // Current (Particle Weight)
135 
136  flux = flux/angleFactor;
137  if ( divideByArea ) flux /= square;
138  //Flux with angle.
139  G4int index = GetIndex(aStep);
140  EvtMap->add(index,flux);
141  return TRUE;
142  }else{
143  return FALSE;
144  }
145  }else{
146  return FALSE;
147  }
148 }
149 
151 
152  G4TouchableHandle theTouchable =
155 
156  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
157  // Entering Geometry
158  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
159  G4ThreeVector localpos1 =
160  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
161  if ( std::fabs(localpos1.z()) > tubsSolid->GetZHalfLength() ) return -1;
162  //if(std::fabs( localpos1.x()*localpos1.x()+localpos1.y()*localpos1.y()
163  // - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius()))
164  // <kCarTolerance ){
165  G4double localR2 = localpos1.x()*localpos1.x()+localpos1.y()*localpos1.y();
166  G4double InsideRadius = tubsSolid->GetInnerRadius();
167  if (localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
168  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
169  return fFlux_In;
170  }
171  }
172 
173  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
174  // Exiting Geometry
175  G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
176  G4ThreeVector localpos2 =
177  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
178  if ( std::fabs(localpos2.z()) > tubsSolid->GetZHalfLength() ) return -1;
179  //if(std::fabs( localpos2.x()*localpos2.x()+localpos2.y()*localpos2.y()
180  // - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius()))
181  // <kCarTolerance ){
182  G4double localR2 = localpos2.x()*localpos2.x()+localpos2.y()*localpos2.y();
183  G4double InsideRadius = tubsSolid->GetInnerRadius();
184  if (localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
185  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
186  return fFlux_Out;
187  }
188  }
189 
190  return -1;
191 }
192 
194 {
196  GetName());
197  if ( HCID < 0 ) HCID = GetCollectionID(0);
198  HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
199 }
200 
202 {;}
203 
205  EvtMap->clear();
206 }
207 
209 {;}
210 
212 {
213  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
214  G4cout << " PrimitiveScorer" << GetName() <<G4endl;
215  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
216  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
217  for(; itr != EvtMap->GetMap()->end(); itr++) {
218  G4cout << " copy no.: " << itr->first
219  << " flux : " << *(itr->second)/GetUnitValue()
220  << " ["<<GetUnit()<<"]"
221  << G4endl;
222  }
223 }
224 
226 {
227  if ( divideByArea ) {
228  CheckAndSetUnit(unit,"Per Unit Surface");
229  } else {
230  if (unit == "" ){
231  unitName = unit;
232  unitValue = 1.0;
233  }else{
234  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
235  G4Exception("G4PSCylinderSurfaceFlux::SetUnit","DetPS0003",JustWarning,msg);
236  }
237  }
238 }
239 
241  // Per Unit Surface
242  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
243  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
244  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
245 }
246 
247 
const XML_Char * name
Definition: expat.h:151
void clear()
Definition: G4THitsMap.hh:267
G4double GetWeight() const
G4String GetName() const
virtual void EndOfEvent(G4HCofThisEvent *)
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
Definition: G4Tubs.hh:85
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
double z() const
static constexpr double mm2
Definition: G4SIunits.hh:116
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
G4GLOB_DLL std::ostream G4cout
G4double GetDeltaPhiAngle() const
G4VPhysicalVolume * GetPhysicalVolume() const
const G4ThreeVector & GetPosition() const
G4double GetUnitValue() const
bool G4bool
Definition: G4Types.hh:79
#define FALSE
Definition: globals.hh:52
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4int GetIndex(G4Step *)
#define TRUE
Definition: globals.hh:55
const G4double kCarTolerance
Definition: G4Step.hh:76
virtual void Initialize(G4HCofThisEvent *)
G4int GetCollectionID(G4int)
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
virtual void SetUnit(const G4String &unit)
G4StepPoint * GetPostStepPoint() const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
double y() const
G4double GetZHalfLength() const
G4PSCylinderSurfaceFlux(G4String name, G4int direction, G4int depth=0)
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:99
#define G4endl
Definition: G4ios.hh:61
G4MultiFunctionalDetector * detector
const G4AffineTransform & GetTopTransform() const
G4int IsSelectedSurface(G4Step *, G4Tubs *)
double G4double
Definition: G4Types.hh:76
const G4String & GetUnit() const
static constexpr double m2
Definition: G4SIunits.hh:130
const G4TouchableHandle & GetTouchableHandle() const
static G4GeometryTolerance * GetInstance()
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)