| [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 | //
|
|---|
| [1055] | 26 | // $Id: G4UniversalFluctuation.cc,v 1.22 2009/03/20 18:11:23 urban Exp $
|
|---|
| [1196] | 27 | // GEANT4 tag $Name: geant4-09-03-cand-01 $
|
|---|
| [819] | 28 | //
|
|---|
| 29 | // -------------------------------------------------------------------
|
|---|
| 30 | //
|
|---|
| 31 | // GEANT4 Class file
|
|---|
| 32 | //
|
|---|
| 33 | //
|
|---|
| 34 | // File name: G4UniversalFluctuation
|
|---|
| 35 | //
|
|---|
| 36 | // Author: Vladimir Ivanchenko
|
|---|
| 37 | //
|
|---|
| 38 | // Creation date: 03.01.2002
|
|---|
| 39 | //
|
|---|
| 40 | // Modifications:
|
|---|
| 41 | //
|
|---|
| 42 | // 28-12-02 add method Dispersion (V.Ivanchenko)
|
|---|
| 43 | // 07-02-03 change signature (V.Ivanchenko)
|
|---|
| 44 | // 13-02-03 Add name (V.Ivanchenko)
|
|---|
| 45 | // 16-10-03 Changed interface to Initialisation (V.Ivanchenko)
|
|---|
| 46 | // 07-11-03 Fix problem of rounding of double in G4UniversalFluctuations
|
|---|
| 47 | // 06-02-04 Add control on big sigma > 2*meanLoss (V.Ivanchenko)
|
|---|
| 48 | // 26-04-04 Comment out the case of very small step (V.Ivanchenko)
|
|---|
| 49 | // 07-02-05 define problim = 5.e-3 (mma)
|
|---|
| 50 | // 03-05-05 conditions of Gaussian fluctuation changed (bugfix)
|
|---|
| 51 | // + smearing for very small loss (L.Urban)
|
|---|
| 52 | // 03-10-05 energy dependent rate -> cut dependence of the
|
|---|
| 53 | // distribution is much weaker (L.Urban)
|
|---|
| 54 | // 17-10-05 correction for very small loss (L.Urban)
|
|---|
| 55 | // 20-03-07 'GLANDZ' part rewritten completely, no 'very small loss'
|
|---|
| 56 | // regime any more (L.Urban)
|
|---|
| 57 | // 03-04-07 correction to get better width of eloss distr.(L.Urban)
|
|---|
| 58 | // 13-07-07 add protection for very small step or low-density material (VI)
|
|---|
| [1055] | 59 | // 19-03-09 new width correction (does not depend on previous steps) (L.Urban)
|
|---|
| 60 | // 20-03-09 modification in the width correction (L.Urban)
|
|---|
| 61 | //
|
|---|
| [819] | 62 |
|
|---|
| 63 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 64 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 65 |
|
|---|
| 66 | #include "G4UniversalFluctuation.hh"
|
|---|
| 67 | #include "Randomize.hh"
|
|---|
| 68 | #include "G4Poisson.hh"
|
|---|
| 69 | #include "G4Step.hh"
|
|---|
| 70 | #include "G4Material.hh"
|
|---|
| 71 | #include "G4DynamicParticle.hh"
|
|---|
| 72 | #include "G4ParticleDefinition.hh"
|
|---|
| 73 |
|
|---|
| 74 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 75 |
|
|---|
| 76 | using namespace std;
|
|---|
| 77 |
|
|---|
| 78 | G4UniversalFluctuation::G4UniversalFluctuation(const G4String& nam)
|
|---|
| 79 | :G4VEmFluctuationModel(nam),
|
|---|
| 80 | particle(0),
|
|---|
| 81 | minNumberInteractionsBohr(10.0),
|
|---|
| 82 | theBohrBeta2(50.0*keV/proton_mass_c2),
|
|---|
| 83 | minLoss(10.*eV),
|
|---|
| 84 | nmaxCont1(4.),
|
|---|
| 85 | nmaxCont2(16.)
|
|---|
| 86 | {
|
|---|
| 87 | lastMaterial = 0;
|
|---|
| 88 | }
|
|---|
| 89 |
|
|---|
| 90 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 91 |
|
|---|
| 92 | G4UniversalFluctuation::~G4UniversalFluctuation()
|
|---|
| 93 | {}
|
|---|
| 94 |
|
|---|
| 95 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 96 |
|
|---|
| 97 | void G4UniversalFluctuation::InitialiseMe(const G4ParticleDefinition* part)
|
|---|
| 98 | {
|
|---|
| 99 | particle = part;
|
|---|
| 100 | particleMass = part->GetPDGMass();
|
|---|
| 101 | G4double q = part->GetPDGCharge()/eplus;
|
|---|
| 102 | chargeSquare = q*q;
|
|---|
| 103 | }
|
|---|
| 104 |
|
|---|
| 105 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 106 |
|
|---|
| 107 | G4double G4UniversalFluctuation::SampleFluctuations(const G4Material* material,
|
|---|
| [961] | 108 | const G4DynamicParticle* dp,
|
|---|
| 109 | G4double& tmax,
|
|---|
| 110 | G4double& length,
|
|---|
| 111 | G4double& meanLoss)
|
|---|
| [819] | 112 | {
|
|---|
| 113 | // Calculate actual loss from the mean loss.
|
|---|
| 114 | // The model used to get the fluctuations is essentially the same
|
|---|
| 115 | // as in Glandz in Geant3 (Cern program library W5013, phys332).
|
|---|
| 116 | // L. Urban et al. NIM A362, p.416 (1995) and Geant4 Physics Reference Manual
|
|---|
| 117 |
|
|---|
| 118 | // shortcut for very very small loss (out of validity of the model)
|
|---|
| 119 | //
|
|---|
| 120 | if (meanLoss < minLoss)
|
|---|
| 121 | return meanLoss;
|
|---|
| 122 |
|
|---|
| 123 | if(!particle) InitialiseMe(dp->GetDefinition());
|
|---|
| 124 |
|
|---|
| 125 | G4double tau = dp->GetKineticEnergy()/particleMass;
|
|---|
| 126 | G4double gam = tau + 1.0;
|
|---|
| 127 | G4double gam2 = gam*gam;
|
|---|
| 128 | G4double beta2 = tau*(tau + 2.0)/gam2;
|
|---|
| 129 |
|
|---|
| 130 | G4double loss(0.), siga(0.);
|
|---|
| 131 |
|
|---|
| 132 | // Gaussian regime
|
|---|
| 133 | // for heavy particles only and conditions
|
|---|
| 134 | // for Gauusian fluct. has been changed
|
|---|
| 135 | //
|
|---|
| 136 | if ((particleMass > electron_mass_c2) &&
|
|---|
| 137 | (meanLoss >= minNumberInteractionsBohr*tmax))
|
|---|
| 138 | {
|
|---|
| 139 | G4double massrate = electron_mass_c2/particleMass ;
|
|---|
| 140 | G4double tmaxkine = 2.*electron_mass_c2*beta2*gam2/
|
|---|
| 141 | (1.+massrate*(2.*gam+massrate)) ;
|
|---|
| 142 | if (tmaxkine <= 2.*tmax)
|
|---|
| 143 | {
|
|---|
| 144 | electronDensity = material->GetElectronDensity();
|
|---|
| 145 | siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
|
|---|
| 146 | * electronDensity * chargeSquare;
|
|---|
| 147 | siga = sqrt(siga);
|
|---|
| 148 | G4double twomeanLoss = meanLoss + meanLoss;
|
|---|
| 149 | if (twomeanLoss < siga) {
|
|---|
| 150 | G4double x;
|
|---|
| 151 | do {
|
|---|
| 152 | loss = twomeanLoss*G4UniformRand();
|
|---|
| 153 | x = (loss - meanLoss)/siga;
|
|---|
| 154 | } while (1.0 - 0.5*x*x < G4UniformRand());
|
|---|
| 155 | } else {
|
|---|
| 156 | do {
|
|---|
| 157 | loss = G4RandGauss::shoot(meanLoss,siga);
|
|---|
| 158 | } while (loss < 0. || loss > twomeanLoss);
|
|---|
| 159 | }
|
|---|
| 160 | return loss;
|
|---|
| 161 | }
|
|---|
| 162 | }
|
|---|
| 163 |
|
|---|
| 164 | // Glandz regime : initialisation
|
|---|
| 165 | //
|
|---|
| 166 | if (material != lastMaterial) {
|
|---|
| 167 | f1Fluct = material->GetIonisation()->GetF1fluct();
|
|---|
| 168 | f2Fluct = material->GetIonisation()->GetF2fluct();
|
|---|
| 169 | e1Fluct = material->GetIonisation()->GetEnergy1fluct();
|
|---|
| 170 | e2Fluct = material->GetIonisation()->GetEnergy2fluct();
|
|---|
| 171 | e1LogFluct = material->GetIonisation()->GetLogEnergy1fluct();
|
|---|
| 172 | e2LogFluct = material->GetIonisation()->GetLogEnergy2fluct();
|
|---|
| 173 | ipotFluct = material->GetIonisation()->GetMeanExcitationEnergy();
|
|---|
| 174 | ipotLogFluct = material->GetIonisation()->GetLogMeanExcEnergy();
|
|---|
| 175 | e0 = material->GetIonisation()->GetEnergy0fluct();
|
|---|
| 176 | lastMaterial = material;
|
|---|
| [1055] | 177 |
|
|---|
| 178 | // modification of some model parameters
|
|---|
| 179 | // (this part should go to materials later)
|
|---|
| 180 | G4double p = 1.40;
|
|---|
| 181 | f2Fluct *= p;
|
|---|
| 182 | f1Fluct = 1.-f2Fluct;
|
|---|
| 183 | G4double q = 1.00;
|
|---|
| 184 | e2Fluct *= q;
|
|---|
| 185 | e2LogFluct = log(e2Fluct);
|
|---|
| 186 | e1LogFluct = (ipotLogFluct-f2Fluct*e2LogFluct)/f1Fluct;
|
|---|
| 187 | e1Fluct = exp(e1LogFluct);
|
|---|
| [819] | 188 | }
|
|---|
| 189 |
|
|---|
| 190 | // very small step or low-density material
|
|---|
| 191 | if(tmax <= e0) return meanLoss;
|
|---|
| 192 |
|
|---|
| 193 | G4double a1 = 0. , a2 = 0., a3 = 0. ;
|
|---|
| 194 |
|
|---|
| 195 | // cut and material dependent rate
|
|---|
| 196 | G4double rate = 1.0;
|
|---|
| 197 | if(tmax > ipotFluct) {
|
|---|
| 198 | G4double w2 = log(2.*electron_mass_c2*beta2*gam2)-beta2;
|
|---|
| 199 |
|
|---|
| 200 | if(w2 > ipotLogFluct && w2 > e2LogFluct) {
|
|---|
| 201 |
|
|---|
| 202 | rate = 0.03+0.23*log(log(tmax/ipotFluct));
|
|---|
| 203 | G4double C = meanLoss*(1.-rate)/(w2-ipotLogFluct);
|
|---|
| [1055] | 204 | a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
|
|---|
| 205 | a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
|
|---|
| 206 | // correction in order to get better FWHM values
|
|---|
| 207 | // ( scale parameters a1 and e1)
|
|---|
| 208 | G4double width = 1.;
|
|---|
| 209 | if(meanLoss > 10.*e1Fluct)
|
|---|
| 210 | {
|
|---|
| 211 | width = 3.1623/sqrt(meanLoss/e1Fluct);
|
|---|
| 212 | if(width < a2/a1)
|
|---|
| 213 | width = a2/a1;
|
|---|
| 214 | }
|
|---|
| 215 | a1 *= width;
|
|---|
| 216 | e1 = e1Fluct/width;
|
|---|
| 217 | e2 = e2Fluct;
|
|---|
| [819] | 218 | }
|
|---|
| 219 | }
|
|---|
| 220 |
|
|---|
| 221 | G4double w1 = tmax/e0;
|
|---|
| 222 | if(tmax > e0)
|
|---|
| 223 | a3 = rate*meanLoss*(tmax-e0)/(e0*tmax*log(w1));
|
|---|
| 224 |
|
|---|
| 225 | //'nearly' Gaussian fluctuation if a1>nmaxCont2&&a2>nmaxCont2&&a3>nmaxCont2
|
|---|
| 226 | G4double emean = 0.;
|
|---|
| 227 | G4double sig2e = 0., sige = 0.;
|
|---|
| 228 | G4double p1 = 0., p2 = 0., p3 = 0.;
|
|---|
| 229 |
|
|---|
| 230 | // excitation of type 1
|
|---|
| 231 | if(a1 > nmaxCont2)
|
|---|
| 232 | {
|
|---|
| 233 | emean += a1*e1;
|
|---|
| 234 | sig2e += a1*e1*e1;
|
|---|
| 235 | }
|
|---|
| 236 | else if(a1 > 0.)
|
|---|
| 237 | {
|
|---|
| 238 | p1 = G4double(G4Poisson(a1));
|
|---|
| 239 | loss += p1*e1;
|
|---|
| 240 | if(p1 > 0.)
|
|---|
| 241 | loss += (1.-2.*G4UniformRand())*e1;
|
|---|
| 242 | }
|
|---|
| 243 |
|
|---|
| 244 | // excitation of type 2
|
|---|
| 245 | if(a2 > nmaxCont2)
|
|---|
| 246 | {
|
|---|
| 247 | emean += a2*e2;
|
|---|
| 248 | sig2e += a2*e2*e2;
|
|---|
| 249 | }
|
|---|
| 250 | else if(a2 > 0.)
|
|---|
| 251 | {
|
|---|
| 252 | p2 = G4double(G4Poisson(a2));
|
|---|
| 253 | loss += p2*e2;
|
|---|
| 254 | if(p2 > 0.)
|
|---|
| 255 | loss += (1.-2.*G4UniformRand())*e2;
|
|---|
| 256 | }
|
|---|
| 257 |
|
|---|
| 258 | // ionisation
|
|---|
| 259 | G4double lossc = 0.;
|
|---|
| 260 | if(a3 > 0.)
|
|---|
| 261 | {
|
|---|
| 262 | p3 = a3;
|
|---|
| 263 | G4double alfa = 1.;
|
|---|
| 264 | if(a3 > nmaxCont2)
|
|---|
| 265 | {
|
|---|
| 266 | alfa = w1*(nmaxCont2+a3)/(w1*nmaxCont2+a3);
|
|---|
| 267 | G4double alfa1 = alfa*log(alfa)/(alfa-1.);
|
|---|
| 268 | G4double namean = a3*w1*(alfa-1.)/((w1-1.)*alfa);
|
|---|
| 269 | emean += namean*e0*alfa1;
|
|---|
| 270 | sig2e += e0*e0*namean*(alfa-alfa1*alfa1);
|
|---|
| 271 | p3 = a3-namean;
|
|---|
| 272 | }
|
|---|
| 273 |
|
|---|
| 274 | G4double w2 = alfa*e0;
|
|---|
| 275 | G4double w = (tmax-w2)/tmax;
|
|---|
| 276 | G4int nb = G4Poisson(p3);
|
|---|
| 277 | if(nb > 0)
|
|---|
| 278 | for (G4int k=0; k<nb; k++) lossc += w2/(1.-w*G4UniformRand());
|
|---|
| 279 | }
|
|---|
| 280 |
|
|---|
| 281 | if(emean > 0.)
|
|---|
| 282 | {
|
|---|
| 283 | sige = sqrt(sig2e);
|
|---|
| 284 | loss += max(0.,G4RandGauss::shoot(emean,sige));
|
|---|
| 285 | }
|
|---|
| 286 |
|
|---|
| 287 | loss += lossc;
|
|---|
| 288 |
|
|---|
| 289 | return loss;
|
|---|
| 290 |
|
|---|
| 291 | }
|
|---|
| 292 |
|
|---|
| 293 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|
| 294 |
|
|---|
| 295 |
|
|---|
| 296 | G4double G4UniversalFluctuation::Dispersion(
|
|---|
| 297 | const G4Material* material,
|
|---|
| 298 | const G4DynamicParticle* dp,
|
|---|
| 299 | G4double& tmax,
|
|---|
| 300 | G4double& length)
|
|---|
| 301 | {
|
|---|
| 302 | if(!particle) InitialiseMe(dp->GetDefinition());
|
|---|
| 303 |
|
|---|
| 304 | electronDensity = material->GetElectronDensity();
|
|---|
| 305 |
|
|---|
| 306 | G4double gam = (dp->GetKineticEnergy())/particleMass + 1.0;
|
|---|
| 307 | G4double beta2 = 1.0 - 1.0/(gam*gam);
|
|---|
| 308 |
|
|---|
| 309 | G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
|
|---|
| 310 | * electronDensity * chargeSquare;
|
|---|
| 311 |
|
|---|
| 312 | return siga;
|
|---|
| 313 | }
|
|---|
| 314 |
|
|---|
| [1055] | 315 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 316 |
|
|---|
| 317 | void
|
|---|
| 318 | G4UniversalFluctuation::SetParticleAndCharge(const G4ParticleDefinition* part,
|
|---|
| 319 | G4double q2)
|
|---|
| 320 | {
|
|---|
| 321 | if(part != particle) {
|
|---|
| 322 | particle = part;
|
|---|
| 323 | particleMass = part->GetPDGMass();
|
|---|
| 324 | }
|
|---|
| 325 | chargeSquare = q2;
|
|---|
| 326 | }
|
|---|
| 327 |
|
|---|
| [819] | 328 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
|
|---|