source: trunk/source/digits_hits/scorer/src/G4PSPassageCellFlux.cc@ 1350

Last change on this file since 1350 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 6.4 KB
Line 
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: G4PSPassageCellFlux.cc,v 1.4 2010/07/22 23:42:01 taso Exp $
28// GEANT4 tag $Name: $
29//
30// G4PSPassageCellFlux
31#include "G4PSPassageCellFlux.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////////////////////////////////////////////////////////////////////////////////
39// (Description)
40// This is a primitive scorer class for scoring cell flux.
41// The Cell Flux is defined by a track length divided by a geometry
42// volume, where only tracks passing through the geometry are taken
43// into account. e.g. the unit of Cell Flux is mm/mm3.
44//
45// If you want to score all tracks in the geometry volume,
46// please use G4PSCellFlux.
47//
48// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
49// 2010-07-22 Introduce Unit specification.
50// 2010-07-22 Add weighted option
51//
52///////////////////////////////////////////////////////////////////////////////
53
54G4PSPassageCellFlux::G4PSPassageCellFlux(G4String name, G4int depth)
55 :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
56 weighted(true)
57{
58 DefineUnitAndCategory();
59 SetUnit("percm2");
60}
61
62G4PSPassageCellFlux::G4PSPassageCellFlux(G4String name, const G4String& unit,
63 G4int depth)
64 :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0)
65{
66 DefineUnitAndCategory();
67 SetUnit(unit);
68}
69
70G4PSPassageCellFlux::~G4PSPassageCellFlux()
71{;}
72
73G4bool G4PSPassageCellFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
74{
75
76 if ( IsPassed(aStep) ) {
77 G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
78 G4VPVParameterisation* physParam = physVol->GetParameterisation();
79 G4VSolid* solid = 0;
80 if(physParam)
81 { // for parameterized volume
82 G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
83 ->GetReplicaNumber(indexDepth);
84 solid = physParam->ComputeSolid(idx, physVol);
85 solid->ComputeDimensions(physParam,idx,physVol);
86 }
87 else
88 { // for ordinary volume
89 solid = physVol->GetLogicalVolume()->GetSolid();
90 }
91
92 fCellFlux /= solid->GetCubicVolume();
93 G4int index = GetIndex(aStep);
94 EvtMap->add(index,fCellFlux);
95 }
96
97 return TRUE;
98}
99
100G4bool G4PSPassageCellFlux::IsPassed(G4Step* aStep){
101 G4bool Passed = FALSE;
102
103 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
104 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
105
106 G4int trkid = aStep->GetTrack()->GetTrackID();
107 G4double trklength = aStep->GetStepLength();
108 if ( weighted ) trklength *= aStep->GetPreStepPoint()->GetWeight();
109
110 if ( IsEnter &&IsExit ){ // Passed at one step
111 fCellFlux = trklength; // Track length is absolutely given.
112 Passed = TRUE;
113 }else if ( IsEnter ){ // Enter a new geometry
114 fCurrentTrkID = trkid; // Resetting the current track.
115 fCellFlux = trklength;
116 }else if ( IsExit ){ // Exit a current geometry
117 if ( fCurrentTrkID == trkid ) {// if the track is same as entered,
118 fCellFlux += trklength; // add the track length to current one.
119 Passed = TRUE;
120 }
121 }else{ // Inside geometry
122 if ( fCurrentTrkID == trkid ){ // if the track is same as entered,
123 fCellFlux += trklength; // adding the track length to current one.
124 }
125 }
126
127 return Passed;
128}
129
130void G4PSPassageCellFlux::Initialize(G4HCofThisEvent* HCE)
131{
132 fCurrentTrkID = -1;
133
134 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
135 GetName());
136 if ( HCID < 0 ) HCID = GetCollectionID(0);
137 HCE->AddHitsCollection(HCID,EvtMap);
138
139}
140
141void G4PSPassageCellFlux::EndOfEvent(G4HCofThisEvent*)
142{;}
143
144void G4PSPassageCellFlux::clear(){
145 EvtMap->clear();
146}
147
148void G4PSPassageCellFlux::DrawAll()
149{;}
150
151void G4PSPassageCellFlux::PrintAll()
152{
153 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
154 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
155 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
156 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
157 for(; itr != EvtMap->GetMap()->end(); itr++) {
158 G4cout << " copy no.: " << itr->first
159 << " cell flux : " << *(itr->second)/GetUnitValue()
160 << " [" << GetUnit()
161 << G4endl;
162 }
163}
164
165void G4PSPassageCellFlux::SetUnit(const G4String& unit)
166{
167 CheckAndSetUnit(unit,"Per Unit Surface");
168}
169
170void G4PSPassageCellFlux::DefineUnitAndCategory(){
171 // Per Unit Surface
172 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
173 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
174 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
175}
176
177
Note: See TracBrowser for help on using the repository browser.