[819] | 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 | // |
---|
| 29 | // GEANT4 Class file |
---|
| 30 | // |
---|
| 31 | // |
---|
| 32 | // File name: G4hBetheBlochModel |
---|
| 33 | // |
---|
| 34 | // Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch) |
---|
| 35 | // |
---|
| 36 | // Creation date: 20 July 2000 |
---|
| 37 | // |
---|
| 38 | // Modifications: |
---|
| 39 | // 20/07/2000 V.Ivanchenko First implementation |
---|
| 40 | // 03/10/2000 V.Ivanchenko clean up accoding to CodeWizard |
---|
| 41 | // |
---|
| 42 | // Class Description: |
---|
| 43 | // |
---|
| 44 | // Bethe-Bloch ionisation model |
---|
| 45 | // |
---|
| 46 | // Class Description: End |
---|
| 47 | // |
---|
| 48 | // ------------------------------------------------------------------- |
---|
| 49 | // |
---|
| 50 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 51 | |
---|
| 52 | #include "G4hBetheBlochModel.hh" |
---|
| 53 | #include "G4DynamicParticle.hh" |
---|
| 54 | #include "G4ParticleDefinition.hh" |
---|
| 55 | #include "G4Material.hh" |
---|
| 56 | #include "globals.hh" |
---|
| 57 | |
---|
| 58 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 59 | |
---|
| 60 | G4hBetheBlochModel::G4hBetheBlochModel(const G4String& name) |
---|
| 61 | : G4VLowEnergyModel(name), |
---|
| 62 | lowEnergyLimit(1.*MeV), |
---|
| 63 | highEnergyLimit(100.*GeV), |
---|
| 64 | twoln10(2.*std::log(10.)), |
---|
| 65 | bg2lim(0.0169), |
---|
| 66 | taulim(8.4146e-3) |
---|
| 67 | {;} |
---|
| 68 | |
---|
| 69 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 70 | |
---|
| 71 | G4hBetheBlochModel::~G4hBetheBlochModel() |
---|
| 72 | {;} |
---|
| 73 | |
---|
| 74 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 75 | |
---|
| 76 | G4double G4hBetheBlochModel::TheValue(const G4DynamicParticle* particle, |
---|
| 77 | const G4Material* material) |
---|
| 78 | { |
---|
| 79 | G4double energy = particle->GetKineticEnergy() ; |
---|
| 80 | G4double particleMass = particle->GetMass() ; |
---|
| 81 | |
---|
| 82 | G4double eloss = BetheBlochFormula(material,energy,particleMass) ; |
---|
| 83 | |
---|
| 84 | return eloss ; |
---|
| 85 | } |
---|
| 86 | |
---|
| 87 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 88 | |
---|
| 89 | G4double G4hBetheBlochModel::TheValue(const G4ParticleDefinition* aParticle, |
---|
| 90 | const G4Material* material, |
---|
| 91 | G4double kineticEnergy) |
---|
| 92 | { |
---|
| 93 | G4double particleMass = aParticle->GetPDGMass() ; |
---|
| 94 | G4double eloss = BetheBlochFormula(material,kineticEnergy,particleMass) ; |
---|
| 95 | |
---|
| 96 | return eloss ; |
---|
| 97 | } |
---|
| 98 | |
---|
| 99 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 100 | |
---|
| 101 | G4double G4hBetheBlochModel::HighEnergyLimit( |
---|
| 102 | const G4ParticleDefinition* , |
---|
| 103 | const G4Material* ) const |
---|
| 104 | { |
---|
| 105 | return highEnergyLimit ; |
---|
| 106 | } |
---|
| 107 | |
---|
| 108 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 109 | |
---|
| 110 | G4double G4hBetheBlochModel::LowEnergyLimit( |
---|
| 111 | const G4ParticleDefinition* aParticle, |
---|
| 112 | const G4Material* material) const |
---|
| 113 | { |
---|
| 114 | G4double taul = (material->GetIonisation()->GetTaul())* |
---|
| 115 | (aParticle->GetPDGMass()) ; |
---|
| 116 | return taul ; |
---|
| 117 | } |
---|
| 118 | |
---|
| 119 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 120 | |
---|
| 121 | G4double G4hBetheBlochModel::HighEnergyLimit( |
---|
| 122 | const G4ParticleDefinition* ) const |
---|
| 123 | { |
---|
| 124 | return highEnergyLimit ; |
---|
| 125 | } |
---|
| 126 | |
---|
| 127 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 128 | |
---|
| 129 | G4double G4hBetheBlochModel::LowEnergyLimit( |
---|
| 130 | const G4ParticleDefinition* ) const |
---|
| 131 | { |
---|
| 132 | return lowEnergyLimit ; |
---|
| 133 | } |
---|
| 134 | |
---|
| 135 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 136 | |
---|
| 137 | G4bool G4hBetheBlochModel::IsInCharge(const G4DynamicParticle* , |
---|
| 138 | const G4Material* ) const |
---|
| 139 | { |
---|
| 140 | return true ; |
---|
| 141 | } |
---|
| 142 | |
---|
| 143 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 144 | |
---|
| 145 | G4bool G4hBetheBlochModel::IsInCharge(const G4ParticleDefinition* , |
---|
| 146 | const G4Material* ) const |
---|
| 147 | { |
---|
| 148 | return true ; |
---|
| 149 | } |
---|
| 150 | |
---|
| 151 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 152 | |
---|
| 153 | G4double G4hBetheBlochModel::BetheBlochFormula( |
---|
| 154 | const G4Material* material, |
---|
| 155 | G4double kineticEnergy, |
---|
| 156 | G4double particleMass) const |
---|
| 157 | { |
---|
| 158 | // This member function is applied normally to proton/antiproton |
---|
| 159 | G4double ionloss ; |
---|
| 160 | |
---|
| 161 | G4double rateMass = electron_mass_c2/particleMass ; |
---|
| 162 | |
---|
| 163 | G4double taul = material->GetIonisation()->GetTaul() ; |
---|
| 164 | G4double tau = kineticEnergy/particleMass ; // tau is relative energy |
---|
| 165 | |
---|
| 166 | // It is not normal case for this function |
---|
| 167 | // for low energy parametrisation have to be applied |
---|
| 168 | if ( tau < taul ) tau = taul ; |
---|
| 169 | |
---|
| 170 | // some local variables |
---|
| 171 | |
---|
| 172 | G4double gamma,bg2,beta2,tmax,x,delta,sh ; |
---|
| 173 | G4double electronDensity = material->GetElectronDensity(); |
---|
| 174 | G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); |
---|
| 175 | G4double eexc2 = eexc*eexc ; |
---|
| 176 | G4double cden = material->GetIonisation()->GetCdensity(); |
---|
| 177 | G4double mden = material->GetIonisation()->GetMdensity(); |
---|
| 178 | G4double aden = material->GetIonisation()->GetAdensity(); |
---|
| 179 | G4double x0den = material->GetIonisation()->GetX0density(); |
---|
| 180 | G4double x1den = material->GetIonisation()->GetX1density(); |
---|
| 181 | G4double* shellCorrectionVector = |
---|
| 182 | material->GetIonisation()->GetShellCorrectionVector(); |
---|
| 183 | |
---|
| 184 | gamma = tau + 1.0 ; |
---|
| 185 | bg2 = tau*(tau+2.0) ; |
---|
| 186 | beta2 = bg2/(gamma*gamma) ; |
---|
| 187 | tmax = 2.*electron_mass_c2*bg2/(1.+2.*gamma*rateMass+rateMass*rateMass) ; |
---|
| 188 | |
---|
| 189 | ionloss = std::log(2.0*electron_mass_c2*bg2*tmax/eexc2)-2.0*beta2 ; |
---|
| 190 | |
---|
| 191 | // density correction |
---|
| 192 | x = std::log(bg2)/twoln10 ; |
---|
| 193 | if ( x < x0den ) { |
---|
| 194 | delta = 0.0 ; |
---|
| 195 | |
---|
| 196 | } else { |
---|
| 197 | delta = twoln10*x - cden ; |
---|
| 198 | if ( x < x1den ) delta += aden*std::pow((x1den-x),mden) ; |
---|
| 199 | } |
---|
| 200 | |
---|
| 201 | // shell correction |
---|
| 202 | sh = 0.0 ; |
---|
| 203 | x = 1.0 ; |
---|
| 204 | |
---|
| 205 | if ( bg2 > bg2lim ) { |
---|
| 206 | for (G4int k=0; k<=2; k++) { |
---|
| 207 | x *= bg2 ; |
---|
| 208 | sh += shellCorrectionVector[k]/x; |
---|
| 209 | } |
---|
| 210 | |
---|
| 211 | } else { |
---|
| 212 | for (G4int k=0; k<=2; k++) { |
---|
| 213 | x *= bg2lim ; |
---|
| 214 | sh += shellCorrectionVector[k]/x; |
---|
| 215 | } |
---|
| 216 | sh *= std::log(tau/taul)/std::log(taulim/taul) ; |
---|
| 217 | } |
---|
| 218 | |
---|
| 219 | // now compute the total ionization loss |
---|
| 220 | |
---|
| 221 | ionloss -= delta + sh ; |
---|
| 222 | ionloss *= twopi_mc2_rcl2*electronDensity/beta2 ; |
---|
| 223 | |
---|
| 224 | if ( ionloss < 0.0) ionloss = 0.0 ; |
---|
| 225 | |
---|
| 226 | return ionloss; |
---|
| 227 | } |
---|
| 228 | |
---|
| 229 | |
---|