// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // // $Id: G4IonisParamElm.cc,v 1.15 2008/06/03 14:30:10 vnivanch Exp $ // GEANT4 tag $Name: geant4-09-02-ref-02 $ // // //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... // // 06-09-04, Update calculated values after any change of ionisation // potential change. V.Ivanchenko // 17-10-02, Fix excitation energy interpolation. V.Ivanchenko // 08-03-01, correct handling of fShellCorrectionVector. M.Maire // 22-11-00, tabulation of ionisation potential from // the ICRU Report N#37. V.Ivanchenko // 09-07-98, data moved from G4Element. M.Maire // //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... #include "G4IonisParamElm.hh" //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... G4IonisParamElm::G4IonisParamElm(G4double Z) { if (Z < 1.) G4Exception (" ERROR! It is not allowed to create an Element with Z < 1" ); // some basic functions of the atomic number fZ = Z; fZ3 = std::pow(fZ, 1./3.); fZZ3 = std::pow(fZ*(fZ+1.), 1./3.); flogZ3 = std::log(fZ)/3.; // Parameters for energy loss by ionisation // Mean excitation energy // from "Stopping Powers for Electrons and Positrons" // ICRU Report N#37, 1984 (energy in eV) static double exc[100] = { 21.8, 20.9, 13.3, 15.9, 15.2, 13.0, 11.7, 11.9, 11.5, 13.7, 13.6, 13.0, 12.8, 12.4, 11.5, 11.3, 10.2, 10.4, 10.0, 9.6, 10.3, 10.6, 10.7, 10.7, 10.9, 11.0, 11.0, 11.1, 11.1, 11.0, 10.8, 10.4, 10.5, 10.2, 9.8, 9.8, 9.8, 9.6, 9.7, 9.8, 10.2, 10.1, 10.0, 10.0, 10.0, 10.2, 10.0, 9.8, 10.0, 9.8, 9.5, 9.3, 9.3, 8.9, 8.9, 8.8, 8.8, 8.8, 9.1, 9.1, 9.2, 9.3, 9.2, 9.2, 9.4, 9.5, 9.7, 9.7, 9.8, 9.8, 9.8, 9.8, 9.8, 9.8, 9.8, 9.8, 9.8, 10.1, 10.0, 10.0, 10.0, 10.0, 9.9, 9.9, 9.7, 9.2, 9.5, 9.4, 9.4, 9.4, 9.6, 9.7, 9.7, 9.8, 9.8, 9.8, 9.8, 9.9, 9.9, 9.9 }; G4int iz = (G4int)Z - 1 ; if(0 > iz) iz = 0; else if(99 < iz) iz = 99 ; fMeanExcitationEnergy = fZ * exc[iz] * eV ; // compute parameters for ion transport // The aproximation from: // J.F.Ziegler, J.P. Biersack, U. Littmark // The Stopping and Range of Ions in Matter, // Vol.1, Pergamon Press, 1985 // Fast ions or hadrons if(91 < iz) iz = 91; static G4double vFermi[92] = { 1.0309, 0.15976, 0.59782, 1.0781, 1.0486, 1.0, 1.058, 0.93942, 0.74562, 0.3424, 0.45259, 0.71074, 0.90519, 0.97411, 0.97184, 0.89852, 0.70827, 0.39816, 0.36552, 0.62712, 0.81707, 0.9943, 1.1423, 1.2381, 1.1222, 0.92705, 1.0047, 1.2, 1.0661, 0.97411, 0.84912, 0.95, 1.0903, 1.0429, 0.49715, 0.37755, 0.35211, 0.57801, 0.77773, 1.0207, 1.029, 1.2542, 1.122, 1.1241, 1.0882, 1.2709, 1.2542, 0.90094, 0.74093, 0.86054, 0.93155, 1.0047, 0.55379, 0.43289, 0.32636, 0.5131, 0.695, 0.72591, 0.71202, 0.67413, 0.71418, 0.71453, 0.5911, 0.70263, 0.68049, 0.68203, 0.68121, 0.68532, 0.68715, 0.61884, 0.71801, 0.83048, 1.1222, 1.2381, 1.045, 1.0733, 1.0953, 1.2381, 1.2879, 0.78654, 0.66401, 0.84912, 0.88433, 0.80746, 0.43357, 0.41923, 0.43638, 0.51464, 0.73087, 0.81065, 1.9578, 1.0257} ; static G4double lFactor[92] = { 1.0, 1.0, 1.1, 1.06, 1.01, 1.03, 1.04, 0.99, 0.95, 0.9, 0.82, 0.81, 0.83, 0.88, 1.0, 0.95, 0.97, 0.99, 0.98, 0.97, 0.98, 0.97, 0.96, 0.93, 0.91, 0.9, 0.88, 0.9, 0.9, 0.9, 0.9, 0.85, 0.9, 0.9, 0.91, 0.92, 0.9, 0.9, 0.9, 0.9, 0.9, 0.88, 0.9, 0.88, 0.88, 0.9, 0.9, 0.88, 0.9, 0.9, 0.9, 0.9, 0.96, 1.2, 0.9, 0.88, 0.88, 0.85, 0.9, 0.9, 0.92, 0.95, 0.99, 1.03, 1.05, 1.07, 1.08, 1.1, 1.08, 1.08, 1.08, 1.08, 1.09, 1.09, 1.1, 1.11, 1.12, 1.13, 1.14, 1.15, 1.17, 1.2, 1.18, 1.17, 1.17, 1.16, 1.16, 1.16, 1.16, 1.16, 1.16, 1.16} ; fVFermi = vFermi[iz]; fLFactor = lFactor[iz]; // obsolete parameters for ionisation fTau0 = 0.1*fZ3*MeV/proton_mass_c2; fTaul = 2.*MeV/proton_mass_c2; // compute the Bethe-Bloch formula for energy = fTaul*particle mass G4double rate = fMeanExcitationEnergy/electron_mass_c2 ; G4double w = fTaul*(fTaul+2.) ; fBetheBlochLow = (fTaul+1.)*(fTaul+1.)*std::log(2.*w/rate)/w - 1. ; fBetheBlochLow = 2.*fZ*twopi_mc2_rcl2*fBetheBlochLow ; fClow = std::sqrt(fTaul)*fBetheBlochLow; fAlow = 6.458040 * fClow/fTau0; G4double Taum = 0.035*fZ3*MeV/proton_mass_c2; fBlow =-3.229020*fClow/(fTau0*std::sqrt(Taum)); // Shell correction parameterization fShellCorrectionVector = new G4double[3]; //[3] rate = 0.001*fMeanExcitationEnergy/eV; G4double rate2 = rate*rate; /* fShellCorrectionVector[0] = ( 1.10289e5 + 5.14781e8*rate)*rate2 ; fShellCorrectionVector[1] = ( 7.93805e3 - 2.22565e7*rate)*rate2 ; fShellCorrectionVector[2] = (-9.92256e1 + 2.10823e5*rate)*rate2 ; */ fShellCorrectionVector[0] = ( 0.422377 + 3.858019*rate)*rate2 ; fShellCorrectionVector[1] = ( 0.0304043 - 0.1667989*rate)*rate2 ; fShellCorrectionVector[2] = (-0.00038106 + 0.00157955*rate)*rate2 ; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... // Fake default constructor - sets only member data and allocates memory // for usage restricted to object persistency G4IonisParamElm::G4IonisParamElm(__void__&) : fShellCorrectionVector(0) { } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... G4IonisParamElm::~G4IonisParamElm() { if (fShellCorrectionVector) delete [] fShellCorrectionVector; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... G4IonisParamElm::G4IonisParamElm(G4IonisParamElm& right) { *this = right; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... const G4IonisParamElm& G4IonisParamElm::operator=(const G4IonisParamElm& right) { if (this != &right) { fZ = right.fZ; fZ3 = right.fZ3; fZZ3 = right.fZZ3; flogZ3 = right.flogZ3; fTau0 = right.fTau0; fTaul = right.fTaul; fBetheBlochLow = right.fBetheBlochLow; fAlow = right.fAlow; fBlow = right.fBlow; fClow = right.fClow; fMeanExcitationEnergy = right.fMeanExcitationEnergy; if (fShellCorrectionVector) delete [] fShellCorrectionVector; fShellCorrectionVector = new G4double[3]; fShellCorrectionVector[0] = right.fShellCorrectionVector[0]; fShellCorrectionVector[1] = right.fShellCorrectionVector[1]; fShellCorrectionVector[2] = right.fShellCorrectionVector[2]; } return *this; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... G4int G4IonisParamElm::operator==(const G4IonisParamElm& right) const { return (this == (G4IonisParamElm *) &right); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo.... G4int G4IonisParamElm::operator!=(const G4IonisParamElm& right) const { return (this != (G4IonisParamElm *) &right); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....