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
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// $Id$
27//
28// G4PSCellCharge
29#include "G4PSCellCharge.hh"
30#include "G4Track.hh"
31
32
33///////////////////////////////////////////////////////////////////////////////
34// (Description)
35// This is a primitive scorer class for scoring cell charge.
36// The Cell Charge is defined by a sum of deposited charge inside the cell.
37//
38// Created: 2006-08-20 Tsukasa ASO
39// 2010-07-22 Introduce Unit specification.
40//
41///////////////////////////////////////////////////////////////////////////////
42
44 :G4VPrimitiveScorer(name,depth),HCID(-1)
45{
46 SetUnit("e+");
47}
48
50 G4int depth)
51 :G4VPrimitiveScorer(name,depth),HCID(-1)
52{
53 SetUnit(unit);
54}
55
57{;}
58
60{
61
62 // Enter or First step of primary.
64 || ( aStep->GetTrack()->GetParentID() == 0 &&
65 aStep->GetTrack()->GetCurrentStepNumber() == 1 ) ){
66 G4double CellCharge = aStep->GetPreStepPoint()->GetCharge();
67 CellCharge *= aStep->GetPreStepPoint()->GetWeight();
68 G4int index = GetIndex(aStep);
69 EvtMap->add(index,CellCharge);
70 }
71
72 // Exit
74 G4double CellCharge = aStep->GetPreStepPoint()->GetCharge();
75 CellCharge *= aStep->GetPreStepPoint()->GetWeight();
76 G4int index = GetIndex(aStep);
77 CellCharge *= -1.0;
78 EvtMap->add(index,CellCharge);
79 }
80
81 return TRUE;
82}
83
85{
87 GetName());
88 if ( HCID < 0 ) HCID = GetCollectionID(0);
89 HCE->AddHitsCollection(HCID,EvtMap);
90}
91
93{;}
94
96 EvtMap->clear();
97}
98
100{;}
101
103{
104 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
105 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
106 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
107 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
108 for(; itr != EvtMap->GetMap()->end(); itr++) {
109 G4cout << " copy no.: " << itr->first
110 << " cell charge : " << *(itr->second)/GetUnitValue()
111 << " [" << GetUnit() << "]"
112 << G4endl;
113 }
114}
115
117{
118 CheckAndSetUnit(unit,"Electric charge");
119}
120
121
@ 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 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: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 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
#define TRUE
Definition: globals.hh:55