Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PSSphereSurfaceCurrent.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: G4PSSphereSurfaceCurrent.cc 81087 2014-05-20 15:44:27Z gcosmo $
28 //
29 // G4PSSphereSurfaceCurrent
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 Current.
43 // Current 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 // 2010-07-22 Introduce Unit specification.
53 //
55 
57  G4int direction, G4int depth)
58  :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
59  EvtMap(0),weighted(true),divideByArea(true)
60 {
62  SetUnit("percm2");
63 }
64 
66  G4int direction,
67  const G4String& unit,
68  G4int depth)
69  :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
70  EvtMap(0),weighted(true),divideByArea(true)
71 {
73  SetUnit(unit);
74 }
75 
77 {;}
78 
80 {
81  G4StepPoint* preStep = aStep->GetPreStepPoint();
82  G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
83  G4VPVParameterisation* physParam = physVol->GetParameterisation();
84  G4VSolid * solid = 0;
85  if(physParam)
86  { // for parameterized volume
88  ->GetReplicaNumber(indexDepth);
89  solid = physParam->ComputeSolid(idx, physVol);
90  solid->ComputeDimensions(physParam,idx,physVol);
91  }
92  else
93  { // for ordinary volume
94  solid = physVol->GetLogicalVolume()->GetSolid();
95  }
96 
97  G4Sphere* sphereSolid = (G4Sphere*)(solid);
98 
99  G4int dirFlag =IsSelectedSurface(aStep,sphereSolid);
100  if ( dirFlag > 0 ) {
101  if ( fDirection == fCurrent_InOut || fDirection == dirFlag ){
102  G4double radi = sphereSolid->GetInnerRadius();
103  G4double dph = sphereSolid->GetDeltaPhiAngle()/radian;
104  G4double stth = sphereSolid->GetStartThetaAngle()/radian;
105  G4double enth = stth+sphereSolid->GetDeltaThetaAngle()/radian;
106  G4double current = 1.0;
107  if ( weighted) current = preStep->GetWeight(); // Current (Particle Weight)
108  if ( divideByArea ){
109  G4double square = radi*radi*dph*( -std::cos(enth) + std::cos(stth) );
110  current = current/square; // Current with angle.
111  }
112 
113  G4int index = GetIndex(aStep);
114  EvtMap->add(index,current);
115  }
116  }
117 
118  return TRUE;
119 }
120 
122 
123  G4TouchableHandle theTouchable =
126 
127  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
128  // Entering Geometry
129  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
130  G4ThreeVector localpos1 =
131  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
132  G4double localR2 = localpos1.x()*localpos1.x()
133  +localpos1.y()*localpos1.y()
134  +localpos1.z()*localpos1.z();
135  //G4double InsideRadius2 =
136  // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
137  //if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
138  G4double InsideRadius = sphereSolid->GetInnerRadius();
139  if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
140  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
141  return fCurrent_In;
142  }
143  }
144 
145  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
146  // Exiting Geometry
147  G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
148  G4ThreeVector localpos2 =
149  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
150  G4double localR2 = localpos2.x()*localpos2.x()
151  +localpos2.y()*localpos2.y()
152  +localpos2.z()*localpos2.z();
153  //G4double InsideRadius2 =
154  // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
155  //if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
156  G4double InsideRadius = sphereSolid->GetInnerRadius();
157  if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
158  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
159  return fCurrent_Out;
160  }
161  }
162 
163  return -1;
164 }
165 
167 {
168  EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
169  if ( HCID < 0 ) HCID = GetCollectionID(0);
170  HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
171 }
172 
174 {;}
175 
177  EvtMap->clear();
178 }
179 
181 {;}
182 
184 {
185  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
186  G4cout << " PrimitiveScorer " << GetName() <<G4endl;
187  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
188  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
189  for(; itr != EvtMap->GetMap()->end(); itr++) {
190  G4cout << " copy no.: " << itr->first << " current : " ;
191  if ( divideByArea ) {
192  G4cout << *(itr->second)/GetUnitValue()
193  << " [" <<GetUnit()<<"]";
194  }else {
195  G4cout << *(itr->second) << " [tracks]" ;
196  }
197  G4cout << G4endl;
198  }
199 }
200 
201 
203 {
204  if ( divideByArea ) {
205  CheckAndSetUnit(unit,"Per Unit Surface");
206  } else {
207  if (unit == "" ){
208  unitName = unit;
209  unitValue = 1.0;
210  }else{
211  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
212  G4Exception("G4PSSphereSurfaceCurrent::SetUnit","DetPS0015",JustWarning,msg);
213  }
214  }
215 }
216 
218  // Per Unit Surface
219  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
220  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
221  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
222 }
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
virtual void EndOfEvent(G4HCofThisEvent *)
G4PSSphereSurfaceCurrent(G4String name, G4int direction, G4int depth=0)
const G4VTouchable * GetTouchable() const
virtual void SetUnit(const G4String &unit)
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
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetStartThetaAngle() const
const G4ThreeVector & GetPosition() const
G4double GetUnitValue() const
bool G4bool
Definition: G4Types.hh:79
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4int GetIndex(G4Step *)
virtual void Initialize(G4HCofThisEvent *)
#define TRUE
Definition: globals.hh:55
const G4double kCarTolerance
Definition: G4Step.hh:76
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
G4StepPoint * GetPostStepPoint() const
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
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
G4int IsSelectedSurface(G4Step *, G4Sphere *)
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()