Geant4
11.4.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
// (Description)
33
// This is a primitive scorer class for scoring cell charge.
34
// The Cell Charge is defined by a sum of deposited charge inside the cell.
35
//
36
// Created: 2006-08-20 Tsukasa ASO
37
// 2010-07-22 Introduce Unit specification.
38
//
39
///////////////////////////////////////////////////////////////////////////////
40
41
G4PSCellCharge::G4PSCellCharge
(
const
G4String
& name,
G4int
depth)
42
:
G4PSCellCharge
(name,
"e+"
, depth)
43
{}
44
45
G4PSCellCharge::G4PSCellCharge
(
const
G4String
& name,
const
G4String
& unit,
G4int
depth)
46
:
G4VPrimitiveScorer
(name, depth)
47
, HCID(-1)
48
, EvtMap(nullptr)
49
{
50
SetUnit
(unit);
51
}
52
53
G4bool
G4PSCellCharge::ProcessHits
(
G4Step
* aStep,
G4TouchableHistory
*)
54
{
55
// Enter or First step of primary.
56
if
(aStep->
GetPreStepPoint
()->
GetStepStatus
() ==
fGeomBoundary
||
57
(aStep->
GetTrack
()->
GetParentID
() == 0 &&
58
aStep->
GetTrack
()->
GetCurrentStepNumber
() == 1))
59
{
60
G4double
CellCharge = aStep->
GetPreStepPoint
()->
GetCharge
();
61
CellCharge *= aStep->
GetPreStepPoint
()->
GetWeight
();
62
G4int
index =
GetIndex
(aStep);
63
EvtMap->add(index, CellCharge);
64
}
65
66
// Exit
67
if
(aStep->
GetPostStepPoint
()->
GetStepStatus
() ==
fGeomBoundary
)
68
{
69
G4double
CellCharge = aStep->
GetPreStepPoint
()->
GetCharge
();
70
CellCharge *= aStep->
GetPreStepPoint
()->
GetWeight
();
71
G4int
index =
GetIndex
(aStep);
72
CellCharge *= -1.0;
73
EvtMap->add(index, CellCharge);
74
}
75
76
return
true
;
77
}
78
79
void
G4PSCellCharge::Initialize
(
G4HCofThisEvent
* HCE)
80
{
81
EvtMap =
new
G4THitsMap<G4double>
(
detector
->GetName(),
GetName
());
82
if
(HCID < 0)
83
HCID =
GetCollectionID
(0);
84
HCE->
AddHitsCollection
(HCID, EvtMap);
85
}
86
87
void
G4PSCellCharge::clear
() { EvtMap->clear(); }
88
89
void
G4PSCellCharge::PrintAll
()
90
{
91
G4cout
<<
" MultiFunctionalDet "
<<
detector
->GetName() <<
G4endl
;
92
G4cout
<<
" PrimitiveScorer "
<<
GetName
() <<
G4endl
;
93
G4cout
<<
" Number of entries "
<< EvtMap->entries() <<
G4endl
;
94
for
(
const
auto
& [copy, charge] : *(EvtMap->GetMap()))
95
{
96
G4cout
<<
" copy no.: "
<< copy
97
<<
" cell charge : "
<< *(charge) /
GetUnitValue
() <<
" ["
98
<<
GetUnit
() <<
"]"
<<
G4endl
;
99
}
100
}
101
102
void
G4PSCellCharge::SetUnit
(
const
G4String
& unit)
103
{
104
CheckAndSetUnit
(unit,
"Electric charge"
);
105
}
G4PSCellCharge.hh
fGeomBoundary
@ fGeomBoundary
Definition
G4StepStatus.hh:43
G4Track.hh
G4double
double G4double
Definition
G4Types.hh:83
G4bool
bool G4bool
Definition
G4Types.hh:86
G4int
int G4int
Definition
G4Types.hh:85
G4endl
#define G4endl
Definition
G4ios.hh:67
G4cout
G4GLOB_DLL std::ostream G4cout
G4HCofThisEvent
Definition
G4HCofThisEvent.hh:51
G4HCofThisEvent::AddHitsCollection
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
Definition
G4HCofThisEvent.cc:52
G4PSCellCharge::G4PSCellCharge
G4PSCellCharge(const G4String &name, G4int depth=0)
Definition
G4PSCellCharge.cc:41
G4PSCellCharge::ProcessHits
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition
G4PSCellCharge.cc:53
G4PSCellCharge::Initialize
void Initialize(G4HCofThisEvent *) override
Definition
G4PSCellCharge.cc:79
G4PSCellCharge::SetUnit
virtual void SetUnit(const G4String &unit)
Definition
G4PSCellCharge.cc:102
G4PSCellCharge::clear
void clear() override
Definition
G4PSCellCharge.cc:87
G4PSCellCharge::PrintAll
void PrintAll() override
Definition
G4PSCellCharge.cc:89
G4StepPoint::GetStepStatus
G4StepStatus GetStepStatus() const
G4StepPoint::GetCharge
G4double GetCharge() const
G4StepPoint::GetWeight
G4double GetWeight() const
G4Step
Definition
G4Step.hh:61
G4Step::GetTrack
G4Track * GetTrack() const
G4Step::GetPreStepPoint
G4StepPoint * GetPreStepPoint() const
G4Step::GetPostStepPoint
G4StepPoint * GetPostStepPoint() const
G4String
Definition
G4String.hh:62
G4THitsMap
Definition
G4THitsMap.hh:521
G4TouchableHistory
G4TouchableHistory is an object representing a touchable detector element, and its history in the geo...
Definition
G4TouchableHistory.hh:107
G4Track::GetCurrentStepNumber
G4int GetCurrentStepNumber() const
G4Track::GetParentID
G4int GetParentID() const
G4VPrimitiveScorer::GetIndex
virtual G4int GetIndex(G4Step *)
Definition
G4VPrimitiveScorer.cc:60
G4VPrimitiveScorer::GetName
const G4String & GetName() const
Definition
G4VPrimitiveScorer.hh:89
G4VPrimitiveScorer::GetUnit
const G4String & GetUnit() const
Definition
G4VPrimitiveScorer.hh:77
G4VPrimitiveScorer::G4VPrimitiveScorer
G4VPrimitiveScorer(const G4String &name, G4int depth=0)
Definition
G4VPrimitiveScorer.cc:39
G4VPrimitiveScorer::detector
G4MultiFunctionalDetector * detector
Definition
G4VPrimitiveScorer.hh:126
G4VPrimitiveScorer::GetCollectionID
G4int GetCollectionID(G4int)
Definition
G4VPrimitiveScorer.cc:43
G4VPrimitiveScorer::CheckAndSetUnit
void CheckAndSetUnit(const G4String &unit, const G4String &category)
Definition
G4VPrimitiveScorer.cc:67
G4VPrimitiveScorer::GetUnitValue
G4double GetUnitValue() const
Definition
G4VPrimitiveScorer.hh:78
geant4-v11.4.0
source
digits_hits
scorer
src
G4PSCellCharge.cc
Generated by
1.16.1