Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PSCylinderSurfaceCurrent.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 // $Id: G4PSCylinderSurfaceCurrent.cc 81087 2014-05-20 15:44:27Z gcosmo $
28 //
29 // G4PSCylinderSurfaceCurrent
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 G4Tubs shape.
44 //
45 // Surface is defined at the inner surface of the tube.
46 // Direction R R+dR
47 // 0 IN || OUT ->|<- |
48 // 1 IN ->| |
49 // 2 OUT |<- |
50 //
51 // Created: 2007-03-21 Tsukasa ASO
52 // 2010-07-22 Introduce Unit specification.
53 //
55 
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  G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
84  G4VPVParameterisation* physParam = physVol->GetParameterisation();
85  G4VSolid * solid = 0;
86  if(physParam)
87  { // for parameterized volume
89  ->GetReplicaNumber(indexDepth);
90  solid = physParam->ComputeSolid(idx, physVol);
91  solid->ComputeDimensions(physParam,idx,physVol);
92  }
93  else
94  { // for ordinary volume
95  solid = physVol->GetLogicalVolume()->GetSolid();
96  }
97 
98  G4Tubs* tubsSolid = (G4Tubs*)(solid);
99 
100  G4int dirFlag =IsSelectedSurface(aStep,tubsSolid);
101  // G4cout << " pos " << preStep->GetPosition() <<" dirFlag " << G4endl;
102  if ( dirFlag > 0 ) {
103  if ( fDirection == fCurrent_InOut || fDirection == dirFlag ){
104  G4TouchableHandle theTouchable = preStep->GetTouchableHandle();
105  //
106  G4double current = 1.0;
107  if ( weighted ) current = preStep->GetWeight(); // Current (Particle Weight)
108  //
109  if ( divideByArea ){
110  G4double square = 2.*tubsSolid->GetZHalfLength()
111  *tubsSolid->GetInnerRadius()* tubsSolid->GetDeltaPhiAngle()/radian;
112  current = current/square; // Current normalized by Area
113  }
114 
115  G4int index = GetIndex(aStep);
116  EvtMap->add(index,current);
117  }
118 
119  }
120 
121  return TRUE;
122 }
123 
125 
126  G4TouchableHandle theTouchable =
129 
130  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
131  // Entering Geometry
132  G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
133  G4ThreeVector localpos1 =
134  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
135  if ( std::fabs(localpos1.z()) > tubsSolid->GetZHalfLength() ) return -1;
136  G4double localR2 = localpos1.x()*localpos1.x()+localpos1.y()*localpos1.y();
137  G4double InsideRadius = tubsSolid->GetInnerRadius();
138  if (localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
139  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
140  return fCurrent_In;
141  }
142  }
143 
144  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
145  // Exiting Geometry
146  G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
147  G4ThreeVector localpos2 =
148  theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
149  if ( std::fabs(localpos2.z()) > tubsSolid->GetZHalfLength() ) return -1;
150  G4double localR2 = localpos2.x()*localpos2.x()+localpos2.y()*localpos2.y();
151  G4double InsideRadius = tubsSolid->GetInnerRadius();
152  if (localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
153  &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
154  return fCurrent_Out;
155  }
156  }
157 
158  return -1;
159 }
160 
162 {
164  if ( HCID < 0 ) HCID = GetCollectionID(0);
166 }
167 
169 {;}
170 
172  EvtMap->clear();
173 }
174 
176 {;}
177 
179 {
180  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
181  G4cout << " PrimitiveScorer " << GetName() <<G4endl;
182  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
183  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
184  for(; itr != EvtMap->GetMap()->end(); itr++) {
185  G4cout << " copy no.: " << itr->first
186  << " current : " ;
187  if ( divideByArea ) {
188  G4cout << *(itr->second)/GetUnitValue()
189  << " ["<<GetUnit()<<"]";
190  } else {
191  G4cout << *(itr->second) << " [tracks]";
192  }
193  G4cout << G4endl;
194  }
195 }
196 
198 {
199  if ( divideByArea ) {
200  CheckAndSetUnit(unit,"Per Unit Surface");
201  } else {
202  if (unit == "" ){
203  unitName = unit;
204  unitValue = 1.0;
205  }else{
206  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
207  G4Exception("G4PSCylinderSurfaceCurrent::SetUnit","DetPS0002",
208  JustWarning,msg);
209  }
210  }
211 }
212 
214  // Per Unit Surface
215  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
216  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
217  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
218 }
219 
const G4String & GetUnit() const
const XML_Char * name
Definition: expat.h:151
const G4VTouchable * GetTouchable() const
G4int GetCollectionID(G4int)
G4LogicalVolume * GetLogicalVolume() const
G4String GetName() const
G4StepPoint * GetPreStepPoint() const
static constexpr double m2
Definition: G4SIunits.hh:130
virtual void Initialize(G4HCofThisEvent *)
Definition: G4Tubs.hh:85
G4int add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:150
virtual G4VPVParameterisation * GetParameterisation() const =0
#define G4endl
Definition: G4ios.hh:61
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
double z() const
G4VPhysicalVolume * GetPhysicalVolume() const
const G4double kCarTolerance
G4double GetDeltaPhiAngle() const
virtual G4int GetIndex(G4Step *)
const G4TouchableHandle & GetTouchableHandle() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4StepStatus GetStepStatus() const
static constexpr double cm2
Definition: G4SIunits.hh:120
virtual void SetUnit(const G4String &unit)
G4double GetInnerRadius() const
const G4ThreeVector & GetPosition() const
static constexpr double mm2
Definition: G4SIunits.hh:116
G4int IsSelectedSurface(G4Step *, G4Tubs *)
const G4AffineTransform & GetTopTransform() const
void clear()
Definition: G4THitsMap.hh:510
#define TRUE
Definition: globals.hh:55
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
G4double GetZHalfLength() const
virtual void EndOfEvent(G4HCofThisEvent *)
void CheckAndSetUnit(const G4String &unit, const G4String &category)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
G4PSCylinderSurfaceCurrent(G4String name, G4int direction, G4int depth=0)
static constexpr double radian
Definition: G4SIunits.hh:142
G4GLOB_DLL std::ostream G4cout
double x() const
G4VSolid * GetSolid() const
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
G4double GetWeight() const
static G4GeometryTolerance * GetInstance()
G4MultiFunctionalDetector * detector
G4double GetUnitValue() const
double y() const
G4int entries() const
Definition: G4THitsMap.hh:121
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:86
G4double GetSurfaceTolerance() const
Map_t * GetMap() const
Definition: G4THitsMap.hh:117