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

Last change on this file since 1110 was 1012, checked in by garnier, 17 years ago

update

File size: 4.3 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.2 2008/12/28 20:32:00 asaim Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
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//
51///////////////////////////////////////////////////////////////////////////////
52
53G4PSCellFlux::G4PSCellFlux(G4String name, G4int depth)
54 :G4VPrimitiveScorer(name,depth),HCID(-1)
55{;}
56
57G4PSCellFlux::~G4PSCellFlux()
58{;}
59
60G4bool G4PSCellFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
61{
62 G4double stepLength = aStep->GetStepLength();
63 if ( stepLength == 0. ) return FALSE;
64
65 G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
66 G4VPVParameterisation* physParam = physVol->GetParameterisation();
67 G4VSolid* solid = 0;
68 if(physParam)
69 { // for parameterized volume
70 G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
71 ->GetReplicaNumber(indexDepth);
72 solid = physParam->ComputeSolid(idx, physVol);
73 solid->ComputeDimensions(physParam,idx,physVol);
74 }
75 else
76 { // for ordinary volume
77 solid = physVol->GetLogicalVolume()->GetSolid();
78 }
79
80 G4double CellFlux = stepLength / (solid->GetCubicVolume());
81 CellFlux *= aStep->GetPreStepPoint()->GetWeight();
82 G4int index = GetIndex(aStep);
83 EvtMap->add(index,CellFlux);
84
85 return TRUE;
86}
87
88void G4PSCellFlux::Initialize(G4HCofThisEvent* HCE)
89{
90 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
91 GetName());
92 if ( HCID < 0 ) HCID = GetCollectionID(0);
93 HCE->AddHitsCollection(HCID,EvtMap);
94}
95
96void G4PSCellFlux::EndOfEvent(G4HCofThisEvent*)
97{;}
98
99void G4PSCellFlux::clear(){
100 EvtMap->clear();
101}
102
103void G4PSCellFlux::DrawAll()
104{;}
105
106void G4PSCellFlux::PrintAll()
107{
108 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
109 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
110 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
111 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
112 for(; itr != EvtMap->GetMap()->end(); itr++) {
113 G4cout << " copy no.: " << itr->first
114 << " cell flux : " << *(itr->second)*cm*cm << " [cm^-2]"
115 << G4endl;
116 }
117}
118
Note: See TracBrowser for help on using the repository browser.