Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0)
50{
51 SetUnit("Gy");
52}
53
55 G4int depth)
56 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0)
57{
58 SetUnit(unit);
59}
60
62{;}
63
65{
66 G4double edep = aStep->GetTotalEnergyDeposit();
67 if ( edep == 0. ) return FALSE;
68
70 (aStep->GetPreStepPoint()->GetTouchable()))
71 ->GetReplicaNumber(indexDepth);
72 G4double cubicVolume = ComputeVolume(aStep, idx);
73
74 G4double density = aStep->GetTrack()->GetStep()->GetPreStepPoint()->GetMaterial()->GetDensity();
75 G4double dose = edep / ( density * cubicVolume );
76 G4double wei = aStep->GetPreStepPoint()->GetWeight();
77 G4int index = GetIndex(aStep);
78 G4double dosew = dose*wei;
79 EvtMap->add(index,dosew);
80
81 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
82 {
83 auto filler = G4VScoreHistFiller::Instance();
84 if(!filler)
85 {
86 G4Exception("G4PSDoseDeposit::ProcessHits","SCORER0123",JustWarning,
87 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
88 }
89 else
90 {
91 filler->FillH1(hitIDMap[index],dose,wei);
92 }
93 }
94
95 return TRUE;
96}
97
99{
101 GetName());
102 if(HCID < 0) {HCID = GetCollectionID(0);}
103 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
104}
105
107{;}
108
110{
111 EvtMap->clear();
112}
113
115{;}
116
118{
119 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
120 G4cout << " PrimitiveScorer " << GetName() << G4endl;
121 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
122 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
123 for(; itr != EvtMap->GetMap()->end(); itr++) {
124 G4cout << " copy no.: " << itr->first
125 << " dose deposit: "
126 << *(itr->second)/GetUnitValue()
127 << " ["<<GetUnit() <<"]"
128 << G4endl;
129 }
130}
131
133{
134 CheckAndSetUnit(unit,"Dose");
135}
136
138{
139 G4VSolid* solid = ComputeSolid(aStep, idx);
140 return solid->GetCubicVolume();
141}
142
143
144
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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
#define FALSE
Definition: Globals.hh:23
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4double GetDensity() const
Definition: G4Material.hh:178
virtual G4double ComputeVolume(G4Step *, G4int idx)
virtual void DrawAll()
virtual void SetUnit(const G4String &unit)
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual void PrintAll()
virtual void EndOfEvent(G4HCofThisEvent *)
G4PSDoseDeposit(G4String name, G4int depth=0)
virtual void clear()
virtual ~G4PSDoseDeposit()
virtual void Initialize(G4HCofThisEvent *)
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:176
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