[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 | // |
---|
[1196] | 26 | // $Id: G4MuBetheBlochModel.cc,v 1.27 2009/11/09 19:18:01 vnivanch Exp $ |
---|
[1337] | 27 | // GEANT4 tag $Name: geant4-09-04-beta-01 $ |
---|
[819] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4MuBetheBlochModel |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code |
---|
| 37 | // |
---|
| 38 | // Creation date: 09.08.2002 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko) |
---|
| 43 | // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko) |
---|
| 44 | // 27-01-03 Make models region aware (V.Ivanchenko) |
---|
| 45 | // 13-02-03 Add name (V.Ivanchenko) |
---|
| 46 | // 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I) |
---|
| 47 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) |
---|
| 48 | // 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko) |
---|
| 49 | // 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma) |
---|
| 50 | // |
---|
| 51 | |
---|
| 52 | // |
---|
| 53 | // ------------------------------------------------------------------- |
---|
| 54 | // |
---|
| 55 | |
---|
| 56 | |
---|
| 57 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 58 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 59 | |
---|
| 60 | #include "G4MuBetheBlochModel.hh" |
---|
| 61 | #include "Randomize.hh" |
---|
| 62 | #include "G4Electron.hh" |
---|
| 63 | #include "G4LossTableManager.hh" |
---|
| 64 | #include "G4EmCorrections.hh" |
---|
| 65 | #include "G4ParticleChangeForLoss.hh" |
---|
| 66 | |
---|
| 67 | G4double G4MuBetheBlochModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917, |
---|
| 68 | 0.7628, 0.8983, 0.9801 }; |
---|
| 69 | |
---|
| 70 | G4double G4MuBetheBlochModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813, |
---|
| 71 | 0.1569, 0.1112, 0.0506 }; |
---|
| 72 | |
---|
| 73 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 74 | |
---|
| 75 | using namespace std; |
---|
| 76 | |
---|
| 77 | G4MuBetheBlochModel::G4MuBetheBlochModel(const G4ParticleDefinition* p, |
---|
| 78 | const G4String& nam) |
---|
| 79 | : G4VEmModel(nam), |
---|
| 80 | particle(0), |
---|
| 81 | limitKinEnergy(100.*keV), |
---|
| 82 | logLimitKinEnergy(log(limitKinEnergy)), |
---|
| 83 | twoln10(2.0*log(10.0)), |
---|
| 84 | bg2lim(0.0169), |
---|
| 85 | taulim(8.4146e-3), |
---|
| 86 | alphaprime(fine_structure_const/twopi) |
---|
| 87 | { |
---|
| 88 | theElectron = G4Electron::Electron(); |
---|
| 89 | corr = G4LossTableManager::Instance()->EmCorrections(); |
---|
[1055] | 90 | fParticleChange = 0; |
---|
[819] | 91 | |
---|
| 92 | if(p) SetParticle(p); |
---|
| 93 | } |
---|
| 94 | |
---|
| 95 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 96 | |
---|
| 97 | G4MuBetheBlochModel::~G4MuBetheBlochModel() |
---|
| 98 | {} |
---|
| 99 | |
---|
| 100 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 101 | |
---|
[1055] | 102 | G4double G4MuBetheBlochModel::MinEnergyCut(const G4ParticleDefinition*, |
---|
| 103 | const G4MaterialCutsCouple* couple) |
---|
[819] | 104 | { |
---|
[1055] | 105 | return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy(); |
---|
[819] | 106 | } |
---|
| 107 | |
---|
[1055] | 108 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... |
---|
[819] | 109 | |
---|
[1055] | 110 | G4double G4MuBetheBlochModel::MaxSecondaryEnergy(const G4ParticleDefinition*, |
---|
| 111 | G4double kinEnergy) |
---|
[819] | 112 | { |
---|
[1055] | 113 | G4double tau = kinEnergy/mass; |
---|
| 114 | G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) / |
---|
| 115 | (1. + 2.0*(tau + 1.)*ratio + ratio*ratio); |
---|
| 116 | return tmax; |
---|
[819] | 117 | } |
---|
| 118 | |
---|
| 119 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 120 | |
---|
| 121 | void G4MuBetheBlochModel::Initialise(const G4ParticleDefinition* p, |
---|
| 122 | const G4DataVector&) |
---|
| 123 | { |
---|
| 124 | if(p) SetParticle(p); |
---|
[1055] | 125 | if(!fParticleChange) fParticleChange = GetParticleChangeForLoss(); |
---|
[819] | 126 | } |
---|
| 127 | |
---|
| 128 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 129 | |
---|
| 130 | G4double G4MuBetheBlochModel::ComputeCrossSectionPerElectron( |
---|
| 131 | const G4ParticleDefinition* p, |
---|
| 132 | G4double kineticEnergy, |
---|
| 133 | G4double cutEnergy, |
---|
| 134 | G4double maxKinEnergy) |
---|
| 135 | { |
---|
| 136 | G4double cross = 0.0; |
---|
| 137 | G4double tmax = MaxSecondaryEnergy(p, kineticEnergy); |
---|
| 138 | G4double maxEnergy = min(tmax,maxKinEnergy); |
---|
| 139 | if(cutEnergy < maxEnergy) { |
---|
| 140 | |
---|
| 141 | G4double totEnergy = kineticEnergy + mass; |
---|
| 142 | G4double energy2 = totEnergy*totEnergy; |
---|
| 143 | G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2; |
---|
| 144 | |
---|
| 145 | cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax |
---|
| 146 | + 0.5*(maxEnergy - cutEnergy)/energy2; |
---|
| 147 | |
---|
| 148 | // radiative corrections of R. Kokoulin |
---|
| 149 | if (maxEnergy > limitKinEnergy) { |
---|
| 150 | |
---|
| 151 | G4double logtmax = log(maxEnergy); |
---|
| 152 | G4double logtmin = log(max(cutEnergy,limitKinEnergy)); |
---|
| 153 | G4double logstep = logtmax - logtmin; |
---|
| 154 | G4double dcross = 0.0; |
---|
| 155 | |
---|
| 156 | for (G4int ll=0; ll<8; ll++) |
---|
| 157 | { |
---|
| 158 | G4double ep = exp(logtmin + xgi[ll]*logstep); |
---|
| 159 | G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2); |
---|
| 160 | G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare); |
---|
| 161 | dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1); |
---|
| 162 | } |
---|
| 163 | |
---|
| 164 | cross += dcross*logstep*alphaprime; |
---|
| 165 | } |
---|
| 166 | |
---|
| 167 | cross *= twopi_mc2_rcl2/beta2; |
---|
| 168 | |
---|
| 169 | } |
---|
| 170 | |
---|
| 171 | // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax |
---|
| 172 | // << " cross= " << cross << G4endl; |
---|
| 173 | |
---|
| 174 | return cross; |
---|
| 175 | } |
---|
| 176 | |
---|
| 177 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 178 | |
---|
| 179 | G4double G4MuBetheBlochModel::ComputeCrossSectionPerAtom( |
---|
| 180 | const G4ParticleDefinition* p, |
---|
| 181 | G4double kineticEnergy, |
---|
| 182 | G4double Z, G4double, |
---|
| 183 | G4double cutEnergy, |
---|
| 184 | G4double maxEnergy) |
---|
| 185 | { |
---|
| 186 | G4double cross = Z*ComputeCrossSectionPerElectron |
---|
| 187 | (p,kineticEnergy,cutEnergy,maxEnergy); |
---|
| 188 | return cross; |
---|
| 189 | } |
---|
| 190 | |
---|
| 191 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 192 | |
---|
| 193 | G4double G4MuBetheBlochModel::CrossSectionPerVolume( |
---|
| 194 | const G4Material* material, |
---|
| 195 | const G4ParticleDefinition* p, |
---|
| 196 | G4double kineticEnergy, |
---|
| 197 | G4double cutEnergy, |
---|
| 198 | G4double maxEnergy) |
---|
| 199 | { |
---|
| 200 | G4double eDensity = material->GetElectronDensity(); |
---|
| 201 | G4double cross = eDensity*ComputeCrossSectionPerElectron |
---|
| 202 | (p,kineticEnergy,cutEnergy,maxEnergy); |
---|
| 203 | return cross; |
---|
| 204 | } |
---|
| 205 | |
---|
| 206 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 207 | |
---|
| 208 | G4double G4MuBetheBlochModel::ComputeDEDXPerVolume(const G4Material* material, |
---|
| 209 | const G4ParticleDefinition* p, |
---|
| 210 | G4double kineticEnergy, |
---|
| 211 | G4double cut) |
---|
| 212 | { |
---|
| 213 | G4double tmax = MaxSecondaryEnergy(p, kineticEnergy); |
---|
| 214 | G4double tau = kineticEnergy/mass; |
---|
| 215 | G4double cutEnergy = min(cut,tmax); |
---|
| 216 | G4double gam = tau + 1.0; |
---|
| 217 | G4double bg2 = tau * (tau+2.0); |
---|
| 218 | G4double beta2 = bg2/(gam*gam); |
---|
| 219 | |
---|
| 220 | G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); |
---|
| 221 | G4double eexc2 = eexc*eexc; |
---|
[1196] | 222 | //G4double cden = material->GetIonisation()->GetCdensity(); |
---|
| 223 | //G4double mden = material->GetIonisation()->GetMdensity(); |
---|
| 224 | //G4double aden = material->GetIonisation()->GetAdensity(); |
---|
| 225 | //G4double x0den = material->GetIonisation()->GetX0density(); |
---|
| 226 | //G4double x1den = material->GetIonisation()->GetX1density(); |
---|
[819] | 227 | |
---|
| 228 | G4double eDensity = material->GetElectronDensity(); |
---|
| 229 | |
---|
| 230 | G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2) |
---|
| 231 | -(1.0 + cutEnergy/tmax)*beta2; |
---|
| 232 | |
---|
| 233 | G4double totEnergy = kineticEnergy + mass; |
---|
| 234 | G4double del = 0.5*cutEnergy/totEnergy; |
---|
| 235 | dedx += del*del; |
---|
| 236 | |
---|
| 237 | // density correction |
---|
| 238 | G4double x = log(bg2)/twoln10; |
---|
[1196] | 239 | //if ( x >= x0den ) { |
---|
| 240 | // dedx -= twoln10*x - cden ; |
---|
| 241 | // if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ; |
---|
| 242 | //} |
---|
| 243 | dedx -= material->GetIonisation()->DensityCorrection(x); |
---|
[819] | 244 | |
---|
| 245 | // shell correction |
---|
| 246 | dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy); |
---|
| 247 | |
---|
| 248 | // now compute the total ionization loss |
---|
| 249 | |
---|
| 250 | if (dedx < 0.0) dedx = 0.0 ; |
---|
| 251 | |
---|
| 252 | // radiative corrections of R. Kokoulin |
---|
| 253 | if (cutEnergy > limitKinEnergy) { |
---|
| 254 | |
---|
| 255 | G4double logtmax = log(cutEnergy); |
---|
| 256 | G4double logstep = logtmax - logLimitKinEnergy; |
---|
| 257 | G4double dloss = 0.0; |
---|
| 258 | G4double ftot2= 0.5/(totEnergy*totEnergy); |
---|
| 259 | |
---|
| 260 | for (G4int ll=0; ll<8; ll++) |
---|
| 261 | { |
---|
| 262 | G4double ep = exp(logLimitKinEnergy + xgi[ll]*logstep); |
---|
| 263 | G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2); |
---|
| 264 | G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare); |
---|
| 265 | dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1); |
---|
| 266 | } |
---|
| 267 | dedx += dloss*logstep*alphaprime; |
---|
| 268 | } |
---|
| 269 | |
---|
| 270 | dedx *= twopi_mc2_rcl2*eDensity/beta2; |
---|
| 271 | |
---|
| 272 | //High order corrections |
---|
[961] | 273 | dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy); |
---|
[819] | 274 | |
---|
| 275 | return dedx; |
---|
| 276 | } |
---|
| 277 | |
---|
| 278 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 279 | |
---|
| 280 | void G4MuBetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp, |
---|
| 281 | const G4MaterialCutsCouple*, |
---|
| 282 | const G4DynamicParticle* dp, |
---|
| 283 | G4double minKinEnergy, |
---|
| 284 | G4double maxEnergy) |
---|
| 285 | { |
---|
| 286 | G4double tmax = MaxSecondaryKinEnergy(dp); |
---|
| 287 | G4double maxKinEnergy = min(maxEnergy,tmax); |
---|
| 288 | if(minKinEnergy >= maxKinEnergy) return; |
---|
| 289 | |
---|
| 290 | G4double kineticEnergy = dp->GetKineticEnergy(); |
---|
| 291 | G4double totEnergy = kineticEnergy + mass; |
---|
| 292 | G4double etot2 = totEnergy*totEnergy; |
---|
| 293 | G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2; |
---|
| 294 | |
---|
| 295 | G4double grej = 1.; |
---|
| 296 | if(tmax > limitKinEnergy) { |
---|
| 297 | G4double a0 = log(2.*totEnergy/mass); |
---|
| 298 | grej += alphaprime*a0*a0; |
---|
| 299 | } |
---|
| 300 | |
---|
| 301 | G4double deltaKinEnergy, f; |
---|
| 302 | |
---|
| 303 | // sampling follows ... |
---|
| 304 | do { |
---|
| 305 | G4double q = G4UniformRand(); |
---|
| 306 | deltaKinEnergy = minKinEnergy*maxKinEnergy |
---|
| 307 | /(minKinEnergy*(1.0 - q) + maxKinEnergy*q); |
---|
| 308 | |
---|
| 309 | |
---|
| 310 | f = 1.0 - beta2*deltaKinEnergy/tmax |
---|
| 311 | + 0.5*deltaKinEnergy*deltaKinEnergy/etot2; |
---|
| 312 | |
---|
| 313 | if(deltaKinEnergy > limitKinEnergy) { |
---|
| 314 | G4double a1 = log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2); |
---|
| 315 | G4double a3 = log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare); |
---|
| 316 | f *= (1. + alphaprime*a1*(a3 - a1)); |
---|
| 317 | } |
---|
| 318 | |
---|
| 319 | if(f > grej) { |
---|
| 320 | G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! " |
---|
| 321 | << "Majorant " << grej << " < " |
---|
| 322 | << f << " for edelta= " << deltaKinEnergy |
---|
| 323 | << " tmin= " << minKinEnergy << " max= " << maxKinEnergy |
---|
| 324 | << G4endl; |
---|
| 325 | } |
---|
| 326 | |
---|
| 327 | |
---|
| 328 | } while( grej*G4UniformRand() > f ); |
---|
| 329 | |
---|
| 330 | G4double deltaMomentum = |
---|
| 331 | sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2)); |
---|
| 332 | G4double totalMomentum = totEnergy*sqrt(beta2); |
---|
| 333 | G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) / |
---|
| 334 | (deltaMomentum * totalMomentum); |
---|
| 335 | |
---|
| 336 | G4double sint = sqrt(1.0 - cost*cost); |
---|
| 337 | |
---|
| 338 | G4double phi = twopi * G4UniformRand() ; |
---|
| 339 | |
---|
| 340 | G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ; |
---|
| 341 | G4ThreeVector direction = dp->GetMomentumDirection(); |
---|
| 342 | deltaDirection.rotateUz(direction); |
---|
| 343 | |
---|
| 344 | // primary change |
---|
| 345 | kineticEnergy -= deltaKinEnergy; |
---|
| 346 | G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection; |
---|
| 347 | direction = dir.unit(); |
---|
| 348 | fParticleChange->SetProposedKineticEnergy(kineticEnergy); |
---|
| 349 | fParticleChange->SetProposedMomentumDirection(direction); |
---|
| 350 | |
---|
| 351 | // create G4DynamicParticle object for delta ray |
---|
| 352 | G4DynamicParticle* delta = new G4DynamicParticle(theElectron, |
---|
| 353 | deltaDirection,deltaKinEnergy); |
---|
| 354 | vdp->push_back(delta); |
---|
| 355 | } |
---|
| 356 | |
---|
| 357 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|