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

Last change on this file since 936 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 5.0 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.1 2007/07/11 01:31:03 asaim Exp $
28// GEANT4 tag $Name: HEAD $
29//
30// G4PSPassageCellFlux
31#include "G4PSPassageCellFlux.hh"
32#include "G4StepStatus.hh"
33#include "G4Track.hh"
34#include "G4VSolid.hh"
35#include "G4UnitsTable.hh"
36////////////////////////////////////////////////////////////////////////////////
37// (Description)
38// This is a primitive scorer class for scoring cell flux.
39// The Cell Flux is defined by a track length divided by a geometry
40// volume, where only tracks passing through the geometry are taken
41// into account. e.g. the unit of Cell Flux is mm/mm3.
42//
43// If you want to score all tracks in the geometry volume,
44// please use G4PSCellFlux.
45//
46// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
47//
48///////////////////////////////////////////////////////////////////////////////
49
50G4PSPassageCellFlux::G4PSPassageCellFlux(G4String name, G4int depth)
51 :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0)
52{;}
53
54G4PSPassageCellFlux::~G4PSPassageCellFlux()
55{;}
56
57G4bool G4PSPassageCellFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
58{
59
60 if ( IsPassed(aStep) ) {
61 fCellFlux /=
62 aStep->GetPreStepPoint()->GetPhysicalVolume()
63 ->GetLogicalVolume()->GetSolid()->GetCubicVolume();
64 G4int index = GetIndex(aStep);
65 EvtMap->add(index,fCellFlux);
66 }
67
68 return TRUE;
69}
70
71G4bool G4PSPassageCellFlux::IsPassed(G4Step* aStep){
72 G4bool Passed = FALSE;
73
74 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
75 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
76
77 G4int trkid = aStep->GetTrack()->GetTrackID();
78 G4double trklength = aStep->GetStepLength();
79 trklength *= aStep->GetPreStepPoint()->GetWeight();
80
81 if ( IsEnter &&IsExit ){ // Passed at one step
82 fCellFlux = trklength; // Track length is absolutely given.
83 Passed = TRUE;
84 }else if ( IsEnter ){ // Enter a new geometry
85 fCurrentTrkID = trkid; // Resetting the current track.
86 fCellFlux = trklength;
87 }else if ( IsExit ){ // Exit a current geometry
88 if ( fCurrentTrkID == trkid ) {// if the track is same as entered,
89 fCellFlux += trklength; // add the track length to current one.
90 Passed = TRUE;
91 }
92 }else{ // Inside geometry
93 if ( fCurrentTrkID == trkid ){ // if the track is same as entered,
94 fCellFlux += trklength; // adding the track length to current one.
95 }
96 }
97
98 return Passed;
99}
100
101void G4PSPassageCellFlux::Initialize(G4HCofThisEvent* HCE)
102{
103 fCurrentTrkID = -1;
104
105 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
106 GetName());
107 if ( HCID < 0 ) HCID = GetCollectionID(0);
108 HCE->AddHitsCollection(HCID,EvtMap);
109
110}
111
112void G4PSPassageCellFlux::EndOfEvent(G4HCofThisEvent*)
113{;}
114
115void G4PSPassageCellFlux::clear(){
116 EvtMap->clear();
117}
118
119void G4PSPassageCellFlux::DrawAll()
120{;}
121
122void G4PSPassageCellFlux::PrintAll()
123{
124 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
125 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
126 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
127 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
128 for(; itr != EvtMap->GetMap()->end(); itr++) {
129 G4cout << " copy no.: " << itr->first
130 << " cell flux : " << *(itr->second)*cm*cm << " [cm^-2]"
131 << G4endl;
132 }
133}
134
Note: See TracBrowser for help on using the repository browser.