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 G4Element. M.Maire 27 // 22-11-00, tabulation of ionisation potential from 28 // the ICRU Report N#37. V.Ivanchenko 29 // 08-03-01, correct handling of fShellCorrectionVector. M.Maire 30 // 17-10-02, Fix excitation energy interpolation. V.Ivanchenko 31 // 06-09-04, Update calculated values after any change of ionisation 32 // potential change. V.Ivanchenko 33 // 29-04-10, Using G4Pow and mean ionisation energy from NIST V.Ivanchenko 34 // 27.10.11: new scheme for G4Exception (mma) 35 36 #include "G4IonisParamElm.hh" 37 38 #include "G4NistManager.hh" 39 #include "G4PhysicalConstants.hh" 40 #include "G4Pow.hh" 41 #include "G4SystemOfUnits.hh" 42 43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 44 45 G4IonisParamElm::G4IonisParamElm(G4double AtomNumber) 46 { 47 G4int Z = G4lrint(AtomNumber); 48 if (Z < 1) { 49 G4Exception("G4IonisParamElm::G4IonisParamElm()", "mat501", FatalException, 50 "It is not allowed to create an Element with Z<1"); 51 } 52 G4Pow* g4pow = G4Pow::GetInstance(); 53 54 // some basic functions of the atomic number 55 fZ = Z; 56 fZ3 = g4pow->Z13(Z); 57 fZZ3 = fZ3 * g4pow->Z13(Z + 1); 58 flogZ3 = g4pow->logZ(Z) / 3.; 59 60 fMeanExcitationEnergy = G4NistManager::Instance()->GetMeanIonisationEnergy(Z); 61 62 // compute parameters for ion transport 63 // The aproximation from: 64 // J.F.Ziegler, J.P. Biersack, U. Littmark 65 // The Stopping and Range of Ions in Matter, 66 // Vol.1, Pergamon Press, 1985 67 // Fast ions or hadrons 68 69 G4int iz = Z - 1; 70 if (91 < iz) { 71 iz = 91; 72 } 73 74 // clang-format off 75 static const G4double vFermi[92] = { 76 1.0309, 0.15976, 0.59782, 1.0781, 1.0486, 1.0, 1.058, 0.93942, 0.74562, 0.3424, 77 0.45259, 0.71074, 0.90519, 0.97411, 0.97184, 0.89852, 0.70827, 0.39816, 0.36552, 0.62712, 78 0.81707, 0.9943, 1.1423, 1.2381, 1.1222, 0.92705, 1.0047, 1.2, 1.0661, 0.97411, 79 0.84912, 0.95, 1.0903, 1.0429, 0.49715, 0.37755, 0.35211, 0.57801, 0.77773, 1.0207, 80 1.029, 1.2542, 1.122, 1.1241, 1.0882, 1.2709, 1.2542, 0.90094, 0.74093, 0.86054, 81 0.93155, 1.0047, 0.55379, 0.43289, 0.32636, 0.5131, 0.695, 0.72591, 0.71202, 0.67413, 82 0.71418, 0.71453, 0.5911, 0.70263, 0.68049, 0.68203, 0.68121, 0.68532, 0.68715, 0.61884, 83 0.71801, 0.83048, 1.1222, 1.2381, 1.045, 1.0733, 1.0953, 1.2381, 1.2879, 0.78654, 84 0.66401, 0.84912, 0.88433, 0.80746, 0.43357, 0.41923, 0.43638, 0.51464, 0.73087, 0.81065, 85 1.9578, 1.0257} ; 86 87 static const G4double lFactor[92] = { 88 1.0, 1.0, 1.1, 1.06, 1.01, 1.03, 1.04, 0.99, 0.95, 0.9, 89 0.82, 0.81, 0.83, 0.88, 1.0, 0.95, 0.97, 0.99, 0.98, 0.97, 90 0.98, 0.97, 0.96, 0.93, 0.91, 0.9, 0.88, 0.9, 0.9, 0.9, 91 0.9, 0.85, 0.9, 0.9, 0.91, 0.92, 0.9, 0.9, 0.9, 0.9, 92 0.9, 0.88, 0.9, 0.88, 0.88, 0.9, 0.9, 0.88, 0.9, 0.9, 93 0.9, 0.9, 0.96, 1.2, 0.9, 0.88, 0.88, 0.85, 0.9, 0.9, 94 0.92, 0.95, 0.99, 1.03, 1.05, 1.07, 1.08, 1.1, 1.08, 1.08, 95 1.08, 1.08, 1.09, 1.09, 1.1, 1.11, 1.12, 1.13, 1.14, 1.15, 96 1.17, 1.2, 1.18, 1.17, 1.17, 1.16, 1.16, 1.16, 1.16, 1.16, 97 1.16, 1.16} ; 98 // clang-format on 99 100 fVFermi = vFermi[iz]; 101 fLFactor = lFactor[iz]; 102 103 // obsolete parameters for ionisation 104 fTau0 = 0.1 * fZ3 * MeV / proton_mass_c2; 105 fTaul = 2. * MeV / proton_mass_c2; 106 107 // compute the Bethe-Bloch formula for energy = fTaul*particle mass 108 G4double rate = fMeanExcitationEnergy / electron_mass_c2; 109 G4double w = fTaul * (fTaul + 2.); 110 fBetheBlochLow = (fTaul + 1.) * (fTaul + 1.) * std::log(2. * w / rate) / w - 1.; 111 fBetheBlochLow = 2. * fZ * twopi_mc2_rcl2 * fBetheBlochLow; 112 113 fClow = std::sqrt(fTaul) * fBetheBlochLow; 114 fAlow = 6.458040 * fClow / fTau0; 115 G4double Taum = 0.035 * fZ3 * MeV / proton_mass_c2; 116 fBlow = -3.229020 * fClow / (fTau0 * std::sqrt(Taum)); 117 118 // Shell correction parameterization 119 fShellCorrectionVector = new G4double[3]; //[3] 120 rate = 0.001 * fMeanExcitationEnergy / eV; 121 G4double rate2 = rate * rate; 122 fShellCorrectionVector[0] = (0.422377 + 3.858019 * rate) * rate2; 123 fShellCorrectionVector[1] = (0.0304043 - 0.1667989 * rate) * rate2; 124 fShellCorrectionVector[2] = (-0.00038106 + 0.00157955 * rate) * rate2; 125 } 126 127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... 128 129 G4IonisParamElm::~G4IonisParamElm() { delete[] fShellCorrectionVector; } 130