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

Last change on this file since 1317 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 5.6 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.2 2008/12/28 20:32:00 asaim Exp $
28// GEANT4 tag $Name: geant4-09-03 $
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//
50///////////////////////////////////////////////////////////////////////////////
51
52G4PSPassageCellFlux::G4PSPassageCellFlux(G4String name, G4int depth)
53 :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0)
54{;}
55
56G4PSPassageCellFlux::~G4PSPassageCellFlux()
57{;}
58
59G4bool G4PSPassageCellFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
60{
61
62 if ( IsPassed(aStep) ) {
63 G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
64 G4VPVParameterisation* physParam = physVol->GetParameterisation();
65 G4VSolid* solid = 0;
66 if(physParam)
67 { // for parameterized volume
68 G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
69 ->GetReplicaNumber(indexDepth);
70 solid = physParam->ComputeSolid(idx, physVol);
71 solid->ComputeDimensions(physParam,idx,physVol);
72 }
73 else
74 { // for ordinary volume
75 solid = physVol->GetLogicalVolume()->GetSolid();
76 }
77
78 fCellFlux /= solid->GetCubicVolume();
79 G4int index = GetIndex(aStep);
80 EvtMap->add(index,fCellFlux);
81 }
82
83 return TRUE;
84}
85
86G4bool G4PSPassageCellFlux::IsPassed(G4Step* aStep){
87 G4bool Passed = FALSE;
88
89 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
90 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
91
92 G4int trkid = aStep->GetTrack()->GetTrackID();
93 G4double trklength = aStep->GetStepLength();
94 trklength *= aStep->GetPreStepPoint()->GetWeight();
95
96 if ( IsEnter &&IsExit ){ // Passed at one step
97 fCellFlux = trklength; // Track length is absolutely given.
98 Passed = TRUE;
99 }else if ( IsEnter ){ // Enter a new geometry
100 fCurrentTrkID = trkid; // Resetting the current track.
101 fCellFlux = trklength;
102 }else if ( IsExit ){ // Exit a current geometry
103 if ( fCurrentTrkID == trkid ) {// if the track is same as entered,
104 fCellFlux += trklength; // add the track length to current one.
105 Passed = TRUE;
106 }
107 }else{ // Inside geometry
108 if ( fCurrentTrkID == trkid ){ // if the track is same as entered,
109 fCellFlux += trklength; // adding the track length to current one.
110 }
111 }
112
113 return Passed;
114}
115
116void G4PSPassageCellFlux::Initialize(G4HCofThisEvent* HCE)
117{
118 fCurrentTrkID = -1;
119
120 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
121 GetName());
122 if ( HCID < 0 ) HCID = GetCollectionID(0);
123 HCE->AddHitsCollection(HCID,EvtMap);
124
125}
126
127void G4PSPassageCellFlux::EndOfEvent(G4HCofThisEvent*)
128{;}
129
130void G4PSPassageCellFlux::clear(){
131 EvtMap->clear();
132}
133
134void G4PSPassageCellFlux::DrawAll()
135{;}
136
137void G4PSPassageCellFlux::PrintAll()
138{
139 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
140 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
141 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
142 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
143 for(; itr != EvtMap->GetMap()->end(); itr++) {
144 G4cout << " copy no.: " << itr->first
145 << " cell flux : " << *(itr->second)*cm*cm << " [cm^-2]"
146 << G4endl;
147 }
148}
149
Note: See TracBrowser for help on using the repository browser.