Geant4 Cross Reference |
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 // 09-07-98, data moved from G4Material, M.Maire 27 // 18-07-98, bug corrected in ComputeDensityEffect() for gas 28 // 16-01-01, bug corrected in ComputeDensityEffect() E100eV (L.Urban) 29 // 08-02-01, fShellCorrectionVector correctly handled (mma) 30 // 28-10-02, add setMeanExcitationEnergy (V.Ivanchenko) 31 // 06-09-04, factor 2 to shell correction term (V.Ivanchenko) 32 // 10-05-05, add a missing coma in FindMeanExcitationEnergy() - Bug#746 (mma) 33 // 27-09-07, add computation of parameters for ions (V.Ivanchenko) 34 // 04-03-08, remove reference to G4NistManager. Add fBirks constant (mma) 35 // 30-10-09, add G4DensityEffectData class and density effect computation (VI) 36 37 #include "G4IonisParamMat.hh" 38 39 #include "G4AtomicShells.hh" 40 #include "G4AutoLock.hh" 41 #include "G4DensityEffectData.hh" 42 #include "G4Exp.hh" 43 #include "G4Log.hh" 44 #include "G4Material.hh" 45 #include "G4NistManager.hh" 46 #include "G4PhysicalConstants.hh" 47 #include "G4Pow.hh" 48 #include "G4SystemOfUnits.hh" 49 50 G4DensityEffectData* G4IonisParamMat::fDensityData = nullptr; 51 52 namespace 53 { 54 G4Mutex ionisMutex = G4MUTEX_INITIALIZER; 55 } 56 57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 58 59 G4IonisParamMat::G4IonisParamMat(const G4Material* material) : fMaterial(material) 60 { 61 fBirks = 0.; 62 fMeanEnergyPerIon = 0.0; 63 twoln10 = 2. * G4Pow::GetInstance()->logZ(10); 64 65 // minimal set of default parameters for density effect 66 fCdensity = 0.0; 67 fD0density = 0.0; 68 fAdjustmentFactor = 1.0; 69 if (fDensityData == nullptr) { 70 fDensityData = new G4DensityEffectData(); 71 } 72 fDensityEffectCalc = nullptr; 73 74 // compute parameters 75 ComputeMeanParameters(); 76 ComputeDensityEffectParameters(); 77 ComputeFluctModel(); 78 ComputeIonParameters(); 79 } 80 81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 82 83 G4IonisParamMat::~G4IonisParamMat() 84 { 85 delete fDensityEffectCalc; 86 delete[] fShellCorrectionVector; 87 delete fDensityData; 88 fDensityData = nullptr; 89 fShellCorrectionVector = nullptr; 90 fDensityEffectCalc = nullptr; 91 } 92 93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 94 95 G4double G4IonisParamMat::GetDensityCorrection(G4double x) const 96 { 97 // x = log10(beta*gamma) 98 G4double y = 0.0; 99 if (x < fX0density) { 100 if (fD0density > 0.0) { 101 y = fD0density * G4Exp(twoln10 * (x - fX0density)); 102 } 103 } 104 else if (x >= fX1density) { 105 y = twoln10 * x - fCdensity; 106 } 107 else { 108 y = twoln10 * x - fCdensity + fAdensity * G4Exp(G4Log(fX1density - x) * fMdensity); 109 } 110 return y; 111 } 112 113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 114 115 void G4IonisParamMat::ComputeMeanParameters() 116 { 117 // compute mean excitation energy and shell correction vector 118 fTaul = (*(fMaterial->GetElementVector()))[0]->GetIonisation()->GetTaul(); 119 120 std::size_t nElements = fMaterial->GetNumberOfElements(); 121 const G4ElementVector* elmVector = fMaterial->GetElementVector(); 122 const G4double* nAtomsPerVolume = fMaterial->GetVecNbOfAtomsPerVolume(); 123 124 fMeanExcitationEnergy = FindMeanExcitationEnergy(fMaterial); 125 fLogMeanExcEnergy = 0.; 126 127 // Chemical formula defines mean excitation energy 128 if (fMeanExcitationEnergy > 0.0) { 129 fLogMeanExcEnergy = G4Log(fMeanExcitationEnergy); 130 131 // Compute average 132 } 133 else { 134 for (std::size_t i = 0; i < nElements; ++i) { 135 const G4Element* elm = (*elmVector)[i]; 136 fLogMeanExcEnergy += 137 nAtomsPerVolume[i] * elm->GetZ() * G4Log(elm->GetIonisation()->GetMeanExcitationEnergy()); 138 } 139 fLogMeanExcEnergy /= fMaterial->GetTotNbOfElectPerVolume(); 140 fMeanExcitationEnergy = G4Exp(fLogMeanExcEnergy); 141 } 142 143 fShellCorrectionVector = new G4double[3]; 144 145 for (G4int j = 0; j <= 2; ++j) { 146 fShellCorrectionVector[j] = 0.; 147 148 for (std::size_t k = 0; k < nElements; ++k) { 149 fShellCorrectionVector[j] += 150 nAtomsPerVolume[k] * (((*elmVector)[k])->GetIonisation()->GetShellCorrectionVector())[j]; 151 } 152 fShellCorrectionVector[j] *= 2.0 / fMaterial->GetTotNbOfElectPerVolume(); 153 } 154 } 155 156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 157 158 G4DensityEffectData* G4IonisParamMat::GetDensityEffectData() { 159 return fDensityData; 160 } 161 162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 163 164 void G4IonisParamMat::ComputeDensityEffectParameters() 165 { 166 G4State State = fMaterial->GetState(); 167 G4double density = fMaterial->GetDensity(); 168 169 // Check if density effect data exist in the table 170 // R.M. Sternheimer, Atomic Data and Nuclear Data Tables, 30: 261 (1984) 171 // or is assign to one of data set in this table 172 G4int idx = fDensityData->GetIndex(fMaterial->GetName()); 173 auto nelm = (G4int)fMaterial->GetNumberOfElements(); 174 G4int Z0 = ((*(fMaterial->GetElementVector()))[0])->GetZasInt(); 175 const G4Material* bmat = fMaterial->GetBaseMaterial(); 176 G4NistManager* nist = G4NistManager::Instance(); 177 178 // arbitrary empirical limits 179 // parameterisation with very different density is not applicable 180 static const G4double corrmax = 1.; 181 static const G4double massfracmax = 0.9; 182 183 // for simple non-NIST materials 184 G4double corr = 0.0; 185 186 if (idx < 0 && 1 == nelm) { 187 G4int z = (1 == Z0 && State == kStateLiquid) ? 0 : Z0; 188 idx = fDensityData->GetElementIndex(z); 189 190 // Correction for base material or for non-nominal density 191 // Except cases of very different density defined in user code 192 if (idx >= 0 && 0 < z) { 193 G4double dens = nist->GetNominalDensity(Z0); 194 if (dens <= 0.0) { 195 idx = -1; 196 } 197 else { 198 corr = G4Log(dens / density); 199 if (std::abs(corr) > corrmax) { 200 idx = -1; 201 } 202 } 203 } 204 } 205 // for base material case 206 if (idx < 0 && nullptr != bmat) { 207 idx = fDensityData->GetIndex(bmat->GetName()); 208 if (idx >= 0) { 209 corr = G4Log(bmat->GetDensity() / density); 210 if (std::abs(corr) > corrmax) { 211 idx = -1; 212 } 213 } 214 } 215 216 // for compound non-NIST materials with one element dominating 217 if (idx < 0 && 1 < nelm) { 218 const G4double tot = fMaterial->GetTotNbOfAtomsPerVolume(); 219 for (G4int i = 0; i < nelm; ++i) { 220 const G4double frac = fMaterial->GetVecNbOfAtomsPerVolume()[i] / tot; 221 if (frac > massfracmax) { 222 Z0 = ((*(fMaterial->GetElementVector()))[i])->GetZasInt(); 223 idx = fDensityData->GetElementIndex(Z0); 224 G4double dens = nist->GetNominalDensity(Z0); 225 if (idx >= 0 && dens > 0.0) { 226 corr = G4Log(dens / density); 227 if (std::abs(corr) > corrmax) { 228 idx = -1; 229 } 230 else { 231 break; 232 } 233 } 234 } 235 } 236 } 237 238 if (idx >= 0) { 239 // Take parameters for the density effect correction from 240 // R.M. Sternheimer et al. Density Effect For The Ionization Loss 241 // of Charged Particles in Various Substances. 242 // Atom. Data Nucl. Data Tabl. 30 (1984) 261-271. 243 244 fCdensity = fDensityData->GetCdensity(idx); 245 fMdensity = fDensityData->GetMdensity(idx); 246 fAdensity = fDensityData->GetAdensity(idx); 247 fX0density = fDensityData->GetX0density(idx); 248 fX1density = fDensityData->GetX1density(idx); 249 fD0density = fDensityData->GetDelta0density(idx); 250 fPlasmaEnergy = fDensityData->GetPlasmaEnergy(idx); 251 fAdjustmentFactor = fDensityData->GetAdjustmentFactor(idx); 252 253 // parameter C is computed and not taken from Sternheimer tables 254 // fCdensity = 1. + 2*G4Log(fMeanExcitationEnergy/fPlasmaEnergy); 255 // G4cout << "IonisParamMat: " << fMaterial->GetName() 256 // << " Cst= " << Cdensity << " C= " << fCdensity << G4endl; 257 258 // correction on nominal density 259 fCdensity += corr; 260 fX0density += corr / twoln10; 261 fX1density += corr / twoln10; 262 } 263 else { 264 static const G4double Cd2 = 4 * CLHEP::pi * CLHEP::hbarc_squared * CLHEP::classic_electr_radius; 265 fPlasmaEnergy = std::sqrt(Cd2 * fMaterial->GetTotNbOfElectPerVolume()); 266 267 // Compute parameters for the density effect correction in DE/Dx formula. 268 // The parametrization is from R.M. Sternheimer, Phys. Rev.B,3:3681 (1971) 269 G4int icase; 270 271 fCdensity = 1. + 2 * G4Log(fMeanExcitationEnergy / fPlasmaEnergy); 272 // 273 // condensed materials 274 // 275 if ((State == kStateSolid) || (State == kStateLiquid)) { 276 static const G4double E100eV = 100. * CLHEP::eV; 277 static const G4double ClimiS[] = {3.681, 5.215}; 278 static const G4double X0valS[] = {1.0, 1.5}; 279 static const G4double X1valS[] = {2.0, 3.0}; 280 281 if (fMeanExcitationEnergy < E100eV) { 282 icase = 0; 283 } 284 else { 285 icase = 1; 286 } 287 288 if (fCdensity < ClimiS[icase]) { 289 fX0density = 0.2; 290 } 291 else { 292 fX0density = 0.326 * fCdensity - X0valS[icase]; 293 } 294 295 fX1density = X1valS[icase]; 296 fMdensity = 3.0; 297 298 // special: Hydrogen 299 if (1 == nelm && 1 == Z0) { 300 fX0density = 0.425; 301 fX1density = 2.0; 302 fMdensity = 5.949; 303 } 304 } 305 else { 306 // 307 // gases 308 // 309 fMdensity = 3.; 310 fX1density = 4.0; 311 312 if (fCdensity <= 10.) { 313 fX0density = 1.6; 314 } 315 else if (fCdensity <= 10.5) { 316 fX0density = 1.7; 317 } 318 else if (fCdensity <= 11.0) { 319 fX0density = 1.8; 320 } 321 else if (fCdensity <= 11.5) { 322 fX0density = 1.9; 323 } 324 else if (fCdensity <= 12.25) { 325 fX0density = 2.0; 326 } 327 else if (fCdensity <= 13.804) { 328 fX0density = 2.0; 329 fX1density = 5.0; 330 } 331 else { 332 fX0density = 0.326 * fCdensity - 2.5; 333 fX1density = 5.0; 334 } 335 336 // special: Hydrogen 337 if (1 == nelm && 1 == Z0) { 338 fX0density = 1.837; 339 fX1density = 3.0; 340 fMdensity = 4.754; 341 } 342 343 // special: Helium 344 if (1 == nelm && 2 == Z0) { 345 fX0density = 2.191; 346 fX1density = 3.0; 347 fMdensity = 3.297; 348 } 349 } 350 } 351 352 // change parameters if the gas is not in STP. 353 // For the correction the density(STP) is needed. 354 // Density(STP) is calculated here : 355 356 if (State == kStateGas) { 357 G4double Pressure = fMaterial->GetPressure(); 358 G4double Temp = fMaterial->GetTemperature(); 359 360 G4double DensitySTP = density * STP_Pressure * Temp / (Pressure * NTP_Temperature); 361 362 G4double ParCorr = G4Log(density / DensitySTP); 363 364 fCdensity -= ParCorr; 365 fX0density -= ParCorr / twoln10; 366 fX1density -= ParCorr / twoln10; 367 } 368 369 // fAdensity parameter can be fixed for not conductive materials 370 if (0.0 == fD0density) { 371 G4double Xa = fCdensity / twoln10; 372 fAdensity = twoln10 * (Xa - fX0density) / std::pow((fX1density - fX0density), fMdensity); 373 } 374 } 375 376 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 377 378 void G4IonisParamMat::ComputeFluctModel() 379 { 380 // compute parameters for the energy loss fluctuation model 381 // needs an 'effective Z' 382 G4double Zeff = 0.; 383 for (std::size_t i = 0; i < fMaterial->GetNumberOfElements(); ++i) { 384 Zeff += (fMaterial->GetFractionVector())[i] * (*(fMaterial->GetElementVector()))[i]->GetZ(); 385 } 386 if (Zeff > 2.1) { 387 fF2fluct = 2.0 / Zeff; 388 fF1fluct = 1. - fF2fluct; 389 fEnergy2fluct = 10. * Zeff * Zeff * CLHEP::eV; 390 fLogEnergy2fluct = G4Log(fEnergy2fluct); 391 fLogEnergy1fluct = (fLogMeanExcEnergy - fF2fluct * fLogEnergy2fluct) / fF1fluct; 392 } else if (Zeff > 1.1) { 393 fF2fluct = 0.0; 394 fF1fluct = 1.0; 395 fEnergy2fluct = 40. * CLHEP::eV; 396 fLogEnergy2fluct = G4Log(fEnergy2fluct); 397 fLogEnergy1fluct = fLogMeanExcEnergy; 398 } else { 399 fF2fluct = 0.0; 400 fF1fluct = 1.0; 401 fEnergy2fluct = 10. * CLHEP::eV; 402 fLogEnergy2fluct = G4Log(fEnergy2fluct); 403 fLogEnergy1fluct = fLogMeanExcEnergy; 404 } 405 fEnergy1fluct = G4Exp(fLogEnergy1fluct); 406 fEnergy0fluct = 10. * CLHEP::eV; 407 fRateionexcfluct = 0.4; 408 } 409 410 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 411 412 void G4IonisParamMat::ComputeIonParameters() 413 { 414 // get elements in the actual material, 415 const G4ElementVector* theElementVector = fMaterial->GetElementVector(); 416 const G4double* theAtomicNumDensityVector = fMaterial->GetAtomicNumDensityVector(); 417 const auto NumberOfElements = (G4int)fMaterial->GetNumberOfElements(); 418 419 // loop for the elements in the material 420 // to find out average values Z, vF, lF 421 G4double z(0.0), vF(0.0), lF(0.0), a23(0.0); 422 423 G4Pow* g4pow = G4Pow::GetInstance(); 424 if (1 == NumberOfElements) { 425 const G4Element* element = (*theElementVector)[0]; 426 z = element->GetZ(); 427 vF = element->GetIonisation()->GetFermiVelocity(); 428 lF = element->GetIonisation()->GetLFactor(); 429 a23 = 1.0 / g4pow->A23(element->GetN()); 430 } 431 else { 432 G4double norm(0.0); 433 for (G4int iel = 0; iel < NumberOfElements; ++iel) { 434 const G4Element* element = (*theElementVector)[iel]; 435 const G4double weight = theAtomicNumDensityVector[iel]; 436 norm += weight; 437 z += element->GetZ() * weight; 438 vF += element->GetIonisation()->GetFermiVelocity() * weight; 439 lF += element->GetIonisation()->GetLFactor() * weight; 440 a23 += weight / g4pow->A23(element->GetN()); 441 } 442 if (norm > 0.0) { norm = 1.0/norm; } 443 z *= norm; 444 vF *= norm; 445 lF *= norm; 446 a23 *= norm; 447 } 448 fZeff = z; 449 fLfactor = lF; 450 fFermiEnergy = 25. * CLHEP::keV * vF * vF; 451 fInvA23 = a23; 452 } 453 454 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 455 456 void G4IonisParamMat::SetMeanExcitationEnergy(G4double value) 457 { 458 if (value == fMeanExcitationEnergy || value <= 0.0) { 459 return; 460 } 461 if (G4NistManager::Instance()->GetVerbose() > 1) { 462 G4cout << "G4Material: Mean excitation energy is changed for " << fMaterial->GetName() 463 << " Iold= " << fMeanExcitationEnergy / CLHEP::eV << "eV; Inew= " << value / eV << " eV;" 464 << G4endl; 465 } 466 467 fMeanExcitationEnergy = value; 468 469 // add corrections to density effect 470 G4double newlog = G4Log(value); 471 G4double corr = 2 * (newlog - fLogMeanExcEnergy); 472 fCdensity += corr; 473 fX0density += corr / twoln10; 474 fX1density += corr / twoln10; 475 476 // recompute parameters of fluctuation model 477 fLogMeanExcEnergy = newlog; 478 ComputeFluctModel(); 479 } 480 481 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 482 483 void G4IonisParamMat::SetDensityEffectParameters( 484 G4double cd, G4double md, G4double ad, G4double x0, G4double x1, G4double d0) 485 { 486 // no check on consistence of user parameters 487 G4AutoLock l(&ionisMutex); 488 fCdensity = cd; 489 fMdensity = md; 490 fAdensity = ad; 491 fX0density = x0; 492 fX1density = x1; 493 fD0density = d0; 494 l.unlock(); 495 } 496 497 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 498 499 void G4IonisParamMat::SetDensityEffectParameters(const G4Material* bmat) 500 { 501 G4AutoLock l(&ionisMutex); 502 const G4IonisParamMat* ipm = bmat->GetIonisation(); 503 fCdensity = ipm->GetCdensity(); 504 fMdensity = ipm->GetMdensity(); 505 fAdensity = ipm->GetAdensity(); 506 fX0density = ipm->GetX0density(); 507 fX1density = ipm->GetX1density(); 508 fD0density = ipm->GetD0density(); 509 510 G4double corr = G4Log(bmat->GetDensity() / fMaterial->GetDensity()); 511 fCdensity += corr; 512 fX0density += corr / twoln10; 513 fX1density += corr / twoln10; 514 l.unlock(); 515 } 516 517 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 518 519 void G4IonisParamMat::ComputeDensityEffectOnFly(G4bool val) 520 { 521 if (val) { 522 if (nullptr == fDensityEffectCalc) { 523 G4int n = 0; 524 for (std::size_t i = 0; i < fMaterial->GetNumberOfElements(); ++i) { 525 const G4int Z = fMaterial->GetElement((G4int)i)->GetZasInt(); 526 n += G4AtomicShells::GetNumberOfShells(Z); 527 } 528 // The last level is the conduction level. If this is *not* a conductor, 529 // make a dummy conductor level with zero electrons in it. 530 fDensityEffectCalc = new G4DensityEffectCalculator(fMaterial, n); 531 } 532 } 533 else { 534 delete fDensityEffectCalc; 535 fDensityEffectCalc = nullptr; 536 } 537 } 538 539 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 540 541 G4double G4IonisParamMat::FindMeanExcitationEnergy(const G4Material* mat) const 542 { 543 G4double res = 0.0; 544 // data from density effect data 545 if (fDensityData != nullptr) { 546 G4int idx = fDensityData->GetIndex(mat->GetName()); 547 if (idx >= 0) { 548 res = fDensityData->GetMeanIonisationPotential(idx); 549 } 550 } 551 552 // The data on mean excitation energy for compaunds 553 // from "Stopping Powers for Electrons and Positrons" 554 // ICRU Report N#37, 1984 (energy in eV) 555 // this value overwrites Density effect data 556 G4String chFormula = mat->GetChemicalFormula(); 557 if (! chFormula.empty()) { 558 static const size_t numberOfMolecula = 54; 559 // clang-format off 560 static const G4String name[numberOfMolecula] = { 561 // gas 0 - 12 562 "NH_3", "C_4H_10", "CO_2", "C_2H_6", "C_7H_16-Gas", 563 // "G4_AMMONIA", "G4_BUTANE","G4_CARBON_DIOXIDE","G4_ETHANE", "G4_N-HEPTANE" 564 "C_6H_14-Gas", "CH_4", "NO", "N_2O", "C_8H_18-Gas", 565 // "G4_N-HEXANE" , "G4_METHANE", "x", "G4_NITROUS_OXIDE", "G4_OCTANE" 566 "C_5H_12-Gas", "C_3H_8", "H_2O-Gas", 567 // "G4_N-PENTANE", "G4_PROPANE", "G4_WATER_VAPOR" 568 569 // liquid 13 - 39 570 "C_3H_6O", "C_6H_5NH_2", "C_6H_6", "C_4H_9OH", "CCl_4", 571 //"G4_ACETONE","G4_ANILINE","G4_BENZENE","G4_N-BUTYL_ALCOHOL","G4_CARBON_TETRACHLORIDE" 572 "C_6H_5Cl", "CHCl_3", "C_6H_12", "C_6H_4Cl_2", "C_4Cl_2H_8O", 573 //"G4_CHLOROBENZENE","G4_CHLOROFORM","G4_CYCLOHEXANE","G4_1,2-DICHLOROBENZENE", 574 //"G4_DICHLORODIETHYL_ETHER" 575 "C_2Cl_2H_4", "(C_2H_5)_2O", "C_2H_5OH", "C_3H_5(OH)_3","C_7H_16", 576 //"G4_1,2-DICHLOROETHANE","G4_DIETHYL_ETHER","G4_ETHYL_ALCOHOL","G4_GLYCEROL","G4_N-HEPTANE" 577 "C_6H_14", "CH_3OH", "C_6H_5NO_2","C_5H_12", "C_3H_7OH", 578 //"G4_N-HEXANE","G4_METHANOL","G4_NITROBENZENE","G4_N-PENTANE","G4_N-PROPYL_ALCOHOL", 579 "C_5H_5N", "C_8H_8", "C_2Cl_4", "C_7H_8", "C_2Cl_3H", 580 //"G4_PYRIDINE","G4_POLYSTYRENE","G4_TETRACHLOROETHYLENE","G4_TOLUENE","G4_TRICHLOROETHYLENE" 581 "H_2O", "C_8H_10", 582 // "G4_WATER", "G4_XYLENE" 583 584 // solid 40 - 53 585 "C_5H_5N_5", "C_5H_5N_5O", "(C_6H_11NO)-nylon", "C_25H_52", 586 // "G4_ADENINE", "G4_GUANINE", "G4_NYLON-6-6", "G4_PARAFFIN" 587 "(C_2H_4)-Polyethylene", "(C_5H_8O_2)-Polymethil_Methacrylate", 588 // "G4_ETHYLENE", "G4_PLEXIGLASS" 589 "(C_8H_8)-Polystyrene", "A-150-tissue", "Al_2O_3", "CaF_2", 590 // "G4_POLYSTYRENE", "G4_A-150_TISSUE", "G4_ALUMINUM_OXIDE", "G4_CALCIUM_FLUORIDE" 591 "LiF", "Photo_Emulsion", "(C_2F_4)-Teflon", "SiO_2" 592 // "G4_LITHIUM_FLUORIDE", "G4_PHOTO_EMULSION", "G4_TEFLON", "G4_SILICON_DIOXIDE" 593 } ; 594 595 static const G4double meanExcitation[numberOfMolecula] = { 596 597 53.7, 48.3, 85.0, 45.4, 49.2, 598 49.1, 41.7, 87.8, 84.9, 49.5, 599 48.2, 47.1, 71.6, 600 601 64.2, 66.2, 63.4, 59.9, 166.3, 602 89.1, 156.0, 56.4, 106.5, 103.3, 603 111.9, 60.0, 62.9, 72.6, 54.4, 604 54.0, 67.6, 75.8, 53.6, 61.1, 605 66.2, 64.0, 159.2, 62.5, 148.1, 606 75.0, 61.8, 607 608 71.4, 75.0, 63.9, 48.3, 57.4, 609 74.0, 68.7, 65.1, 145.2, 166., 610 94.0, 331.0, 99.1, 139.2 611 }; 612 // clang-format on 613 614 for (std::size_t i = 0; i < numberOfMolecula; ++i) { 615 if (chFormula == name[i]) { 616 res = meanExcitation[i] * CLHEP::eV; 617 break; 618 } 619 } 620 } 621 return res; 622 } 623