Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PSPassageCellFlux.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 // GEANT4 tag $Name: geant4-09-04 $
29 //
30 // G4PSPassageCellFlux
31 #include "G4PSPassageCellFlux.hh"
32 
33 #include "G4SystemOfUnits.hh"
34 #include "G4StepStatus.hh"
35 #include "G4Track.hh"
36 #include "G4VSolid.hh"
37 #include "G4VPhysicalVolume.hh"
38 #include "G4VPVParameterisation.hh"
39 #include "G4UnitsTable.hh"
41 // (Description)
42 // This is a primitive scorer class for scoring cell flux.
43 // The Cell Flux is defined by a track length divided by a geometry
44 // volume, where only tracks passing through the geometry are taken
45 // into account. e.g. the unit of Cell Flux is mm/mm3.
46 //
47 // If you want to score all tracks in the geometry volume,
48 // please use G4PSCellFlux.
49 //
50 // Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
51 // 2010-07-22 Introduce Unit specification.
52 // 2010-07-22 Add weighted option
53 //
55 
57  : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
58  weighted(true)
59 {
61  SetUnit("percm2");
62 }
63 
65  G4int depth)
66  : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
67  weighted(true)
68 {
70  SetUnit(unit);
71 }
72 
74 {;}
75 
77 {
78 
79  if ( IsPassed(aStep) ) {
80  G4int idx = ((G4TouchableHistory*)
81  (aStep->GetPreStepPoint()->GetTouchable()))
82  ->GetReplicaNumber(indexDepth);
83  G4double cubicVolume = ComputeVolume(aStep, idx);
84 
85  fCellFlux /= cubicVolume;
86  G4int index = GetIndex(aStep);
87  EvtMap->add(index,fCellFlux);
88  }
89 
90  return TRUE;
91 }
92 
94  G4bool Passed = FALSE;
95 
96  G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
97  G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
98 
99  G4int trkid = aStep->GetTrack()->GetTrackID();
100  G4double trklength = aStep->GetStepLength();
101  if ( weighted ) trklength *= aStep->GetPreStepPoint()->GetWeight();
102 
103  if ( IsEnter &&IsExit ){ // Passed at one step
104  fCellFlux = trklength; // Track length is absolutely given.
105  Passed = TRUE;
106  }else if ( IsEnter ){ // Enter a new geometry
107  fCurrentTrkID = trkid; // Resetting the current track.
108  fCellFlux = trklength;
109  }else if ( IsExit ){ // Exit a current geometry
110  if ( fCurrentTrkID == trkid ) {// if the track is same as entered,
111  fCellFlux += trklength; // add the track length to current one.
112  Passed = TRUE;
113  }
114  }else{ // Inside geometry
115  if ( fCurrentTrkID == trkid ){ // if the track is same as entered,
116  fCellFlux += trklength; // adding the track length to current one.
117  }
118  }
119 
120  return Passed;
121 }
122 
124 {
125  fCurrentTrkID = -1;
126 
127  EvtMap = new G4THitsMap<G4double>(detector->GetName(),
128  GetName());
129  if ( HCID < 0 ) HCID = GetCollectionID(0);
130  HCE->AddHitsCollection(HCID,EvtMap);
131 
132 }
133 
135 {;}
136 
138  EvtMap->clear();
139 }
140 
142 {;}
143 
145 {
146  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
147  G4cout << " PrimitiveScorer " << GetName() <<G4endl;
148  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
149  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
150  for(; itr != EvtMap->GetMap()->end(); itr++) {
151  G4cout << " copy no.: " << itr->first
152  << " cell flux : " << *(itr->second)/GetUnitValue()
153  << " [" << GetUnit()
154  << G4endl;
155  }
156 }
157 
159 {
160  CheckAndSetUnit(unit,"Per Unit Surface");
161 }
162 
164  // Per Unit Surface
165  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
166  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
167  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
168 }
169 
170 
172 
173  G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
174  G4VPVParameterisation* physParam = physVol->GetParameterisation();
175  G4VSolid* solid = 0;
176  if(physParam)
177  { // for parameterized volume
178  if(idx<0)
179  {
181  ED << "Incorrect replica number --- GetReplicaNumber : " << idx << G4endl;
182  G4Exception("G4PSPassageCellFlux::ComputeVolume","DetPS0013",JustWarning,ED);
183  }
184  solid = physParam->ComputeSolid(idx, physVol);
185  solid->ComputeDimensions(physParam,idx,physVol);
186  }
187  else
188  { // for ordinary volume
189  solid = physVol->GetLogicalVolume()->GetSolid();
190  }
191 
192  return solid->GetCubicVolume();
193 }