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

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

update ti head

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: G4PSCellFlux.cc,v 1.4 2010/07/22 23:42:01 taso Exp $
28// GEANT4 tag $Name: $
29//
30// G4PSCellFlux
31#include "G4PSCellFlux.hh"
32#include "G4Track.hh"
33#include "G4VSolid.hh"
34#include "G4VPhysicalVolume.hh"
35#include "G4VPVParameterisation.hh"
36#include "G4UnitsTable.hh"
37
38///////////////////////////////////////////////////////////////////////////////
39// (Description)
40// This is a primitive scorer class for scoring cell flux.
41// The Cell Flux is defined by a sum of track length divided
42// by the geometry volume, where all of the tracks in the geometry
43// are taken into account.
44//
45// If you want to score only tracks passing through the geometry volume,
46// please use G4PSPassageCellFlux.
47//
48//
49// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
50// 2010-07-22 Introduce Unit specification.
51// 2010-07-22 Add weighted option
52//
53///////////////////////////////////////////////////////////////////////////////
54
55G4PSCellFlux::G4PSCellFlux(G4String name, G4int depth)
56 :G4VPrimitiveScorer(name,depth),HCID(-1),weighted(true)
57{
58 DefineUnitAndCategory();
59 SetUnit("percm2");
60}
61
62G4PSCellFlux::G4PSCellFlux(G4String name, const G4String& unit, G4int depth)
63 :G4VPrimitiveScorer(name,depth),HCID(-1)
64{
65 DefineUnitAndCategory();
66 SetUnit(unit);
67}
68
69G4PSCellFlux::~G4PSCellFlux()
70{;}
71
72G4bool G4PSCellFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
73{
74 G4double stepLength = aStep->GetStepLength();
75 if ( stepLength == 0. ) return FALSE;
76
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 G4double CellFlux = stepLength / (solid->GetCubicVolume());
93 if (weighted) CellFlux *= aStep->GetPreStepPoint()->GetWeight();
94 G4int index = GetIndex(aStep);
95 EvtMap->add(index,CellFlux);
96
97 return TRUE;
98}
99
100void G4PSCellFlux::Initialize(G4HCofThisEvent* HCE)
101{
102 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
103 GetName());
104 if ( HCID < 0 ) HCID = GetCollectionID(0);
105 HCE->AddHitsCollection(HCID,EvtMap);
106}
107
108void G4PSCellFlux::EndOfEvent(G4HCofThisEvent*)
109{;}
110
111void G4PSCellFlux::clear(){
112 EvtMap->clear();
113}
114
115void G4PSCellFlux::DrawAll()
116{;}
117
118void G4PSCellFlux::PrintAll()
119{
120 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
121 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
122 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
123 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
124 for(; itr != EvtMap->GetMap()->end(); itr++) {
125 G4cout << " copy no.: " << itr->first
126 << " cell flux : " << *(itr->second)/GetUnitValue()
127 << " [" << GetUnit() << "]"
128 << G4endl;
129 }
130}
131
132void G4PSCellFlux::SetUnit(const G4String& unit)
133{
134 CheckAndSetUnit(unit,"Per Unit Surface");
135}
136
137void G4PSCellFlux::DefineUnitAndCategory(){
138 // Per Unit Surface
139 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
140 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
141 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
142}
Note: See TracBrowser for help on using the repository browser.