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
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// $Id$
28//
29// G4PSPassageCellCurrent
31#include "G4StepStatus.hh"
32#include "G4Track.hh"
33#include "G4VSolid.hh"
34#include "G4UnitsTable.hh"
35////////////////////////////////////////////////////////////////////////////////
36// (Description)
37// This is a primitive scorer class for scoring number of tracks,
38// where only tracks passing through the geometry are taken
39// into account.
40//
41// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
42// 2010-07-22 Introduce Unit specification.
43// 2010-07-22 Add weighted option
44//
45///////////////////////////////////////////////////////////////////////////////
46
48 :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCurrent(0),
49 weighted(true)
50{
51 SetUnit("");
52}
53
55{;}
56
58{
59
60 if ( IsPassed(aStep) ) {
61 if(weighted) fCurrent = aStep->GetPreStepPoint()->GetWeight();
62 G4int index = GetIndex(aStep);
63 EvtMap->add(index,fCurrent);
64 }
65
66 return TRUE;
67}
68
70 G4bool Passed = FALSE;
71
72 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
73 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
74
75 G4int trkid = aStep->GetTrack()->GetTrackID();
76
77 if ( IsEnter &&IsExit ){ // Passed at one step
78 Passed = TRUE;
79 }else if ( IsEnter ){ // Enter a new geometry
80 fCurrentTrkID = trkid; // Resetting the current track.
81 }else if ( IsExit ){ // Exit a current geometry
82 if ( fCurrentTrkID == trkid ) {
83 Passed = TRUE; // if the track is same as entered.
84 }
85 }else{ // Inside geometry
86 if ( fCurrentTrkID == trkid ){ // Adding the track length to current one ,
87 }
88 }
89 return Passed;
90}
91
93{
94 fCurrentTrkID = -1;
95
97 GetName());
98 if ( HCID < 0 ) HCID = GetCollectionID(0);
99 HCE->AddHitsCollection(HCID,EvtMap);
100
101}
102
104{
105}
106
108 EvtMap->clear();
109}
110
112{;}
113
115{
116 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
117 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
118 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
119 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
120 for(; itr != EvtMap->GetMap()->end(); itr++) {
121 G4cout << " copy no.: " << itr->first
122 << " cell current : " << *(itr->second)
123 << " [tracks] "
124 << G4endl;
125 }
126}
127
129{
130 if (unit == "" ){
131 unitName = unit;
132 unitValue = 1.0;
133 }else{
134 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
135 G4Exception("G4PSPassageCellCurrent::SetUnit","DetPS0012",JustWarning,msg);
136 }
137}
138
139
@ JustWarning
@ fGeomBoundary
Definition: G4StepStatus.hh:54
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 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 GetWeight() const
Definition: G4Step.hh:78
G4Track * GetTrack() const
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
G4int GetTrackID() const
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