Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PSDoseDeposit.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// G4PSDoseDeposit
29#include "G4PSDoseDeposit.hh"
30#include "G4VSolid.hh"
31#include "G4VPhysicalVolume.hh"
33#include "G4UnitsTable.hh"
34#include "G4VScoreHistFiller.hh"
35
36////////////////////////////////////////////////////////////////////////////////
37// (Description)
38// This is a primitive scorer class for scoring only energy deposit.
39//
40//
41// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
42// 2010-07-22 Introduce Unit specification.
43// 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of dose deposit (x)
44// vs. weight (y) (Makoto Asai)
45//
46///////////////////////////////////////////////////////////////////////////////
47
49 : G4PSDoseDeposit(name, "Gy", depth)
50{}
51
53 G4int depth)
54 : G4VPrimitivePlotter(name, depth)
55 , HCID(-1)
56 , EvtMap(nullptr)
57{
58 SetUnit(unit);
59}
60
62{
63 G4double edep = aStep->GetTotalEnergyDeposit();
64 if(edep == 0.)
65 return false;
66
68 ->GetReplicaNumber(indexDepth);
69 G4double cubicVolume = ComputeVolume(aStep, idx);
70
71 G4double density = aStep->GetTrack()
72 ->GetStep()
74 ->GetMaterial()
75 ->GetDensity();
76 G4double dose = edep / (density * cubicVolume);
77 G4double wei = aStep->GetPreStepPoint()->GetWeight();
78 G4int index = GetIndex(aStep);
79 G4double dosew = dose * wei;
80 EvtMap->add(index, dosew);
81
82 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend())
83 {
84 auto filler = G4VScoreHistFiller::Instance();
85 if(filler == nullptr)
86 {
88 "G4PSDoseDeposit::ProcessHits", "SCORER0123", JustWarning,
89 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
90 }
91 else
92 {
93 filler->FillH1(hitIDMap[index], dose, wei);
94 }
95 }
96
97 return true;
98}
99
101{
103 GetName());
104 if(HCID < 0)
105 {
106 HCID = GetCollectionID(0);
107 }
108 HCE->AddHitsCollection(HCID, (G4VHitsCollection*) EvtMap);
109}
110
111void G4PSDoseDeposit::clear() { EvtMap->clear(); }
112
114{
115 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
116 G4cout << " PrimitiveScorer " << GetName() << G4endl;
117 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
118 for(const auto& [copy, dose] : *(EvtMap->GetMap()))
119 {
120 G4cout << " copy no.: " << copy
121 << " dose deposit: " << *(dose) / GetUnitValue() << " ["
122 << GetUnit() << "]" << G4endl;
123 }
124}
125
127{
128 CheckAndSetUnit(unit, "Dose");
129}
130
132{
133 G4VSolid* solid = ComputeSolid(aStep, idx);
134 return solid->GetCubicVolume();
135}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
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
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4double GetDensity() const
Definition: G4Material.hh:175
virtual G4double ComputeVolume(G4Step *, G4int idx)
virtual void SetUnit(const G4String &unit)
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
void Initialize(G4HCofThisEvent *) override
G4PSDoseDeposit(G4String name, G4int depth=0)
void PrintAll() override
void clear() override
const G4VTouchable * GetTouchable() const
G4Material * GetMaterial() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
const G4Step * GetStep() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
G4String GetName() const
G4MultiFunctionalDetector * GetMultiFunctionalDetector() const
G4VSolid * ComputeSolid(G4Step *aStep, G4int replicaIdx)
const G4String & GetUnit() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)
void CheckAndSetUnit(const G4String &unit, const G4String &category)
G4double GetUnitValue() const
static G4VScoreHistFiller * Instance()
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
size_t entries() const
Definition: G4THitsMap.hh:158
void clear()
Definition: G4THitsMap.hh:524
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177