Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PSCellFlux.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 // G4PSCellFlux
31 #include "G4PSCellFlux.hh"
32 
33 #include "G4SystemOfUnits.hh"
34 #include "G4Track.hh"
35 #include "G4VSolid.hh"
36 #include "G4VPhysicalVolume.hh"
37 #include "G4VPVParameterisation.hh"
38 #include "G4UnitsTable.hh"
40 // (Description)
41 // This is a primitive scorer class for scoring cell flux.
42 // The Cell Flux is defined by a sum of track length divided
43 // by the geometry volume, where all of the tracks in the geometry
44 // are taken into account.
45 //
46 // If you want to score only tracks passing through the geometry volume,
47 // please use G4PSPassageCellFlux.
48 //
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),weighted(true)
58 {
60  SetUnit("percm2");
61  //verboseLevel = 10;
62 }
63 
65  :G4VPrimitiveScorer(name,depth),HCID(-1),weighted(true)
66 {
68  SetUnit(unit);
69 }
70 
72 {;}
73 
75 {
76  G4double stepLength = aStep->GetStepLength();
77  if ( stepLength == 0. ) return FALSE;
78 
79  G4int idx = ((G4TouchableHistory*)
80  (aStep->GetPreStepPoint()->GetTouchable()))
81  ->GetReplicaNumber(indexDepth);
82  G4double cubicVolume = ComputeVolume(aStep, idx);
83 
84  G4double CellFlux = stepLength / cubicVolume;
85  if (weighted) CellFlux *= aStep->GetPreStepPoint()->GetWeight();
86  G4int index = GetIndex(aStep);
87  EvtMap->add(index,CellFlux);
88 
89  return TRUE;
90 }
91 
93 {
94  EvtMap = new G4THitsMap<G4double>(detector->GetName(),
95  GetName());
96  if ( HCID < 0 ) HCID = GetCollectionID(0);
97  HCE->AddHitsCollection(HCID,EvtMap);
98 }
99 
101 {;}
102 
104  EvtMap->clear();
105 }
106 
108 {;}
109 
111 {
112  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
113  G4cout << " PrimitiveScorer " << GetName() <<G4endl;
114  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
115  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
116  for(; itr != EvtMap->GetMap()->end(); itr++) {
117  G4cout << " copy no.: " << itr->first
118  << " cell flux : " << *(itr->second)/GetUnitValue()
119  << " [" << GetUnit() << "]"
120  << G4endl;
121  }
122 }
123 
125 {
126  CheckAndSetUnit(unit,"Per Unit Surface");
127 }
128 
130  // Per Unit Surface
131  new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
132  new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
133  new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
134 }
135 
137 
138  G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
139  G4VPVParameterisation* physParam = physVol->GetParameterisation();
140  G4VSolid* solid = 0;
141  if(physParam)
142  { // for parameterized volume
143  if(idx<0)
144  {
146  ED << "Incorrect replica number --- GetReplicaNumber : " << idx << G4endl;
147  G4Exception("G4PSCellFlux::ComputeVolume","DetPS0001",JustWarning,ED);
148  }
149  solid = physParam->ComputeSolid(idx, physVol);
150  solid->ComputeDimensions(physParam,idx,physVol);
151  }
152  else
153  { // for ordinary volume
154  solid = physVol->GetLogicalVolume()->GetSolid();
155  }
156 
157  return solid->GetCubicVolume();
158 }