Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PSPassageTrackLength.cc
Go to the documentation of this file.
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$
28//
29// G4PSPassageTrackLength
31#include "G4StepStatus.hh"
32#include "G4Track.hh"
33#include "G4UnitsTable.hh"
34
35////////////////////////////////////////////////////////////////////////////////
36// (Description)
37// This is a primitive scorer class for scoring only track length.
38// The tracks which passed a geometry is taken into account.
39//
40//
41// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
42// 2010-07-22 Introduce Unit specification.
43//
44///////////////////////////////////////////////////////////////////////////////
45
47 :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
48 weighted(false)
49{
50 SetUnit("mm");
51}
52
54 const G4String& unit,
55 G4int depth)
56 :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
57 weighted(false)
58{
59 SetUnit(unit);
60}
61
62
64{;}
65
67{
68
69 if ( IsPassed(aStep) ) {
70 G4int index = GetIndex(aStep);
71 EvtMap->add(index,fTrackLength);
72 }
73
74 return TRUE;
75}
76
78 G4bool Passed = FALSE;
79
80 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
81 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
82
83 G4int trkid = aStep->GetTrack()->GetTrackID();
84 G4double trklength = aStep->GetStepLength();
85 if(weighted) trklength *= aStep->GetPreStepPoint()->GetWeight();
86
87 if ( IsEnter &&IsExit ){ // Passed at one step
88 fTrackLength = trklength; // Track length is absolutely given.
89 Passed = TRUE;
90 }else if ( IsEnter ){ // Enter a new geometry
91 fCurrentTrkID = trkid; // Resetting the current track.
92 fTrackLength = trklength;
93 }else if ( IsExit ){ // Exit a current geometry
94 if ( fCurrentTrkID == trkid ) {
95 fTrackLength += trklength; // Adding the track length to current one,
96 Passed = TRUE; // if the track is same as entered.
97 }
98 }else{ // Inside geometry
99 if ( fCurrentTrkID == trkid ){ // Adding the track length to current one ,
100 fTrackLength += trklength; // if the track is same as entered.
101 }
102 }
103
104 return Passed;
105}
106
108{
109 fCurrentTrkID = -1;
110
111 EvtMap = new G4THitsMap<G4double>(detector->GetName(),GetName());
112 if ( HCID < 0 ) HCID = GetCollectionID(0);
113 HCE->AddHitsCollection(HCID,EvtMap);
114}
115
117{;}
118
120 EvtMap->clear();
121}
122
124{;}
125
127{
128 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
129 G4cout << " PrimitiveSenstivity " << GetName() << G4endl;
130 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
131 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
132 for(; itr != EvtMap->GetMap()->end(); itr++) {
133 G4cout << " copy no.: " << itr->first
134 << " track length : "
135 << *(itr->second)/GetUnitValue()
136 << " [" << GetUnit()<< "]"
137 << G4endl;
138 }
139}
140
142{
143 CheckAndSetUnit(unit,"Length");
144}
@ fGeomBoundary
Definition: G4StepStatus.hh:54
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
virtual void SetUnit(const G4String &unit)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void EndOfEvent(G4HCofThisEvent *)
virtual void Initialize(G4HCofThisEvent *)
G4PSPassageTrackLength(G4String name, G4int depth=0)
G4StepStatus GetStepStatus() const
G4double GetWeight() const
Definition: G4Step.hh:78
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:138
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
G4int entries() const
Definition: G4THitsMap.hh:79
void clear()
Definition: G4THitsMap.hh:209
G4int GetTrackID() const
virtual G4int GetIndex(G4Step *)
G4String GetName() const
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
void CheckAndSetUnit(const G4String &unit, const G4String &category)
G4double GetUnitValue() const
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52