Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EnergySplitter Class Reference

#include <G4EnergySplitter.hh>

Public Member Functions

 G4EnergySplitter ()
 
virtual ~G4EnergySplitter ()
 
G4int SplitEnergyInVolumes (const G4Step *aStep)
 
void GetLastVoxelID (G4int &voxelID)
 
void GetFirstVoxelID (G4int &voxelID)
 
void GetVoxelID (G4int stepNo, G4int &voxelID)
 
void GetVoxelIDAndLength (G4int stepNo, G4int &voxelID, G4double &stepLength)
 
void GetLengthAndEnergyDeposited (G4int stepNo, G4int &voxelID, G4double &stepLength, G4double &energyLoss)
 
void GetLengthAndInitialEnergy (G4double &preStepEnergy, G4int stepNo, G4int &voxelID, G4double &stepLength, G4double &initialEnergy)
 
void SetNIterations (G4int niter)
 
G4MaterialGetVoxelMaterial (G4int stepNo)
 

Detailed Description

Definition at line 46 of file G4EnergySplitter.hh.

Constructor & Destructor Documentation

◆ G4EnergySplitter()

G4EnergySplitter::G4EnergySplitter ( )

Definition at line 43 of file G4EnergySplitter.cc.

44{
45 theElossExt = new G4EnergyLossForExtrapolator(0);
46 thePhantomParam = 0;
47 theNIterations = 2;
48}

◆ ~G4EnergySplitter()

G4EnergySplitter::~G4EnergySplitter ( )
virtual

Definition at line 50 of file G4EnergySplitter.cc.

51{;}

Member Function Documentation

◆ GetFirstVoxelID()

void G4EnergySplitter::GetFirstVoxelID ( G4int voxelID)

Definition at line 306 of file G4EnergySplitter.cc.

307{
308 voxelID = (*(G4RegularNavigationHelper::Instance()->GetStepLengths().rbegin())).first;
309}
const std::vector< std::pair< G4int, G4double > > & GetStepLengths()
static G4RegularNavigationHelper * Instance()

◆ GetLastVoxelID()

void G4EnergySplitter::GetLastVoxelID ( G4int voxelID)

Definition at line 300 of file G4EnergySplitter.cc.

301{
302 voxelID = (*(G4RegularNavigationHelper::Instance()->GetStepLengths().begin())).first;
303}

◆ GetLengthAndEnergyDeposited()

void G4EnergySplitter::GetLengthAndEnergyDeposited ( G4int  stepNo,
G4int voxelID,
G4double stepLength,
G4double energyLoss 
)
inline

◆ GetLengthAndInitialEnergy()

void G4EnergySplitter::GetLengthAndInitialEnergy ( G4double preStepEnergy,
G4int  stepNo,
G4int voxelID,
G4double stepLength,
G4double initialEnergy 
)
inline

◆ GetVoxelID()

void G4EnergySplitter::GetVoxelID ( G4int  stepNo,
G4int voxelID 
)

Definition at line 312 of file G4EnergySplitter.cc.

313{
314 if( stepNo < 0 || stepNo >= G4int(G4RegularNavigationHelper::Instance()->GetStepLengths().size()) ) {
315 G4Exception("G4EnergySplitter::GetVoxelID",
316 "Invalid stepNo, smaller than 0 or bigger or equal to number of voxels traversed",
318 G4String("stepNo = " + G4UIcommand::ConvertToString(stepNo) + ", number of voxels = " + G4UIcommand::ConvertToString(G4int(G4RegularNavigationHelper::Instance()->GetStepLengths().size())) ).c_str());
319 }
320 std::vector< std::pair<G4int,G4double> >::const_iterator ite = G4RegularNavigationHelper::Instance()->GetStepLengths().begin();
321 advance( ite, stepNo );
322 voxelID = (*ite).first;
323
324}
@ FatalErrorInArgument
int G4int
Definition: G4Types.hh:66
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:349
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4ScoreSplittingProcess::PostStepDoIt().

◆ GetVoxelIDAndLength()

void G4EnergySplitter::GetVoxelIDAndLength ( G4int  stepNo,
G4int voxelID,
G4double stepLength 
)
inline

◆ GetVoxelMaterial()

G4Material * G4EnergySplitter::GetVoxelMaterial ( G4int  stepNo)
inline

◆ SetNIterations()

void G4EnergySplitter::SetNIterations ( G4int  niter)
inline

◆ SplitEnergyInVolumes()

G4int G4EnergySplitter::SplitEnergyInVolumes ( const G4Step aStep)

Definition at line 53 of file G4EnergySplitter.cc.

