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
G4DensityEffectCalculator.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 * Implements calculation of the Fermi density effect as per the method
29 * described in:
30 *
31 * R. M. Sternheimer, M. J. Berger, and S. M. Seltzer. Density
32 * effect for the ionization loss of charged particles in various sub-
33 * stances. Atom. Data Nucl. Data Tabl., 30:261, 1984.
34 *
35 * Which (among other Sternheimer references) builds on:
36 *
37 * R. M. Sternheimer. The density effect for ionization loss in
38 * materials. Phys. Rev., 88:851­859, 1952.
39 *
40 * The returned values of delta are directly from the Sternheimer calculation,
41 * and not Sternheimer's popular three-part approximate parameterization
42 * introduced in the same paper.
43 *
44 * Author: Matthew Strait <straitm@umn.edu> 2019
45 */
46
48#include "G4AtomicShells.hh"
49#include "G4NistManager.hh"
50#include "G4Pow.hh"
51
52static G4Pow * gpow = G4Pow::GetInstance();
53
54const G4int maxWarnings = 20;
55
57 : fMaterial(mat), fVerbose(0), fWarnings(0), nlev(n)
58{
59 fVerbose = std::max(fVerbose, G4NistManager::Instance()->GetVerbose());
60
61 sternf = new G4double [nlev];
62 levE = new G4double [nlev];
63 sternl = new G4double [nlev];
64 sternEbar = new G4double [nlev];
65 for(G4int i=0; i<nlev; ++i) {
66 sternf[i] = 0.0;
67 levE[i] = 0.0;
68 sternl[i] = 0.0;
69 sternEbar[i] = 0.0;
70 }
71
72 fConductivity = sternx = 0.0;
73 G4bool conductor = (fMaterial->GetFreeElectronDensity() > 0.0);
74
75 G4int sh = 0;
76 G4double sum = 0.;
77 const G4double tot = fMaterial->GetTotNbOfAtomsPerVolume();
78 for(size_t j = 0; j < fMaterial->GetNumberOfElements(); ++j) {
79 // The last subshell is considered to contain the conduction
80 // electrons. Sternheimer 1984 says "the lowest chemical valance of
81 // the element" is used to set the number of conduction electrons.
82 // I'm not sure if that means the highest subshell or the whole
83 // shell, but in any case, he also says that the choice is arbitrary
84 // and offers a possible alternative. This is one of the sources of
85 // uncertainty in the model.
86 const G4double frac = fMaterial->GetVecNbOfAtomsPerVolume()[j]/tot;
87 const G4int Z = fMaterial->GetElement((G4int)j)->GetZasInt();
89 for(G4int i = 0; i < nshell; ++i) {
90 // For conductors, put *all* top shell electrons into the conduction
91 // band, regardless of element.
93 if(i < nshell-1 || !conductor) {
94 sternf[sh] += xx;
95 } else {
96 fConductivity += xx;
97 }
98 levE[sh] = G4AtomicShells::GetBindingEnergy(Z, i)/CLHEP::eV;
99 ++sh;
100 }
101 }
102 for(G4int i=0; i<nlev; ++i) {
103 sum += sternf[i];
104 }
105 sum += fConductivity;
106
107 const G4double invsum = (sum > 0.0) ? 1./sum : 0.0;
108 for(G4int i=0; i<nlev; ++i) {
109 sternf[i] *= invsum;
110 }
111 fConductivity *= invsum;
112 plasmaE = fMaterial->GetIonisation()->GetPlasmaEnergy()/CLHEP::eV;
113 meanexcite = fMaterial->GetIonisation()->GetMeanExcitationEnergy()/CLHEP::eV;
114}
115
117{
118 delete [] sternf;
119 delete [] levE;
120 delete [] sternl;
121 delete [] sternEbar;
122}
123
125{
126 if(fVerbose > 1) {
127 G4cout << "G4DensityEffectCalculator::ComputeDensityCorrection for "
128 << fMaterial->GetName() << ", x= " << x << G4endl;
129 }
130 const G4double approx = fMaterial->GetIonisation()->GetDensityCorrection(x);
131 const G4double exact = FermiDeltaCalculation(x);
132
133 if(fVerbose > 1) {
134 G4cout << " Delta: computed= " << exact
135 << ", parametrized= " << approx << G4endl;
136 }
137 if(approx >= 0. && exact < 0.) {
138 if(fVerbose > 0) {
139 ++fWarnings;
140 if(fWarnings < maxWarnings) {
142 ed << "Sternheimer fit failed for " << fMaterial->GetName()
143 << ", x = " << x << ": Delta exact= "
144 << exact << ", approx= " << approx;
145 G4Exception("G4DensityEffectCalculator::DensityCorrection", "mat008",
146 JustWarning, ed);
147 }
148 }
149 return approx;
150 }
151 // Fall back to approx if exact and approx are very different, under the
152 // assumption that this means the exact calculation has gone haywire
153 // somehow, with the exception of the case where approx is negative. I
154 // have seen this clearly-wrong result occur for substances with extremely
155 // low density (1e-25 g/cc).
156 if(approx >= 0. && std::abs(exact - approx) > 1.) {
157 if(fVerbose > 0) {
158 ++fWarnings;
159 if(fWarnings < maxWarnings) {
161 ed << "Sternheimer exact= " << exact << " and approx= "
162 << approx << " are too different for "
163 << fMaterial->GetName() << ", x = " << x;
164 G4Exception("G4DensityEffectCalculator::DensityCorrection", "mat008",
165 JustWarning, ed);
166 }
167 }
168 return approx;
169 }
170 return exact;
171}
172
173G4double G4DensityEffectCalculator::FermiDeltaCalculation(G4double x)
174{
175 // Above beta*gamma of 10^10, the exact treatment is within machine
176 // precision of the limiting case, for ordinary solids, at least. The
177 // convergence goes up as the density goes down, but even in a pretty
178 // hard vacuum it converges by 10^20. Also, it's hard to imagine how
179 // this energy is relevant (x = 20 -> 10^19 GeV for muons). So this
180 // is mostly not here for physical reasons, but rather to avoid ugly
181 // discontinuities in the return value.
182 if(x > 20.) { return -1.; }
183
184 sternx = x;
185 G4double sternrho = Newton(1.5, true);
186
187 // Negative values, and values much larger than unity are non-physical.
188 // Values between zero and one are also suspect, but not as clearly wrong.
189 if(sternrho <= 0. || sternrho > 100.) {
190 if(fVerbose > 0) {
191 ++fWarnings;
192 if(fWarnings < maxWarnings) {
194 ed << "Sternheimer computation failed for " << fMaterial->GetName()
195 << ", x = " << x << ":\n"
196 << "Could not solve for Sternheimer rho. Probably you have a \n"
197 << "mean ionization energy which is incompatible with your\n"
198 << "distribution of energy levels, or an unusually dense material.\n"
199 << "Number of levels: " << nlev
200 << " Mean ionization energy(eV): " << meanexcite
201 << " Plasma energy(eV): " << plasmaE << "\n";
202 for(G4int i = 0; i < nlev; ++i) {
203 ed << "Level " << i << ": strength " << sternf[i]
204 << ": energy(eV)= " << levE[i] << "\n";
205 }
206 G4Exception("G4DensityEffectCalculator::SetupFermiDeltaCalc", "mat008",
207 JustWarning, ed);
208 }
209 }
210 return -1.;
211 }
212
213 // Calculate the Sternheimer adjusted energy levels and parameters l_i given
214 // the Sternheimer parameter rho.
215 for(G4int i=0; i<nlev; ++i) {
216 sternEbar[i] = levE[i] * (sternrho/plasmaE);
217 sternl[i] = std::sqrt(gpow->powN(sternEbar[i], 2) + (2./3.)*sternf[i]);
218 }
219 // The derivative of the function we are solving for is strictly
220 // negative for positive (physical) values, so if the value at
221 // zero is less than zero, it has no solution, and there is no
222 // density effect in the Sternheimer "exact" treatment (which is
223 // still an approximation).
224 //
225 // For conductors, this test is not needed, because Ell(L) contains
226 // the term fConductivity/(L*L), so the value at L=0 is always
227 // positive infinity. In the code we don't return inf, though, but
228 // rather set that term to zero, which means that if this test were
229 // used, it would give the wrong result for some materials.
230 if(fConductivity == 0 && Ell(0) <= 0)
231 {
232 return 0;
233 }
234
235 // Attempt to find the root from 40 starting points evenly distributed
236 // in log space. Trying a single starting point is not sufficient for
237 // convergence in most cases.
238 for(G4int startLi = -10; startLi < 30; ++startLi){
239 const G4double sternL = Newton(gpow->powN(2, startLi), false);
240 if(sternL != -1.) {
241 return DeltaOnceSolved(sternL);
242 }
243 }
244 return -1.; // Signal the caller to use the Sternheimer approximation,
245 // because we have been unable to solve the exact form.
246}
247
248/* Newton's method for finding roots. Adapted from G4PolynominalSolver, but
249 * without the assumption that the input is a polynomial. Also, here we
250 * always expect the roots to be positive, so return -1 as an error value. */
251G4double G4DensityEffectCalculator::Newton(G4double start, G4bool first)
252{
253 const G4int maxIter = 100;
254 G4int nbad = 0, ngood = 0;
255
256 G4double lambda(start), value(0.), dvalue(0.);
257
258 if(fVerbose > 2) {
259 G4cout << "G4DensityEffectCalculator::Newton: strat= " << start
260 << " type: " << first << G4endl;
261 }
262 while(true) {
263 if(first) {
264 value = FRho(lambda);
265 dvalue = DFRho(lambda);
266 } else {
267 value = Ell(lambda);
268 dvalue = DEll(lambda);
269 }
270 if(dvalue == 0.0) { break; }
271 const G4double del = value/dvalue;
272 lambda -= del;
273
274 const G4double eps = std::abs(del/lambda);
275 if(eps <= 1.e-12) {
276 ++ngood;
277 if(ngood == 2) {
278 if(fVerbose > 2) {
279 G4cout << " Converged with result= " << lambda << G4endl;
280 }
281 return lambda;
282 }
283 } else {
284 ++nbad;
285 }
286 if(nbad > maxIter || std::isnan(value) || std::isinf(value)) { break; }
287 }
288 if(fVerbose > 2) {
289 G4cout << " Failed to converge last value= " << value
290 << " dvalue= " << dvalue << " lambda= " << lambda << G4endl;
291 }
292 return -1.;
293}
294
295/* Return the derivative of the equation used
296 * to solve for the Sternheimer parameter rho. */
297G4double G4DensityEffectCalculator::DFRho(G4double rho)
298{
299 G4double ans = 0.0;
300 for(G4int i = 0; i < nlev; ++i) {
301 if(sternf[i] > 0.) {
302 ans += sternf[i] * gpow->powN(levE[i], 2) * rho /
303 (gpow->powN(levE[i] * rho, 2)
304 + 2./3. * sternf[i] * gpow->powN(plasmaE, 2));
305 }
306 }
307 return ans;
308}
309
310/* Return the functional value for the equation used
311 * to solve for the Sternheimer parameter rho. */
312G4double G4DensityEffectCalculator::FRho(G4double rho)
313{
314 G4double ans = 0.0;
315 for(G4int i = 0; i<nlev; ++i) {
316 if(sternf[i] > 0.) {
317 ans += sternf[i] * G4Log(gpow->powN(levE[i]*rho, 2) +
318 2./3. * sternf[i]*gpow->powN(plasmaE, 2));
319 }
320 }
321 ans *= 0.5; // pulled out of loop for efficiency
322
323 if(fConductivity > 0.) {
324 ans += fConductivity * G4Log(plasmaE * std::sqrt(fConductivity));
325 }
326 ans -= G4Log(meanexcite);
327 return ans;
328}
329
330/* Return the derivative for the equation used to
331 * solve for the Sternheimer parameter l, called 'L' here. */
332G4double G4DensityEffectCalculator::DEll(G4double L)
333{
334 G4double ans = 0.;
335 for(G4int i=0; i<nlev; ++i) {
336 if(sternf[i] > 0 && (sternEbar[i] > 0. || L != 0.)) {
337 const G4double y = gpow->powN(sternEbar[i], 2);
338 ans += sternf[i]/gpow->powN(y + L*L, 2);
339 }
340 }
341 ans += fConductivity/gpow->powN(L*L, 2);
342 ans *= (-2*L); // pulled out of the loop for efficiency
343 return ans;
344}
345
346/* Return the functional value for the equation used to
347 * solve for the Sternheimer parameter l, called 'L' here. */
348G4double G4DensityEffectCalculator::Ell(G4double L)
349{
350 G4double ans = 0.;
351 for(G4int i=0; i<nlev; ++i) {
352 if(sternf[i] > 0. && (sternEbar[i] > 0. || L != 0.)) {
353 ans += sternf[i]/(gpow->powN(sternEbar[i], 2) + L*L);
354 }
355 }
356 if(fConductivity > 0. && L != 0.) {
357 ans += fConductivity/(L*L);
358 }
359 ans -= gpow->powZ(10, -2 * sternx);
360 return ans;
361}
362
363/**
364 * Given the Sternheimer parameter l (called 'sternL' here), and that
365 * the l_i and adjusted energies have been found with SetupFermiDeltaCalc(),
366 * return the value of delta. Helper function for DoFermiDeltaCalc().
367 */
368G4double G4DensityEffectCalculator::DeltaOnceSolved(G4double sternL)
369{
370 G4double ans = 0.;
371 for(G4int i=0; i<nlev; ++i) {
372 if(sternf[i] > 0.) {
373 ans += sternf[i] * G4Log((gpow->powN(sternl[i], 2)
374 + gpow->powN(sternL, 2))/gpow->powN(sternl[i], 2));
375 }
376 }
377 // sternl for the conduction electrons is sqrt(fConductivity), with
378 // no factor of 2./3 as with the other levels.
379 if(fConductivity > 0) {
380 ans += fConductivity * G4Log((fConductivity
381 + gpow->powN(sternL, 2))/fConductivity);
382 }
383 ans -= gpow->powN(sternL, 2)/(1 + gpow->powZ(10, 2 * sternx));
384 return ans;
385}
const G4int maxWarnings
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
static G4int GetNumberOfShells(G4int Z)
G4DensityEffectCalculator(const G4Material *, G4int)
G4double ComputeDensityCorrection(G4double x)
G4double GetMeanExcitationEnergy() const
G4double GetDensityCorrection(G4double x) const
G4double GetPlasmaEnergy() const
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:197
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
G4double GetFreeElectronDensity() const
Definition: G4Material.hh:174
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
const G4String & GetName() const
Definition: G4Material.hh:172
static G4NistManager * Instance()
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:225
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:162