Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EnergySplitter.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#include "G4EnergySplitter.hh"
27
28#include "G4EmCalculator.hh"
30#include "G4PVParameterised.hh"
33#include "G4Step.hh"
34#include "G4UnitsTable.hh"
35#include "G4VSolid.hh"
36
37////////////////////////////////////////////////////////////////////////////////
38// (Description)
39//
40// Created:
41//
42///////////////////////////////////////////////////////////////////////////////
43
45{
46 theElossExt = new G4EnergyLossForExtrapolator(0);
47 thePhantomParam = nullptr;
48 theNIterations = 2;
49}
50
52{
53 delete theElossExt;
54}
55
57{
58 theEnergies.clear();
59
60 if (aStep == nullptr) return false; // it is 0 when called by GmScoringMgr after last event
61
62 G4double edep = aStep->GetTotalEnergyDeposit();
63
64#ifdef VERBOSE_ENERSPLIT
65 G4bool verbose = 1;
66 if (verbose)
67 G4cout << "G4EnergySplitter::SplitEnergyInVolumes totalEdepo " << aStep->GetTotalEnergyDeposit()
68 << " Nsteps " << G4RegularNavigationHelper::Instance()->GetStepLengths().size()
69 << G4endl;
70#endif
71 if (G4RegularNavigationHelper::Instance()->GetStepLengths().empty()
72 || aStep->GetTrack()->GetDefinition()->GetPDGCharge() == 0)
73 { // we are only counting dose deposit
74 return (G4int)theEnergies.size();
75 }
76 if (G4RegularNavigationHelper::Instance()->GetStepLengths().size() == 1) {
77 theEnergies.push_back(edep);
78 return (G4int)theEnergies.size();
79 }
80
81 if (thePhantomParam == nullptr) GetPhantomParam(true);
82
83 //----- Distribute energy deposited in voxels
84 std::vector<std::pair<G4int, G4double>> rnsl =
86
87 const G4ParticleDefinition* part = aStep->GetTrack()->GetDefinition();
88 G4double kinEnergyPreOrig = aStep->GetPreStepPoint()->GetKineticEnergy();
89 G4double kinEnergyPre = kinEnergyPreOrig;
90
91 G4double stepLength = aStep->GetStepLength();
92 G4double slSum = 0.;
93 unsigned int ii;
94 for (ii = 0; ii < rnsl.size(); ++ii) {
95 G4double sl = rnsl[ii].second;
96 slSum += sl;
97#ifdef VERBOSE_ENERSPLIT
98 if (verbose)
99 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter1 step length geom "
100 << sl << G4endl;
101#endif
102 }
103
104#ifdef VERBOSE_ENERSPLIT
105 if (verbose)
106 G4cout << "G4EnergySplitter RN: step length geom TOTAL " << slSum << " true TOTAL "
107 << stepLength << " ratio " << stepLength / slSum << " Energy "
108 << aStep->GetPreStepPoint()->GetKineticEnergy() << " Material "
109 << aStep->GetPreStepPoint()->GetMaterial()->GetName() << " Number of geom steps "
110 << rnsl.size() << G4endl;
111#endif
112 //----- No iterations to correct elost and msc => distribute energy deposited according to
113 // geometrical step length in each voxel
114 if (theNIterations == 0) {
115 for (ii = 0; ii < rnsl.size(); ++ii) {
116 G4double sl = rnsl[ii].second;
117 G4double edepStep = edep * sl / slSum; // divide edep along steps, proportional to step
118 // length
119#ifdef VERBOSE_ENERSPLIT
120 if (verbose)
121 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " edep " << edepStep << G4endl;
122#endif
123
124 theEnergies.push_back(edepStep);
125 }
126 }
127 else { // 1 or more iterations demanded
128
129#ifdef VERBOSE_ENERSPLIT
130 // print corrected energy at iteration 0
131 if (verbose) {
132 G4double slSum = 0.;
133 for (ii = 0; ii < rnsl.size(); ++ii) {
134 G4double sl = rnsl[ii].second;
135 slSum += sl;
136 }
137 for (ii = 0; ii < rnsl.size(); ii++) {
138 G4cout << "G4EnergySplitter::SplitEnergyInVolumes " << ii
139 << " RN: iter0 corrected energy lost " << edep * rnsl[ii].second / slSum << G4endl;
140 }
141 }
142#endif
143
144 G4double slRatio = stepLength / slSum;
145#ifdef VERBOSE_ENERSPLIT
146 if (verbose)
147 G4cout << "G4EnergySplitter::SplitEnergyInVolumes RN: iter 0, step ratio " << slRatio
148 << G4endl;
149#endif
150
151 //--- energy at each interaction
152 G4EmCalculator emcalc;
153 G4double totalELost = 0.;
154 std::vector<G4double> stepLengths;
155 for (G4int iiter = 1; iiter <= theNIterations; ++iiter) {
156 //--- iter1: distribute true step length in each voxel: geom SL in each voxel is multiplied by
157 // a constant so that the sum gives the total true step length
158 if (iiter == 1) {
159 for (ii = 0; ii < rnsl.size(); ++ii) {
160 G4double sl = rnsl[ii].second;
161 stepLengths.push_back(sl * slRatio);
162#ifdef VERBOSE_ENERSPLIT
163 if (verbose)
164 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter" << iiter
165 << " corrected step length " << sl * slRatio << G4endl;
166#endif
167 }
168
169 for (ii = 0; ii < rnsl.size(); ++ii) {
170 const G4Material* mate = thePhantomParam->GetMaterial(rnsl[ii].first);
171 G4double dEdx = 0.;
172 if (kinEnergyPre > 0.) { // t check this
173 dEdx = emcalc.GetDEDX(kinEnergyPre, part, mate);
174 }
175 G4double elost = stepLengths[ii] * dEdx;
176
177#ifdef VERBOSE_ENERSPLIT
178 if (verbose)
179 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter1 energy lost "
180 << elost << " energy at interaction " << kinEnergyPre << " = stepLength "
181 << stepLengths[ii] << " * dEdx " << dEdx << G4endl;
182#endif
183 kinEnergyPre -= elost;
184 theEnergies.push_back(elost);
185 totalELost += elost;
186 }
187 }
188 else {
189 //------ 2nd and other iterations
190 //----- Get step lengths corrected by changing geom2true correction
191 //-- Get ratios for each energy
192 slSum = 0.;
193 kinEnergyPre = kinEnergyPreOrig;
194 for (ii = 0; ii < rnsl.size(); ++ii) {
195 const G4Material* mate = thePhantomParam->GetMaterial(rnsl[ii].first);
196 stepLengths[ii] = theElossExt->TrueStepLength(kinEnergyPre, rnsl[ii].second, mate, part);
197 kinEnergyPre -= theEnergies[ii];
198
199#ifdef VERBOSE_ENERSPLIT
200 if (verbose)
201 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter" << iiter
202 << " step length geom " << stepLengths[ii] << " geom2true "
203 << rnsl[ii].second / stepLengths[ii] << G4endl;
204#endif
205
206 slSum += stepLengths[ii];
207 }
208
209 // Correct step lengths so that they sum the total step length
210 G4double slratio = aStep->GetStepLength() / slSum;
211#ifdef VERBOSE_ENERSPLIT
212 if (verbose)
213 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter" << iiter
214 << " step ratio " << slRatio << G4endl;
215#endif
216 for (ii = 0; ii < rnsl.size(); ++ii) {
217 stepLengths[ii] *= slratio;
218#ifdef VERBOSE_ENERSPLIT
219 if (verbose)
220 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter" << iiter
221 << " corrected step length " << stepLengths[ii] << G4endl;
222#endif
223 }
224
225 //---- Recalculate energy lost with this new step lengths
226 kinEnergyPre = aStep->GetPreStepPoint()->GetKineticEnergy();
227 totalELost = 0.;
228 for (ii = 0; ii < rnsl.size(); ++ii) {
229 const G4Material* mate = thePhantomParam->GetMaterial(rnsl[ii].first);
230 G4double dEdx = 0.;
231 if (kinEnergyPre > 0.) {
232 dEdx = emcalc.GetDEDX(kinEnergyPre, part, mate);
233 }
234 G4double elost = stepLengths[ii] * dEdx;
235#ifdef VERBOSE_ENERSPLIT
236 if (verbose)
237 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter" << iiter
238 << " energy lost " << elost << " energy at interaction " << kinEnergyPre
239 << " = stepLength " << stepLengths[ii] << " * dEdx " << dEdx << G4endl;
240#endif
241 kinEnergyPre -= elost;
242 theEnergies[ii] = elost;
243 totalELost += elost;
244 }
245 }
246
247 // correct energies so that they reproduce the real step energy lost
248 G4double enerRatio = (edep / totalELost);
249
250#ifdef VERBOSE_ENERSPLIT
251 if (verbose)
252 G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter" << iiter
253 << " energy ratio " << enerRatio << G4endl;
254#endif
255
256#ifdef VERBOSE_ENERSPLIT
257 G4double elostTot = 0.;
258#endif
259 for (ii = 0; ii < theEnergies.size(); ++ii) {
260 theEnergies[ii] *= enerRatio;
261#ifdef VERBOSE_ENERSPLIT
262 elostTot += theEnergies[ii];
263 if (verbose)
264 G4cout << "G4EnergySplitter::SplitEnergyInVolumes " << ii << " RN: iter" << iiter
265 << " corrected energy lost " << theEnergies[ii] << " orig elost "
266 << theEnergies[ii] / enerRatio << " energy before interaction "
267 << kinEnergyPreOrig - elostTot + theEnergies[ii] << " energy after interaction "
268 << kinEnergyPreOrig - elostTot << G4endl;
269#endif
270 }
271 }
272 }
273
274 return (G4int)theEnergies.size();
275}
276
277//-----------------------------------------------------------------------
278void G4EnergySplitter::GetPhantomParam(G4bool mustExist)
279{
281 for (const auto pv : *pvs) {
282 if (IsPhantomVolume(pv)) {
283 const auto pvparam = static_cast<const G4PVParameterised*>(pv);
284 G4VPVParameterisation* param = pvparam->GetParameterisation();
285 thePhantomParam = static_cast<G4PhantomParameterisation*>(param);
286 }
287 }
288
289 if ((thePhantomParam == nullptr) && mustExist)
290 G4Exception("G4EnergySplitter::GetPhantomParam", "PhantomParamError", FatalException,
291 "No G4PhantomParameterisation found !");
292}
293
294//-----------------------------------------------------------------------
295G4bool G4EnergySplitter::IsPhantomVolume(G4VPhysicalVolume* pv)
296{
297 EAxis axis;
298 G4int nReplicas;
299 G4double width, offset;
300 G4bool consuming;
301 pv->GetReplicationData(axis, nReplicas, width, offset, consuming);
302 EVolume type = (consuming) ? kReplica : kParameterised;
303 return type == kParameterised && pv->GetRegularStructureId() == 1;
304}
305
306//-----------------------------------------------------------------------
308{
309 voxelID = (*(G4RegularNavigationHelper::Instance()->GetStepLengths().cbegin())).first;
310}
311
312//-----------------------------------------------------------------------
314{
315 voxelID = (*(G4RegularNavigationHelper::Instance()->GetStepLengths().crbegin())).first;
316}
317
318//-----------------------------------------------------------------------
320{
321 if (stepNo < 0 || stepNo >= G4int(G4RegularNavigationHelper::Instance()->GetStepLengths().size()))
322 {
323 G4Exception("G4EnergySplitter::GetVoxelID",
324 "Invalid stepNo, smaller than 0 or bigger or equal to number of voxels traversed",
326 G4String("stepNo = " + G4UIcommand::ConvertToString(stepNo)
327 + ", number of voxels = "
329 G4int(G4RegularNavigationHelper::Instance()->GetStepLengths().size())))
330 .c_str());
331 }
333 advance(ite, stepNo);
334 voxelID = (*ite).first;
335}
336
337//-----------------------------------------------------------------------
338void G4EnergySplitter::GetStepLength(G4int stepNo, G4double& stepLength)
339{
341 advance(ite, stepNo);
342 stepLength = (*ite).second;
343}
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=nullptr)
G4double TrueStepLength(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
virtual ~G4EnergySplitter()
void GetFirstVoxelID(G4int &voxelID)
G4int SplitEnergyInVolumes(const G4Step *aStep)
void GetLastVoxelID(G4int &voxelID)
void GetVoxelID(G4int stepNo, G4int &voxelID)
const G4String & GetName() const
G4Material * GetMaterial(std::size_t nx, std::size_t ny, std::size_t nz) const
static G4PhysicalVolumeStore * GetInstance()
const std::vector< std::pair< G4int, G4double > > & GetStepLengths()
static G4RegularNavigationHelper * Instance()
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
G4ParticleDefinition * GetDefinition() const
static G4String ConvertToString(G4bool boolVal)
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
virtual G4int GetRegularStructureId() const =0
EAxis
Definition geomdefs.hh:54
EVolume
Definition geomdefs.hh:83
@ kParameterised
Definition geomdefs.hh:86
@ kReplica
Definition geomdefs.hh:85