Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // 26 26 27 /* 27 /* 28 * Implements calculation of the Fermi density 28 * Implements calculation of the Fermi density effect as per the method 29 * described in: 29 * described in: 30 * 30 * 31 * R. M. Sternheimer, M. J. Berger, and S. M 31 * R. M. Sternheimer, M. J. Berger, and S. M. Seltzer. Density 32 * effect for the ionization loss of charged 32 * effect for the ionization loss of charged particles in various sub- 33 * stances. Atom. Data Nucl. Data Tabl., 30: 33 * stances. Atom. Data Nucl. Data Tabl., 30:261, 1984. 34 * 34 * 35 * Which (among other Sternheimer references) 35 * Which (among other Sternheimer references) builds on: 36 * 36 * 37 * R. M. Sternheimer. The density effect for 37 * R. M. Sternheimer. The density effect for ionization loss in 38 * materials. Phys. Rev., 88:851859, 1952. 38 * materials. Phys. Rev., 88:851859, 1952. 39 * 39 * 40 * The returned values of delta are directly f 40 * The returned values of delta are directly from the Sternheimer calculation, 41 * and not Sternheimer's popular three-part ap 41 * and not Sternheimer's popular three-part approximate parameterization 42 * introduced in the same paper. 42 * introduced in the same paper. 43 * 43 * 44 * Author: Matthew Strait <straitm@umn.edu> 20 44 * Author: Matthew Strait <straitm@umn.edu> 2019 45 */ 45 */ 46 46 47 #include "G4DensityEffectCalculator.hh" 47 #include "G4DensityEffectCalculator.hh" 48 << 49 #include "G4AtomicShells.hh" 48 #include "G4AtomicShells.hh" 50 #include "G4NistManager.hh" 49 #include "G4NistManager.hh" 51 #include "G4Pow.hh" 50 #include "G4Pow.hh" 52 51 53 static G4Pow* gpow = G4Pow::GetInstance(); << 52 static G4Pow * gpow = G4Pow::GetInstance(); 54 53 55 const G4int maxWarnings = 20; 54 const G4int maxWarnings = 20; 56 55 57 G4DensityEffectCalculator::G4DensityEffectCalc 56 G4DensityEffectCalculator::G4DensityEffectCalculator(const G4Material* mat, G4int n) 58 : fMaterial(mat), nlev(n) << 57 : fMaterial(mat), fVerbose(0), fWarnings(0), nlev(n) 59 { 58 { 60 fVerbose = std::max(fVerbose, G4NistManager: 59 fVerbose = std::max(fVerbose, G4NistManager::Instance()->GetVerbose()); 61 60 62 sternf = new G4double[nlev]; << 61 sternf = new G4double [nlev]; 63 levE = new G4double[nlev]; << 62 levE = new G4double [nlev]; 64 sternl = new G4double[nlev]; << 63 sternl = new G4double [nlev]; 65 sternEbar = new G4double[nlev]; << 64 sternEbar = new G4double [nlev]; 66 for (G4int i = 0; i < nlev; ++i) { << 65 for(G4int i=0; i<nlev; ++i) { 67 sternf[i] = 0.0; << 66 sternf[i] = 0.0; 68 levE[i] = 0.0; << 67 levE[i] = 0.0; 69 sternl[i] = 0.0; << 68 sternl[i] = 0.0; 70 sternEbar[i] = 0.0; 69 sternEbar[i] = 0.0; 71 } 70 } 72 71 73 fConductivity = sternx = 0.0; 72 fConductivity = sternx = 0.0; 74 G4bool conductor = (fMaterial->GetFreeElectr 73 G4bool conductor = (fMaterial->GetFreeElectronDensity() > 0.0); 75 74 76 G4int sh = 0; 75 G4int sh = 0; 77 G4double sum = 0.; 76 G4double sum = 0.; 78 const G4double tot = fMaterial->GetTotNbOfAt 77 const G4double tot = fMaterial->GetTotNbOfAtomsPerVolume(); 79 for (size_t j = 0; j < fMaterial->GetNumberO << 78 for(size_t j = 0; j < fMaterial->GetNumberOfElements(); ++j) { 80 // The last subshell is considered to cont 79 // The last subshell is considered to contain the conduction 81 // electrons. Sternheimer 1984 says "the l 80 // electrons. Sternheimer 1984 says "the lowest chemical valance of 82 // the element" is used to set the number 81 // the element" is used to set the number of conduction electrons. 83 // I'm not sure if that means the highest 82 // I'm not sure if that means the highest subshell or the whole 84 // shell, but in any case, he also says th 83 // shell, but in any case, he also says that the choice is arbitrary 85 // and offers a possible alternative. This 84 // and offers a possible alternative. This is one of the sources of 86 // uncertainty in the model. 85 // uncertainty in the model. 87 const G4double frac = fMaterial->GetVecNbO << 86 const G4double frac = fMaterial->GetVecNbOfAtomsPerVolume()[j]/tot; 88 const G4int Z = fMaterial->GetElement((G4i << 87 const G4int Z = fMaterial->GetElement(j)->GetZasInt(); 89 const G4int nshell = G4AtomicShells::GetNu 88 const G4int nshell = G4AtomicShells::GetNumberOfShells(Z); 90 for (G4int i = 0; i < nshell; ++i) { << 89 for(G4int i = 0; i < nshell; ++i) { 91 // For conductors, put *all* top shell e 90 // For conductors, put *all* top shell electrons into the conduction 92 // band, regardless of element. 91 // band, regardless of element. 93 const G4double xx = frac * G4AtomicShell << 92 const G4double xx = frac*G4AtomicShells::GetNumberOfElectrons(Z, i); 94 if (i < nshell - 1 || ! conductor) { << 93 if(i < nshell-1 || !conductor) { 95 sternf[sh] += xx; << 94 sternf[sh] += xx; 96 } << 95 } else { 97 else { << 98 fConductivity += xx; 96 fConductivity += xx; 99 } 97 } 100 levE[sh] = G4AtomicShells::GetBindingEne << 98 levE[sh] = G4AtomicShells::GetBindingEnergy(Z, i)/CLHEP::eV; 101 ++sh; 99 ++sh; 102 } 100 } 103 } 101 } 104 for (G4int i = 0; i < nlev; ++i) { << 102 for(G4int i=0; i<nlev; ++i) { 105 sum += sternf[i]; 103 sum += sternf[i]; 106 } 104 } 107 sum += fConductivity; << 105 sum = (sum > 0.0) ? 1./sum : 0.0; 108 << 106 for(G4int i=0; i<nlev; ++i) { 109 const G4double invsum = (sum > 0.0) ? 1. / s << 107 sternf[i] *= sum; 110 for (G4int i = 0; i < nlev; ++i) { << 108 } 111 sternf[i] *= invsum; << 109 plasmaE = fMaterial->GetIonisation()->GetPlasmaEnergy()/CLHEP::eV; 112 } << 110 meanexcite = fMaterial->GetIonisation()->GetMeanExcitationEnergy()/CLHEP::eV; 113 fConductivity *= invsum; << 114 plasmaE = fMaterial->GetIonisation()->GetPla << 115 meanexcite = fMaterial->GetIonisation()->Get << 116 } 111 } 117 112 118 G4DensityEffectCalculator::~G4DensityEffectCal 113 G4DensityEffectCalculator::~G4DensityEffectCalculator() 119 { 114 { 120 delete[] sternf; << 115 delete [] sternf; 121 delete[] levE; << 116 delete [] levE; 122 delete[] sternl; << 117 delete [] sternl; 123 delete[] sternEbar; << 118 delete [] sternEbar; 124 } 119 } 125 120 126 G4double G4DensityEffectCalculator::ComputeDen 121 G4double G4DensityEffectCalculator::ComputeDensityCorrection(G4double x) 127 { 122 { 128 if (fVerbose > 1) { << 123 if(fVerbose > 1) { 129 G4cout << "G4DensityEffectCalculator::Comp << 124 G4cout << "G4DensityEffectCalculator::ComputeDensityCorrection for " 130 << ", x= " << x << G4endl; << 125 << fMaterial->GetName() << ", x= " << x << G4endl; 131 } 126 } 132 const G4double approx = fMaterial->GetIonisa 127 const G4double approx = fMaterial->GetIonisation()->GetDensityCorrection(x); 133 const G4double exact = FermiDeltaCalculation << 128 const G4double exact = FermiDeltaCalculation(x); 134 129 135 if (fVerbose > 1) { << 130 if(fVerbose > 1) { 136 G4cout << " Delta: computed= " << exact << 131 G4cout << " Delta: computed= " << exact >> 132 << ", parametrized= " << approx << G4endl; 137 } 133 } 138 if (approx >= 0. && exact < 0.) { << 134 if(approx > 0. && exact < 0.) { 139 if (fVerbose > 0) { << 135 if(fVerbose > 0) { 140 ++fWarnings; 136 ++fWarnings; 141 if (fWarnings < maxWarnings) { << 137 if(fWarnings < maxWarnings) { 142 G4ExceptionDescription ed; << 138 G4ExceptionDescription ed; 143 ed << "Sternheimer fit failed for " << << 139 ed << "Sternheimer fit failed for " << fMaterial->GetName() 144 << ": Delta exact= " << exact << ", << 140 << ", x = " << x << ": Delta exact= " 145 G4Exception("G4DensityEffectCalculator << 141 << exact << ", approx= " << approx; >> 142 G4Exception("G4DensityEffectCalculator::DensityCorrection", "mat008", >> 143 JustWarning, ed); 146 } 144 } 147 } 145 } 148 return approx; 146 return approx; 149 } 147 } 150 // Fall back to approx if exact and approx a 148 // Fall back to approx if exact and approx are very different, under the 151 // assumption that this means the exact calc 149 // assumption that this means the exact calculation has gone haywire 152 // somehow, with the exception of the case w 150 // somehow, with the exception of the case where approx is negative. I 153 // have seen this clearly-wrong result occur 151 // have seen this clearly-wrong result occur for substances with extremely 154 // low density (1e-25 g/cc). 152 // low density (1e-25 g/cc). 155 if (approx >= 0. && std::abs(exact - approx) << 153 if(approx >= 0. && std::abs(exact - approx) > 1.) { 156 if (fVerbose > 0) { << 154 if(fVerbose > 0) { 157 ++fWarnings; 155 ++fWarnings; 158 if (fWarnings < maxWarnings) { << 156 if(fWarnings < maxWarnings) { 159 G4ExceptionDescription ed; << 157 G4ExceptionDescription ed; 160 ed << "Sternheimer exact= " << exact < << 158 ed << "Sternheimer exact= " << exact << " and approx= " 161 << " are too different for " << fMa << 159 << approx << " are too different for " 162 G4Exception("G4DensityEffectCalculator << 160 << fMaterial->GetName() << ", x = " << x; >> 161 G4Exception("G4DensityEffectCalculator::DensityCorrection", "mat008", >> 162 JustWarning, ed); 163 } 163 } 164 } 164 } 165 return approx; 165 return approx; 166 } 166 } 167 return exact; 167 return exact; 168 } 168 } 169 169 170 G4double G4DensityEffectCalculator::FermiDelta 170 G4double G4DensityEffectCalculator::FermiDeltaCalculation(G4double x) 171 { 171 { 172 // Above beta*gamma of 10^10, the exact trea 172 // Above beta*gamma of 10^10, the exact treatment is within machine 173 // precision of the limiting case, for ordin 173 // precision of the limiting case, for ordinary solids, at least. The 174 // convergence goes up as the density goes d 174 // convergence goes up as the density goes down, but even in a pretty 175 // hard vacuum it converges by 10^20. Also, 175 // hard vacuum it converges by 10^20. Also, it's hard to imagine how 176 // this energy is relevant (x = 20 -> 10^19 176 // this energy is relevant (x = 20 -> 10^19 GeV for muons). So this 177 // is mostly not here for physical reasons, 177 // is mostly not here for physical reasons, but rather to avoid ugly 178 // discontinuities in the return value. 178 // discontinuities in the return value. 179 if (x > 20.) { << 179 if(x > 20.) { return -1.; } 180 return -1.; << 181 } << 182 180 183 sternx = x; 181 sternx = x; 184 G4double sternrho = Newton(1.5, true); 182 G4double sternrho = Newton(1.5, true); 185 183 186 // Negative values, and values much larger t 184 // Negative values, and values much larger than unity are non-physical. 187 // Values between zero and one are also susp 185 // Values between zero and one are also suspect, but not as clearly wrong. 188 if (sternrho <= 0. || sternrho > 100.) { << 186 if(sternrho <= 0. || sternrho > 100.) { 189 if (fVerbose > 0) { << 187 if(fVerbose > 0) { 190 ++fWarnings; 188 ++fWarnings; 191 if (fWarnings < maxWarnings) { << 189 if(fWarnings < maxWarnings) { 192 G4ExceptionDescription ed; << 190 G4ExceptionDescription ed; 193 ed << "Sternheimer computation failed << 191 ed << "Sternheimer computation failed for " << fMaterial->GetName() 194 << ":\n" << 192 << ", x = " << x << ":\n" 195 << "Could not solve for Sternheimer << 193 << "Could not solve for Sternheimer rho. Probably you have a \n" 196 << "mean ionization energy which is << 194 << "mean ionization energy which is incompatible with your\n" 197 << "distribution of energy levels, << 195 << "distribution of energy levels, or an unusually dense material.\n" 198 << "Number of levels: " << nlev << << 196 << "Number of levels: " << nlev 199 << " Plasma energy(eV): " << plasma << 197 << " Mean ionization energy(eV): " << meanexcite 200 for (G4int i = 0; i < nlev; ++i) { << 198 << " Plasma energy(eV): " << plasmaE << "\n"; 201 ed << "Level " << i << ": strength " << 199 for(G4int i = 0; i < nlev; ++i) { 202 } << 200 ed << "Level " << i << ": strength " << sternf[i] 203 G4Exception("G4DensityEffectCalculator << 201 << ": energy(eV)= " << levE[i] << "\n"; >> 202 } >> 203 G4Exception("G4DensityEffectCalculator::SetupFermiDeltaCalc", "mat008", >> 204 JustWarning, ed); 204 } 205 } 205 } 206 } 206 return -1.; 207 return -1.; 207 } 208 } 208 209 209 // Calculate the Sternheimer adjusted energy 210 // Calculate the Sternheimer adjusted energy levels and parameters l_i given 210 // the Sternheimer parameter rho. 211 // the Sternheimer parameter rho. 211 for (G4int i = 0; i < nlev; ++i) { << 212 sternrho /= plasmaE; 212 sternEbar[i] = levE[i] * (sternrho / plasm << 213 for(G4int i=0; i<nlev; ++i) { 213 sternl[i] = std::sqrt(gpow->powN(sternEbar << 214 sternEbar[i] = levE[i] * sternrho; 214 } << 215 sternl[i] = std::sqrt(gpow->powN(sternEbar[i], 2) + 2./3.*sternf[i]); 215 // The derivative of the function we are sol << 216 } 216 // negative for positive (physical) values, << 217 217 // zero is less than zero, it has no solutio << 218 // Make imphirical initial guess 218 // density effect in the Sternheimer "exact" << 219 const G4double sternL = Newton(sternrho, false); 219 // still an approximation). << 220 if(sternL > -1.) { 220 // << 221 return DeltaOnceSolved(sternL); 221 // For conductors, this test is not needed, << 222 // the term fConductivity/(L*L), so the valu << 223 // positive infinity. In the code we don't r << 224 // rather set that term to zero, which means << 225 // used, it would give the wrong result for << 226 if (fConductivity == 0 && Ell(0) <= 0) { << 227 return 0; << 228 } << 229 << 230 // Attempt to find the root from 40 starting << 231 // in log space. Trying a single starting p << 232 // convergence in most cases. << 233 for (G4int startLi = -10; startLi < 30; ++st << 234 const G4double sternL = Newton(gpow->powN( << 235 if (sternL != -1.) { << 236 return DeltaOnceSolved(sternL); << 237 } << 238 } 222 } 239 return -1.; // Signal the caller to use the << 223 240 // because we have been unable << 224 return -1.; // Signal the caller to use the Sternheimer approximation, >> 225 // because we have been unable to solve the exact form. 241 } 226 } 242 227 243 /* Newton's method for finding roots. Adapted 228 /* Newton's method for finding roots. Adapted from G4PolynominalSolver, but 244 * without the assumption that the input is a 229 * without the assumption that the input is a polynomial. Also, here we 245 * always expect the roots to be positive, so 230 * always expect the roots to be positive, so return -1 as an error value. */ 246 G4double G4DensityEffectCalculator::Newton(G4d 231 G4double G4DensityEffectCalculator::Newton(G4double start, G4bool first) 247 { 232 { 248 const G4int maxIter = 100; 233 const G4int maxIter = 100; 249 G4int nbad = 0, ngood = 0; 234 G4int nbad = 0, ngood = 0; 250 235 251 G4double lambda(start), value(0.), dvalue(0. 236 G4double lambda(start), value(0.), dvalue(0.); 252 237 253 if (fVerbose > 2) { << 238 if(fVerbose > 2) { 254 G4cout << "G4DensityEffectCalculator::Newt << 239 G4cout << "G4DensityEffectCalculator::Newton: strat= " << start >> 240 << " type: " << first << G4endl; 255 } 241 } 256 while (true) { << 242 while(true) { 257 if (first) { << 243 if(first) { 258 value = FRho(lambda); 244 value = FRho(lambda); 259 dvalue = DFRho(lambda); 245 dvalue = DFRho(lambda); 260 } << 246 } else { 261 else { << 262 value = Ell(lambda); 247 value = Ell(lambda); 263 dvalue = DEll(lambda); 248 dvalue = DEll(lambda); 264 } 249 } 265 if (dvalue == 0.0) { << 250 if(dvalue == 0.0) { break; } 266 break; << 251 const G4double del = value/dvalue; 267 } << 268 const G4double del = value / dvalue; << 269 lambda -= del; 252 lambda -= del; 270 253 271 const G4double eps = std::abs(del / lambda << 254 const G4double eps = std::abs(del); 272 if (eps <= 1.e-12) { << 255 if(eps <= 1.e-12) { 273 ++ngood; 256 ++ngood; 274 if (ngood == 2) { << 257 if(ngood == 2) { 275 if (fVerbose > 2) { << 258 if(fVerbose > 2) { 276 G4cout << " Converged with result= << 259 G4cout << " Converged with result= " << lambda << G4endl; 277 } << 260 } 278 return lambda; << 261 return lambda; 279 } 262 } 280 } << 263 } else { 281 else { << 282 ++nbad; 264 ++nbad; 283 } 265 } 284 if (nbad > maxIter || std::isnan(value) || << 266 if(nbad > maxIter || eps > 1.) { break; } 285 break; << 286 } << 287 } 267 } 288 if (fVerbose > 2) { << 268 if(fVerbose > 2) { 289 G4cout << " Failed to converge last value << 269 G4cout << " Failed to converge last value= " << value 290 << " lambda= " << lambda << G4endl; << 270 << " dvalue= " << dvalue << " lambda= " << lambda << G4endl; 291 } 271 } 292 return -1.; 272 return -1.; 293 } 273 } 294 274 295 /* Return the derivative of the equation used 275 /* Return the derivative of the equation used 296 * to solve for the Sternheimer parameter rho. 276 * to solve for the Sternheimer parameter rho. */ 297 G4double G4DensityEffectCalculator::DFRho(G4do 277 G4double G4DensityEffectCalculator::DFRho(G4double rho) 298 { 278 { 299 G4double ans = 0.0; 279 G4double ans = 0.0; 300 for (G4int i = 0; i < nlev; ++i) { << 280 for(G4int i = 0; i < nlev; ++i) { 301 if (sternf[i] > 0.) { << 281 if(sternf[i] > 0.) { 302 ans += sternf[i] * gpow->powN(levE[i], 2 282 ans += sternf[i] * gpow->powN(levE[i], 2) * rho / 303 (gpow->powN(levE[i] * rho, 2) + 2 << 283 (gpow->powN(levE[i] * rho, 2) >> 284 + 2./3. * sternf[i] * gpow->powN(plasmaE, 2)); 304 } 285 } 305 } 286 } 306 return ans; 287 return ans; 307 } 288 } 308 289 309 /* Return the functional value for the equatio 290 /* Return the functional value for the equation used 310 * to solve for the Sternheimer parameter rho. 291 * to solve for the Sternheimer parameter rho. */ 311 G4double G4DensityEffectCalculator::FRho(G4dou 292 G4double G4DensityEffectCalculator::FRho(G4double rho) 312 { 293 { 313 G4double ans = 0.0; 294 G4double ans = 0.0; 314 for (G4int i = 0; i < nlev; ++i) { << 295 for(G4int i = 0; i<nlev; ++i) { 315 if (sternf[i] > 0.) { << 296 if(sternf[i] > 0.) { 316 ans += sternf[i] * << 297 ans += sternf[i] * G4Log(gpow->powN(levE[i]*rho, 2) + 317 G4Log(gpow->powN(levE[i] * rho, 2 << 298 2./3. * sternf[i]*gpow->powN(plasmaE, 2)); 318 } << 299 } 319 } 300 } 320 ans *= 0.5; // pulled out of loop for effic << 301 ans *= 0.5; // pulled out of loop for efficiency 321 302 322 if (fConductivity > 0.) { << 303 if(fConductivity > 0.) { 323 ans += fConductivity * G4Log(plasmaE * std 304 ans += fConductivity * G4Log(plasmaE * std::sqrt(fConductivity)); 324 } 305 } 325 ans -= G4Log(meanexcite); 306 ans -= G4Log(meanexcite); 326 return ans; 307 return ans; 327 } 308 } 328 309 329 /* Return the derivative for the equation used 310 /* Return the derivative for the equation used to 330 * solve for the Sternheimer parameter l, call 311 * solve for the Sternheimer parameter l, called 'L' here. */ 331 G4double G4DensityEffectCalculator::DEll(G4dou 312 G4double G4DensityEffectCalculator::DEll(G4double L) 332 { 313 { 333 G4double ans = 0.; 314 G4double ans = 0.; 334 for (G4int i = 0; i < nlev; ++i) { << 315 for(G4int i=0; i<nlev; ++i) { 335 if (sternf[i] > 0 && (sternEbar[i] > 0. || << 316 if(sternf[i] > 0 && (sternEbar[i] > 0. || L != 0.)) { 336 const G4double y = gpow->powN(sternEbar[ 317 const G4double y = gpow->powN(sternEbar[i], 2); 337 ans += sternf[i] / gpow->powN(y + L * L, << 318 ans += sternf[i]/gpow->powN(y + L*L, 2); 338 } 319 } 339 } 320 } 340 ans += fConductivity / gpow->powN(L * L, 2); << 321 ans *= (-2*L); // pulled out of the loop for efficiency 341 ans *= (-2 * L); // pulled out of the loop << 342 return ans; 322 return ans; 343 } 323 } 344 324 345 /* Return the functional value for the equatio 325 /* Return the functional value for the equation used to 346 * solve for the Sternheimer parameter l, call 326 * solve for the Sternheimer parameter l, called 'L' here. */ 347 G4double G4DensityEffectCalculator::Ell(G4doub 327 G4double G4DensityEffectCalculator::Ell(G4double L) 348 { 328 { 349 G4double ans = 0.; 329 G4double ans = 0.; 350 for (G4int i = 0; i < nlev; ++i) { << 330 for(G4int i=0; i<nlev; ++i) { 351 if (sternf[i] > 0. && (sternEbar[i] > 0. | << 331 if(sternf[i] > 0. && (sternEbar[i] > 0. || L != 0.)) { 352 ans += sternf[i] / (gpow->powN(sternEbar << 332 ans += sternf[i]/(gpow->powN(sternEbar[i], 2) + L*L); 353 } 333 } 354 } 334 } 355 if (fConductivity > 0. && L != 0.) { << 356 ans += fConductivity / (L * L); << 357 } << 358 ans -= gpow->powZ(10, -2 * sternx); 335 ans -= gpow->powZ(10, -2 * sternx); 359 return ans; 336 return ans; 360 } 337 } 361 338 362 /** 339 /** 363 * Given the Sternheimer parameter l (called ' << 340 * Given the Sternheimer parameter l^2 (called 'sternL' here), and that 364 * the l_i and adjusted energies have been fou 341 * the l_i and adjusted energies have been found with SetupFermiDeltaCalc(), 365 * return the value of delta. Helper function 342 * return the value of delta. Helper function for DoFermiDeltaCalc(). 366 */ 343 */ 367 G4double G4DensityEffectCalculator::DeltaOnceS 344 G4double G4DensityEffectCalculator::DeltaOnceSolved(G4double sternL) 368 { 345 { 369 G4double ans = 0.; 346 G4double ans = 0.; 370 for (G4int i = 0; i < nlev; ++i) { << 347 for(G4int i=0; i<nlev; ++i) { 371 if (sternf[i] > 0.) { << 348 if(sternf[i] > 0.) { 372 ans += sternf[i] * << 349 ans += sternf[i] * G4Log((gpow->powN(sternl[i], 2) 373 G4Log((gpow->powN(sternl[i], 2) + << 350 + gpow->powN(sternL, 2))/gpow->powN(sternl[i], 2)); 374 } 351 } 375 } 352 } 376 // sternl for the conduction electrons is sq << 353 ans -= gpow->powN(sternL, 2)/(1 + gpow->powZ(10, 2 * sternx)); 377 // no factor of 2./3 as with the other level << 378 if (fConductivity > 0) { << 379 ans += fConductivity * G4Log((fConductivit << 380 } << 381 ans -= gpow->powN(sternL, 2) / (1 + gpow->po << 382 return ans; 354 return ans; 383 } 355 } 384 356