Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4PSTrackCounter.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// G4PSTrackCounter
30#include "G4PSTrackCounter.hh"
31#include "G4UnitsTable.hh"
32
33///////////////////////////////////////////////////////////////////////////////
34// (Description)
35// This is a primitive scorer class for scoring number of tracks in a cell.
36//
37// Created: 2007-02-02 Tsukasa ASO, Akinori Kimura.
38// 2010-07-22 Introduce Unit specification.
39//
40///////////////////////////////////////////////////////////////////////////////
41
43 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),weighted(false)
44{
45 SetUnit("");
46}
47
49{;}
50
52{
53
54 G4StepPoint* preStep = aStep->GetPreStepPoint();
55 G4StepPoint* posStep = aStep->GetPostStepPoint();
56 G4bool IsEnter = preStep->GetStepStatus()==fGeomBoundary;
57 G4bool IsExit = posStep->GetStepStatus()==fGeomBoundary;
58
59 // Regular : count in prestep volume.
60 G4int index = GetIndex(aStep);
61
62 // G4cout << " trk " << aStep->GetTrack()->GetTrackID()
63 // << " index " << index << " In " << IsEnter << " Out " <<IsExit
64 // << " PreStatus " << preStep->GetStepStatus()
65 // << " PostStatus " << posStep->GetStepStatus()
66 // << " w " << aStep->GetPreStepPoint()->GetWeight()
67 // << G4endl;
68
69 G4bool flag = FALSE;
70
71 if ( IsEnter && fDirection == fCurrent_In ) flag = TRUE;
72 else if ( IsExit && fDirection == fCurrent_Out ) flag = TRUE;
73 else if ( (IsExit||IsEnter) && fDirection == fCurrent_InOut ) flag = TRUE;
74
75 if ( flag ){
76 //G4cout << " Count " << G4endl;
77 G4double val = 1.0;
78 if(weighted) val *= aStep->GetPreStepPoint()->GetWeight();
79 EvtMap->add(index,val);
80 }
81
82 return TRUE;
83}
84
86{
88 if(HCID < 0) {HCID = GetCollectionID(0);}
89 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
90 //G4cout << "----- PS Initialize " << G4endl;
91}
92
94{;}
95
97 EvtMap->clear();
98}
99
101{;}
102
104{
105 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
106 G4cout << " PrimitiveScorer " << GetName() << G4endl;
107 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
108 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
109 for(; itr != EvtMap->GetMap()->end(); itr++) {
110 G4cout << " copy no.: " << itr->first
111 << " track count: " << *(itr->second)
112 << " [tracks] "
113 << G4endl;
114 }
115}
116
118{
119 if (unit == "" ){
120 unitName = unit;
121 unitValue = 1.0;
122 }else{
123 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
124 G4Exception("G4PSTrackCounter::SetUnit","DetPS0018",JustWarning,msg);
125 }
126
127}
128
@ JustWarning
@ fCurrent_Out
@ fCurrent_In
@ fCurrent_InOut
@ 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 DrawAll()
G4PSTrackCounter(G4String name, G4int direction, G4int depth=0)
virtual void PrintAll()
virtual ~G4PSTrackCounter()
virtual void clear()
virtual void Initialize(G4HCofThisEvent *)
G4StepStatus GetStepStatus() const
G4double GetWeight() const
Definition: G4Step.hh:78
G4StepPoint * GetPreStepPoint() 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
virtual G4int GetIndex(G4Step *)
G4String GetName() const
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52