Geant4-11
G4PSPassageCellCurrent.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//
28// G4PSPassageCellCurrent
30#include "G4StepStatus.hh"
31#include "G4Track.hh"
32#include "G4VSolid.hh"
33#include "G4UnitsTable.hh"
34#include "G4VScoreHistFiller.hh"
35
37// (Description)
38// This is a primitive scorer class for scoring number of tracks,
39// where only tracks passing through the geometry are taken
40// into account.
41//
42// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
43// 2010-07-22 Introduce Unit specification.
44// 2010-07-22 Add weighted option
45// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x)
46// vs. current * track weight (y) (Makoto Asai)
47//
49
52 , HCID(-1)
53 , fCurrentTrkID(-1)
54 , fCurrent(0)
55 , EvtMap(0)
56 , weighted(true)
57{
58 SetUnit("");
59}
60
62
64{
65 if(IsPassed(aStep))
66 {
67 fCurrent = 1.;
68 if(weighted)
70 G4int index = GetIndex(aStep);
71 EvtMap->add(index, fCurrent);
72
73 if(hitIDMap.size() > 0 && hitIDMap.find(index) != hitIDMap.end())
74 {
75 auto filler = G4VScoreHistFiller::Instance();
76 if(!filler)
77 {
79 "G4PSVolumeFlux::ProcessHits", "SCORER0123", JustWarning,
80 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
81 }
82 else
83 {
84 filler->FillH1(hitIDMap[index],
86 }
87 }
88 }
89
90 return TRUE;
91}
92
94{
95 G4bool Passed = FALSE;
96
97 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
98 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
99
100 G4int trkid = aStep->GetTrack()->GetTrackID();
101
102 if(IsEnter && IsExit)
103 { // Passed at one step
104 Passed = TRUE;
105 }
106 else if(IsEnter)
107 { // Enter a new geometry
108 fCurrentTrkID = trkid; // Resetting the current track.
109 }
110 else if(IsExit)
111 { // Exit a current geometry
112 if(fCurrentTrkID == trkid)
113 {
114 Passed = TRUE; // if the track is same as entered.
115 }
116 }
117 else
118 { // Inside geometry
119 if(fCurrentTrkID == trkid)
120 { // Adding the track length to current one ,
121 }
122 }
123 return Passed;
124}
125
127{
128 fCurrentTrkID = -1;
129
131 if(HCID < 0)
134}
135
137
139
141
143{
144 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
145 G4cout << " PrimitiveScorer " << GetName() << G4endl;
146 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
147 std::map<G4int, G4double*>::iterator itr = EvtMap->GetMap()->begin();
148 for(; itr != EvtMap->GetMap()->end(); itr++)
149 {
150 G4cout << " copy no.: " << itr->first
151 << " cell current : " << *(itr->second) << " [tracks] " << G4endl;
152 }
153}
154
156{
157 if(unit == "")
158 {
159 unitName = unit;
160 unitValue = 1.0;
161 }
162 else
163 {
164 G4String msg = "Invalid unit [" + unit + "] (Current unit is [" +
165 GetUnit() + "] ) for " + GetName();
166 G4Exception("G4PSPassageCellCurrent::SetUnit", "DetPS0012", JustWarning,
167 msg);
168 }
169}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fGeomBoundary
Definition: G4StepStatus.hh:43
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
virtual void SetUnit(const G4String &unit)
virtual G4bool IsPassed(G4Step *)
G4PSPassageCellCurrent(G4String name, G4int depth=0)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
G4THitsMap< G4double > * EvtMap
virtual void EndOfEvent(G4HCofThisEvent *)
virtual void Initialize(G4HCofThisEvent *)
G4StepStatus GetStepStatus() const
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
static G4VScoreHistFiller * Instance()
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t entries() const
Definition: G4THitsMap.hh:158
void clear()
Definition: G4THitsMap.hh:524
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177
const char * name(G4int ptype)