Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
36////////////////////////////////////////////////////////////////////////////////
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//
48///////////////////////////////////////////////////////////////////////////////
49
51 :G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fCurrent(0),
52 EvtMap(0),weighted(true)
53{
54 SetUnit("");
55}
56
58{;}
59
61{
62
63 if ( IsPassed(aStep) ) {
64 fCurrent = 1.;
65 if(weighted) fCurrent = aStep->GetPreStepPoint()->GetWeight();
66 G4int index = GetIndex(aStep);
67 EvtMap->add(index,fCurrent);
68
69 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
70 {
71 auto filler = G4VScoreHistFiller::Instance();
72 if(!filler)
73 {
74 G4Exception("G4PSVolumeFlux::ProcessHits","SCORER0123",JustWarning,
75 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
76 }
77 else
78 {
79 filler->FillH1(hitIDMap[index],aStep->GetPreStepPoint()->GetKineticEnergy(),fCurrent);
80 }
81 }
82
83 }
84
85 return TRUE;
86}
87
89 G4bool Passed = FALSE;
90
91 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
92 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
93
94 G4int trkid = aStep->GetTrack()->GetTrackID();
95
96 if ( IsEnter &&IsExit ){ // Passed at one step
97 Passed = TRUE;
98 }else if ( IsEnter ){ // Enter a new geometry
99 fCurrentTrkID = trkid; // Resetting the current track.
100 }else if ( IsExit ){ // Exit a current geometry
101 if ( fCurrentTrkID == trkid ) {
102 Passed = TRUE; // if the track is same as entered.
103 }
104 }else{ // Inside geometry
105 if ( fCurrentTrkID == trkid ){ // Adding the track length to current one ,
106 }
107 }
108 return Passed;
109}
110
112{
113 fCurrentTrkID = -1;
114
115 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
116 GetName());
117 if ( HCID < 0 ) HCID = GetCollectionID(0);
118 HCE->AddHitsCollection(HCID,EvtMap);
119
120}
121
123{
124}
125
127 EvtMap->clear();
128}
129
131{;}
132
134{
135 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
136 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
137 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
138 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
139 for(; itr != EvtMap->GetMap()->end(); itr++) {
140 G4cout << " copy no.: " << itr->first
141 << " cell current : " << *(itr->second)
142 << " [tracks] "
143 << G4endl;
144 }
145}
146
148{
149 if (unit == "" ){
150 unitName = unit;
151 unitValue = 1.0;
152 }else{
153 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
154 G4Exception("G4PSPassageCellCurrent::SetUnit","DetPS0012",JustWarning,msg);
155 }
156}
157
158
@ 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 *)
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:151
size_t entries() const
Definition: G4THitsMap.hh:155
void clear()
Definition: G4THitsMap.hh:537
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177