ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PSCylinderSurfaceFlux.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4PSCylinderSurfaceFlux.cc
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 //
28 // // G4PSCylinderSurfaceFlux
30 
31 #include "G4SystemOfUnits.hh"
32 #include "G4StepStatus.hh"
33 #include "G4Track.hh"
34 #include "G4VSolid.hh"
35 #include "G4VPhysicalVolume.hh"
36 #include "G4VPVParameterisation.hh"
37 #include "G4UnitsTable.hh"
38 #include "G4GeometryTolerance.hh"
39 // ////////////////////////////////////////////////////////////////////////////////
40 // (Description)
41 // This is a primitive scorer class for scoring Surface Flux.
42 // Current version assumes only for G4Tubs shape, and the surface
43 // is fixed on inner plane of the tube.
44 //
45 // Surface is defined at the innner surface of the tube.
46 // Direction R R+dR
47 // 0 IN || OUT ->|<- |
48 // 1 IN ->| |
49 // 2 OUT |<- |
50 //
51 // Created: 2007-03-29 Tsukasa ASO
52 // 2010-07-22 Introduce Unit specification.
53 // 2010-07-22 Add weighted and divideByArea options
54 // 2011-02-21 Get correct momentum direction in Flux_Out.
56 
58  G4int direction, G4int depth)
59  : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
60  weighted(true),divideByArea(true)
61 {
63  SetUnit("percm2");
64 }
65 
67  G4int direction,
68  const G4String& unit,
69  G4int depth)
70  : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),EvtMap(0),
71  weighted(true),divideByArea(true)
72 {
74  SetUnit(unit);
75 }
76 
78 {;}
79 
81 {
82  G4StepPoint* preStep = aStep->GetPreStepPoint();
83 
84  G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
85  G4VPVParameterisation* physParam = physVol->GetParameterisation();
86  G4VSolid * solid = 0;
87  if(physParam)
88  { // for parameterized volume
90  ->GetReplicaNumber(indexDepth);
91  solid = physParam->ComputeSolid(idx, physVol);
92  solid->ComputeDimensions(physParam,idx,physVol);
93  }
94  else
95  { // for ordinary volume
96  solid = physVol->GetLogicalVolume()->GetSolid();
97  }
98 
99  G4Tubs* tubsSolid = (G4Tubs*)(solid);
100 
101  G4int dirFlag =IsSelectedSurface(aStep,tubsSolid);
102 
103  if ( dirFlag > 0 ){
104  if (fDirection == fFlux_InOut || dirFlag == fDirection ){
105 
106  G4StepPoint* thisStep=0;
107  if ( dirFlag == fFlux_In ){
108  thisStep = preStep;
109  }else if ( dirFlag == fFlux_Out ){
110  thisStep = aStep->GetPostStepPoint();
111  }else{
112  return FALSE;
113  }
114 
115  G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
116  G4ThreeVector pdirection = thisStep->GetMomentumDirection();
117  G4ThreeVector localdir =
118  theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
119  G4ThreeVector position = thisStep->GetPosition();
120  G4ThreeVector localpos =
121  theTouchable->GetHistory()->GetTopTransform().TransformAxis(position);
122  G4double angleFactor = (localdir.x()*localpos.x()+localdir.y()*localpos.y())
123  /std::sqrt(localdir.x()*localdir.x()
124  +localdir.y()*localdir.y()+localdir.z()*localdir.z())
125  /std::sqrt(localpos.x()*localpos.x()+localpos.y()*localpos.y());
126 
127  if ( angleFactor < 0 ) angleFactor *= -1.;
128  G4double square = 2.*tubsSolid->GetZHalfLength()
129  *tubsSolid->GetInnerRadius()* tubsSolid->GetDeltaPhiAngle()/radian;
130 
131  G4double flux = 1.0;
132  if ( weighted ) flux *=preStep->GetWeight();
133  // Current (Particle Weight)
134 
135  flux = flux/angleFactor;
136  if ( divideByArea ) flux /= square;
137  //Flux with angle.
138  G4int index = GetIndex(aStep);
139  EvtMap->add(index,flux);
140  return TRUE;
141  }else{
142  return FALSE;
143  }
144  }else{
145  return FALSE;
146  }
147 }
148 
150 
151  G4TouchableHandle theTouchable =
154 
155  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
156  // Entering Geometry
157  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
158  G4ThreeVector localpos1 =
159  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
160  if ( std::fabs(localpos1.z()) > tubsSolid->GetZHalfLength() ) return -1;
161  //if(std::fabs( localpos1.x()*localpos1.x()+localpos1.y()*localpos1.y()
162  // - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius()))
163  // <kCarTolerance ){
164  G4double localR2 = localpos1.x()*localpos1.x()+localpos1.y()*localpos1.y();
165  G4double InsideRadius = tubsSolid->GetInnerRadius();
166  if (localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
167  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
168  return fFlux_In;
169  }
170  }
171 
172  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
173  // Exiting Geometry
174  G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
175  G4ThreeVector localpos2 =
176  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
177  if ( std::fabs(localpos2.z()) > tubsSolid->GetZHalfLength() ) return -1;
178  //if(std::fabs( localpos2.x()*localpos2.x()+localpos2.y()*localpos2.y()
179  // - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius()))
180  // <kCarTolerance ){
181  G4double localR2 = localpos2.x()*localpos2.x()+localpos2.y()*localpos2.y();
182  G4double InsideRadius = tubsSolid->GetInnerRadius();
183  if (localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
184  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
185  return fFlux_Out;
186  }
187  }
188 
189  return -1;
190 }
191 
193 {
195  GetName());
196  if ( HCID < 0 ) HCID = GetCollectionID(0);
198 }
199 
201 {;}
202 
204  EvtMap->clear();
205 }
206 
208 {;}
209 
211 {
212  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
213  G4cout << " PrimitiveScorer" << GetName() <<G4endl;
214  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
215  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
216  for(; itr != EvtMap->GetMap()->end(); itr++) {
217  G4cout << " copy no.: " << itr->first
218  << " flux : " << *(itr->second)/GetUnitValue()
219  << " ["<<GetUnit()<<"]"
220  << G4endl;
221  }
222 }
223 
225 {
226  if ( divideByArea ) {
227  CheckAndSetUnit(unit,"Per Unit Surface");
228  } else {
229  if (unit == "" ){
230  unitName = unit;
231  unitValue = 1.0;
232  }else{
233  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
234  G4Exception("G4PSCylinderSurfaceFlux::SetUnit","DetPS0003",JustWarning,msg);
235  }
236  }
237 }
238 
240  // Per Unit Surface
241  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
242  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
243  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
244 }
245 
246