Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PSCellCharge.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// G4PSCellCharge
28#include "G4PSCellCharge.hh"
29#include "G4Track.hh"
30
31
32///////////////////////////////////////////////////////////////////////////////
33// (Description)
34// This is a primitive scorer class for scoring cell charge.
35// The Cell Charge is defined by a sum of deposited charge inside the cell.
36//
37// Created: 2006-08-20 Tsukasa ASO
38// 2010-07-22 Introduce Unit specification.
39//
40///////////////////////////////////////////////////////////////////////////////
41
43 :G4VPrimitiveScorer(name,depth),HCID(-1),EvtMap(0)
44{
45 SetUnit("e+");
46}
47
49 G4int depth)
50 :G4VPrimitiveScorer(name,depth),HCID(-1),EvtMap(0)
51{
52 SetUnit(unit);
53}
54
56{;}
57
59{
60
61 // Enter or First step of primary.
63 || ( aStep->GetTrack()->GetParentID() == 0 &&
64 aStep->GetTrack()->GetCurrentStepNumber() == 1 ) ){
65 G4double CellCharge = aStep->GetPreStepPoint()->GetCharge();
66 CellCharge *= aStep->GetPreStepPoint()->GetWeight();
67 G4int index = GetIndex(aStep);
68 EvtMap->add(index,CellCharge);
69 }
70
71 // Exit
73 G4double CellCharge = aStep->GetPreStepPoint()->GetCharge();
74 CellCharge *= aStep->GetPreStepPoint()->GetWeight();
75 G4int index = GetIndex(aStep);
76 CellCharge *= -1.0;
77 EvtMap->add(index,CellCharge);
78 }
79
80 return TRUE;
81}
82
84{
86 GetName());
87 if ( HCID < 0 ) HCID = GetCollectionID(0);
88 HCE->AddHitsCollection(HCID,EvtMap);
89}
90
92{;}
93
95 EvtMap->clear();
96}
97
99{;}
100
102{
103 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
104 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
105 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
106 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
107 for(; itr != EvtMap->GetMap()->end(); itr++) {
108 G4cout << " copy no.: " << itr->first
109 << " cell charge : " << *(itr->second)/GetUnitValue()
110 << " [" << GetUnit() << "]"
111 << G4endl;
112 }
113}
114
116{
117 CheckAndSetUnit(unit,"Electric charge");
118}
119
120
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
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
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
virtual void Initialize(G4HCofThisEvent *)
virtual void EndOfEvent(G4HCofThisEvent *)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
G4PSCellCharge(G4String name, G4int depth=0)
virtual void PrintAll()
virtual ~G4PSCellCharge()
virtual void SetUnit(const G4String &unit)
virtual void clear()
virtual void DrawAll()
G4StepStatus GetStepStatus() const
G4double GetCharge() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int GetCurrentStepNumber() const
G4int GetParentID() 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
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