54{
55 theEnergies.clear();
56
57 G4double edep = aStep->GetTotalEnergyDeposit();
58
59#ifdef VERBOSE_ENERSPLIT
60 G4bool verbose = 1;
61 if( verbose ) G4cout << "G4EnergySplitter::SplitEnergyInVolumes totalEdepo " << aStep->GetTotalEnergyDeposit()
62 << " Nsteps " << G4RegularNavigationHelper::Instance()->GetStepLengths().size() << G4endl;
63#endif
64 if( G4RegularNavigationHelper::Instance()->GetStepLengths().size() == 0 ||
65 aStep->GetTrack()->GetDefinition()->GetPDGCharge() == 0) { // we are only counting dose deposit
66 return theEnergies.size();
67 }
68 if( G4RegularNavigationHelper::Instance()->GetStepLengths().size() == 1 ) {
69 theEnergies.push_back(edep);
70 return theEnergies.size();
71 }
72
73 if( !thePhantomParam ) GetPhantomParam(TRUE);
74
75 if( aStep == 0 ) return FALSE; // it is 0 when called by GmScoringMgr after last event
76
77 //----- Distribute energy deposited in voxels
78 std::vector< std::pair<G4int,G4double> > rnsl = G4RegularNavigationHelper::Instance()->GetStepLengths();
79
80 const G4ParticleDefinition* part = aStep->GetTrack()->GetDefinition();
81 G4double kinEnergyPreOrig = aStep->GetPreStepPoint()->GetKineticEnergy();
82 G4double kinEnergyPre = kinEnergyPreOrig;
83
84 G4double stepLength = aStep->GetStepLength();
85 G4double slSum = 0.;
86 unsigned int ii;
87 for( ii = 0; ii < rnsl.size(); ii++ ){
88 G4double sl = rnsl[ii].second;
89 slSum += sl;
90#ifdef VERBOSE_ENERSPLIT
91 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes"<< ii << " RN: iter1 step length geom " << sl << G4endl;
92#endif
93 }
94
95#ifdef VERBOSE_ENERSPLIT
96 if( verbose )
97 G4cout << "G4EnergySplitter RN: step length geom TOTAL " << slSum
98 << " true TOTAL " << stepLength
99 << " ratio " << stepLength/slSum
100 << " Energy " << aStep->GetPreStepPoint()->GetKineticEnergy()
101 << " Material " << aStep->GetPreStepPoint()->GetMaterial()->GetName()
102 << " Number of geom steps " << rnsl.size() << G4endl;
103#endif
104 //----- No iterations to correct elost and msc => distribute energy deposited according to geometrical step length in each voxel
105 if( theNIterations == 0 ) {
106 for( ii = 0; ii < rnsl.size(); ii++ ){
107 G4double sl = rnsl[ii].second;
108 G4double edepStep = edep * sl/slSum; //divide edep along steps, proportional to step length
109#ifdef VERBOSE_ENERSPLIT
110 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes"<< ii
111 << " edep " << edepStep << G4endl;
112#endif
113
114 theEnergies.push_back(edepStep);
115
116 }
117 } else { // 1 or more iterations demanded
118
119#ifdef VERBOSE_ENERSPLIT
120 // print corrected energy at iteration 0
121 if(verbose) {
122 G4double slSum = 0.;
123 for( ii = 0; ii < rnsl.size(); ii++ ){
124 G4double sl = rnsl[ii].second;
125 slSum += sl;
126 }
127 for( ii = 0; ii < rnsl.size(); ii++ ){
128 G4cout << "G4EnergySplitter::SplitEnergyInVolumes "<< ii
129 << " RN: iter0 corrected energy lost " << edep*rnsl[ii].second/slSum
130 << G4endl;
131 }
132 }
133#endif
134
135 G4double slRatio = stepLength/slSum;
136#ifdef VERBOSE_ENERSPLIT
137 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes RN: iter 0, step ratio " << slRatio << G4endl;
138#endif
139
140 //--- energy at each interaction
141 G4EmCalculator emcalc;
142 G4double totalELost = 0.;
143 std::vector<G4double> stepLengths;
144 for( int iiter = 1; iiter <= theNIterations; iiter++ ) {
145 //--- iter1: distribute true step length in each voxel: geom SL in each voxel is multiplied by a constant so that the sum gives the total true step length
146 if( iiter == 1 ) {
147 for( ii = 0; ii < rnsl.size(); ii++ ){
148 G4double sl = rnsl[ii].second;
149 stepLengths.push_back( sl * slRatio );
150#ifdef VERBOSE_ENERSPLIT
151 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes"<< ii << " RN: iter" << iiter << " corrected step length " << sl*slRatio << G4endl;
152#endif
153 }
154
155 for( ii = 0; ii < rnsl.size(); ii++ ){
156 const G4Material* mate = thePhantomParam->GetMaterial( rnsl[ii].first );
157 G4double dEdx = 0.;
158 if( kinEnergyPre > 0. ) { //t check this
159 dEdx = emcalc.GetDEDX(kinEnergyPre, part, mate);
160 }
161 G4double elost = stepLengths[ii] * dEdx;
162
163#ifdef VERBOSE_ENERSPLIT
164 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes"<< ii << " RN: iter1 energy lost " << elost
165 << " energy at interaction " << kinEnergyPre
166 << " = stepLength " << stepLengths[ii]
167 << " * dEdx " << dEdx << G4endl;
168#endif
169 kinEnergyPre -= elost;
170 theEnergies.push_back( elost );
171 totalELost += elost;
172 }
173
174 } else{
175 //------ 2nd and other iterations
176 //----- Get step lengths corrected by changing geom2true correction
177 //-- Get ratios for each energy
178 slSum = 0.;
179 kinEnergyPre = kinEnergyPreOrig;
180 for( ii = 0; ii < rnsl.size(); ii++ ){
181 const G4Material* mate = thePhantomParam->GetMaterial( rnsl[ii].first );
182 stepLengths[ii] = theElossExt->TrueStepLength( kinEnergyPre, rnsl[ii].second , mate, part );
183 kinEnergyPre -= theEnergies[ii];
184
185#ifdef VERBOSE_ENERSPLIT
186 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii
187 << " RN: iter" << iiter << " step length geom " << stepLengths[ii]
188 << " geom2true " << rnsl[ii].second / stepLengths[ii] << G4endl;
189#endif
190
191 slSum += stepLengths[ii];
192 }
193
194 //Correct step lengths so that they sum the total step length
195 G4double slratio = aStep->GetStepLength()/slSum;
196#ifdef VERBOSE_ENERSPLIT
197 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes" << ii << " RN: iter" << iiter << " step ratio " << slRatio << G4endl;
198#endif
199 for( ii = 0; ii < rnsl.size(); ii++ ){
200 stepLengths[ii] *= slratio;
201#ifdef VERBOSE_ENERSPLIT
202 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes"<< ii << " RN: iter" << iiter << " corrected step length " << stepLengths[ii] << G4endl;
203#endif
204 }
205
206 //---- Recalculate energy lost with this new step lengths
207 kinEnergyPre = aStep->GetPreStepPoint()->GetKineticEnergy();
208 totalELost = 0.;
209 for( ii = 0; ii < rnsl.size(); ii++ ){
210 const G4Material* mate = thePhantomParam->GetMaterial( rnsl[ii].first );
211 G4double dEdx = 0.;
212 if( kinEnergyPre > 0. ) {
213 dEdx = emcalc.GetDEDX(kinEnergyPre, part, mate);
214 }
215 G4double elost = stepLengths[ii] * dEdx;
216#ifdef VERBOSE_ENERSPLIT
217 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes"<< ii << " RN: iter" << iiter << " energy lost " << elost
218 << " energy at interaction " << kinEnergyPre
219 << " = stepLength " << stepLengths[ii]
220 << " * dEdx " << dEdx << G4endl;
221#endif
222 kinEnergyPre -= elost;
223 theEnergies[ii] = elost;
224 totalELost += elost;
225 }
226
227 }
228
229 //correct energies so that they reproduce the real step energy lost
230 G4double enerRatio = (edep/totalELost);
231
232#ifdef VERBOSE_ENERSPLIT
233 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes"<< ii << " RN: iter" << iiter << " energy ratio " << enerRatio << G4endl;
234#endif
235
236#ifdef VERBOSE_ENERSPLIT
237 G4double elostTot = 0.;
238#endif
239 for( ii = 0; ii < theEnergies.size(); ii++ ){
240 theEnergies[ii] *= enerRatio;
241#ifdef VERBOSE_ENERSPLIT
242 elostTot += theEnergies[ii];
243 if(verbose) G4cout << "G4EnergySplitter::SplitEnergyInVolumes "<< ii << " RN: iter" << iiter << " corrected energy lost " << theEnergies[ii]
244 << " orig elost " << theEnergies[ii]/enerRatio
245 << " energy before interaction " << kinEnergyPreOrig-elostTot+theEnergies[ii]
246 << " energy after interaction " << kinEnergyPreOrig-elostTot
247 << G4endl;
248#endif
249 }
250 }
251
252 }
253
254 return theEnergies.size();
255}
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *, const G4Region *r=0)
G4double TrueStepLength(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *part)
const G4String & GetName() const
Definition: G4Material.hh:177
G4double GetPDGCharge() const
G4Material * GetMaterial(size_t nx, size_t ny, size_t nz) const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
G4ParticleDefinition * GetDefinition() const
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52

Referenced by G4ScoreSplittingProcess::PostStepDoIt().


The documentation for this class was generated from the following files: