Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // >> 26 // $Id: G4UniversalFluctuation.cc,v 1.13 2007/03/21 15:23:45 urban Exp $ >> 27 // GEANT4 tag $Name: geant4-08-03 $ 26 // 28 // 27 // ------------------------------------------- 29 // ------------------------------------------------------------------- 28 // 30 // 29 // GEANT4 Class file 31 // GEANT4 Class file 30 // 32 // 31 // 33 // 32 // File name: G4UniversalFluctuation 34 // File name: G4UniversalFluctuation 33 // 35 // 34 // Author: V. Ivanchenko for Laszlo Urb << 36 // Author: Vladimir Ivanchenko 35 // 37 // 36 // Creation date: 03.01.2002 38 // Creation date: 03.01.2002 37 // 39 // 38 // Modifications: 40 // Modifications: 39 // 41 // 40 // << 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 // 21-03-07 optimization in ionization part (L.Urban) >> 58 // 41 59 42 //....oooOO0OOooo........oooOO0OOooo........oo 60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 43 //....oooOO0OOooo........oooOO0OOooo........oo 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 44 62 45 #include "G4UniversalFluctuation.hh" 63 #include "G4UniversalFluctuation.hh" 46 #include "G4PhysicalConstants.hh" << 47 #include "G4SystemOfUnits.hh" << 48 #include "Randomize.hh" 64 #include "Randomize.hh" 49 #include "G4Poisson.hh" 65 #include "G4Poisson.hh" >> 66 #include "G4Step.hh" 50 #include "G4Material.hh" 67 #include "G4Material.hh" 51 #include "G4MaterialCutsCouple.hh" << 52 #include "G4DynamicParticle.hh" 68 #include "G4DynamicParticle.hh" 53 #include "G4ParticleDefinition.hh" 69 #include "G4ParticleDefinition.hh" 54 #include "G4Log.hh" << 55 70 56 //....oooOO0OOooo........oooOO0OOooo........oo 71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 57 72 >> 73 using namespace std; >> 74 58 G4UniversalFluctuation::G4UniversalFluctuation 75 G4UniversalFluctuation::G4UniversalFluctuation(const G4String& nam) 59 :G4VEmFluctuationModel(nam), 76 :G4VEmFluctuationModel(nam), 60 minLoss(10.*CLHEP::eV) << 77 particle(0), >> 78 minNumberInteractionsBohr(10.0), >> 79 theBohrBeta2(50.0*keV/proton_mass_c2), >> 80 minLoss(10.*eV), >> 81 nmaxCont1(4.), >> 82 nmaxCont2(16.) 61 { 83 { 62 rndmarray = new G4double[sizearray]; << 84 lastMaterial = 0; 63 } 85 } 64 86 65 //....oooOO0OOooo........oooOO0OOooo........oo 87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 66 88 67 G4UniversalFluctuation::~G4UniversalFluctuatio 89 G4UniversalFluctuation::~G4UniversalFluctuation() 68 { << 90 {} 69 delete [] rndmarray; << 70 } << 71 91 72 //....oooOO0OOooo........oooOO0OOooo........oo 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 73 93 74 void G4UniversalFluctuation::InitialiseMe(cons 94 void G4UniversalFluctuation::InitialiseMe(const G4ParticleDefinition* part) 75 { 95 { 76 particle = part; << 96 particle = part; 77 particleMass = part->GetPDGMass(); << 97 particleMass = part->GetPDGMass(); 78 const G4double q = part->GetPDGCharge()/CLHE << 98 G4double q = part->GetPDGCharge()/eplus; 79 << 99 chargeSquare = q*q; 80 // Derived quantities << 81 m_Inv_particleMass = 1.0 / particleMass; << 82 m_massrate = CLHEP::electron_mass_c2 * m_Inv << 83 chargeSquare = q*q; << 84 } 100 } 85 101 86 //....oooOO0OOooo........oooOO0OOooo........oo 102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 87 103 88 G4double << 104 G4double G4UniversalFluctuation::SampleFluctuations(const G4Material* material, 89 G4UniversalFluctuation::SampleFluctuations(con << 105 const G4DynamicParticle* dp, 90 con << 106 G4double& tmax, 91 con << 107 G4double& length, 92 con << 108 G4double& meanLoss) 93 con << 94 con << 95 { 109 { 96 // Calculate actual loss from the mean loss. << 110 // Calculate actual loss from the mean loss. 97 // The model used to get the fluctuations is << 111 // The model used to get the fluctuations is essentially the same 98 // as in Glandz in Geant3 (Cern program libr << 112 // as in Glandz in Geant3 (Cern program library W5013, phys332). 99 // L. Urban et al. NIM A362, p.416 (1995) an << 113 // L. Urban et al. NIM A362, p.416 (1995) and Geant4 Physics Reference Manual 100 114 101 // shortcut for very small loss or from a st << 115 // shortcut for very very small loss (out of validity of the model) 102 // (out of validity of the model) << 103 // 116 // 104 if (averageLoss < minLoss) { return averageL << 117 if (meanLoss < minLoss) return meanLoss; 105 meanLoss = averageLoss; << 106 const G4double tkin = dp->GetKineticEnergy( << 107 //G4cout<< "Emean= "<< meanLoss<< " tmax= "< << 108 << 109 if(dp->GetDefinition() != particle) { Initia << 110 << 111 CLHEP::HepRandomEngine* rndmEngineF = G4Rand << 112 << 113 const G4double gam = tkin * m_Inv_particle << 114 const G4double gam2 = gam*gam; << 115 const G4double beta = dp->GetBeta(); << 116 const G4double beta2 = beta*beta; << 117 118 118 G4double loss(0.), siga(0.); << 119 if(!particle) InitialiseMe(dp->GetDefinition()); 119 120 120 const G4Material* material = couple->GetMate << 121 G4double tau = dp->GetKineticEnergy()/particleMass; >> 122 G4double gam = tau + 1.0; >> 123 G4double gam2 = gam*gam; >> 124 G4double beta2 = tau*(tau + 2.0)/gam2; >> 125 >> 126 G4double loss(0.), siga(0.); 121 127 122 // Gaussian regime 128 // Gaussian regime 123 // for heavy particles only and conditions 129 // for heavy particles only and conditions 124 // for Gauusian fluct. has been changed 130 // for Gauusian fluct. has been changed 125 // 131 // 126 if (particleMass > CLHEP::electron_mass_c2 & << 132 if ((particleMass > electron_mass_c2) && 127 meanLoss >= minNumberInteractionsBohr*tc << 133 (meanLoss >= minNumberInteractionsBohr*tmax)) 128 << 134 { 129 siga = std::sqrt((tmax/beta2 - 0.5*tcut)*C << 135 G4double massrate = electron_mass_c2/particleMass ; 130 length*chargeSquare*mate << 136 G4double tmaxkine = 2.*electron_mass_c2*beta2*gam2/ 131 const G4double sn = meanLoss/siga; << 137 (1.+massrate*(2.*gam+massrate)) ; 132 << 138 if (tmaxkine <= 2.*tmax) 133 // thick target case << 139 { 134 if (sn >= 2.0) { << 140 electronDensity = material->GetElectronDensity(); 135 << 141 siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length 136 const G4double twomeanLoss = meanLoss + << 142 * electronDensity * chargeSquare; 137 do { << 143 siga = sqrt(siga); 138 loss = G4RandGauss::shoot(rndmEngineF, meanL << 144 G4double twomeanLoss = meanLoss + meanLoss; 139 // Loop checking, 03-Aug-2015, Vladimir Ivan << 145 if (twomeanLoss < siga) { 140 } while (0.0 > loss || twomeanLoss < lo << 146 G4double x; 141 << 147 do { 142 // Gamma distribution << 148 loss = twomeanLoss*G4UniformRand(); 143 } else { << 149 x = (loss - meanLoss)/siga; 144 << 150 } while (1.0 - 0.5*x*x < G4UniformRand()); 145 const G4double neff = sn*sn; << 151 } else { 146 loss = meanLoss*G4RandGamma::shoot(rndmE << 152 do { >> 153 loss = G4RandGauss::shoot(meanLoss,siga); >> 154 } while (loss < 0. || loss > twomeanLoss); >> 155 } >> 156 return loss; 147 } 157 } 148 //G4cout << "Gauss: " << loss << G4endl; << 149 return loss; << 150 } 158 } 151 159 152 auto ioni = material->GetIonisation(); << 160 // Glandz regime : initialisation 153 e0 = ioni->GetEnergy0fluct(); << 161 // 154 << 162 if (material != lastMaterial) { 155 // very small step or low-density material << 163 f1Fluct = material->GetIonisation()->GetF1fluct(); 156 if(tcut <= e0) { return meanLoss; } << 164 f2Fluct = material->GetIonisation()->GetF2fluct(); 157 << 165 e1Fluct = material->GetIonisation()->GetEnergy1fluct(); 158 ipotFluct = ioni->GetMeanExcitationEnergy(); << 166 e2Fluct = material->GetIonisation()->GetEnergy2fluct(); 159 ipotLogFluct = ioni->GetLogMeanExcEnergy(); << 167 e1LogFluct = material->GetIonisation()->GetLogEnergy1fluct(); >> 168 e2LogFluct = material->GetIonisation()->GetLogEnergy2fluct(); >> 169 ipotFluct = material->GetIonisation()->GetMeanExcitationEnergy(); >> 170 ipotLogFluct = material->GetIonisation()->GetLogMeanExcEnergy(); >> 171 e0 = material->GetIonisation()->GetEnergy0fluct(); >> 172 lastMaterial = material; >> 173 } 160 174 161 // width correction for small cuts << 175 G4double a1 = 0. , a2 = 0., a3 = 0. ; 162 const G4double scaling = std::min(1.+0.5*CLH << 163 meanLoss /= scaling; << 164 176 165 w2 = (tcut > ipotFluct) ? << 177 // cut and material dependent rate 166 G4Log(2.*CLHEP::electron_mass_c2*beta2*gam << 178 G4double rate = 1.0; 167 return SampleGlandz(rndmEngineF, material, t << 179 if(tmax > ipotFluct) { 168 } << 180 G4double w2 = log(2.*electron_mass_c2*beta2*gam2)-beta2; >> 181 >> 182 if(w2 > ipotLogFluct && w2 > e2LogFluct) { >> 183 >> 184 rate = 0.03+0.23*log(log(tmax/ipotFluct)); >> 185 G4double C = meanLoss*(1.-rate)/(w2-ipotLogFluct); >> 186 a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct; >> 187 a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct; >> 188 } >> 189 } 169 190 170 //....oooOO0OOooo........oooOO0OOooo........oo << 191 G4double w1 = tmax/e0; >> 192 if(tmax > e0) >> 193 a3 = rate*meanLoss*(tmax-e0)/(e0*tmax*log(w1)); 171 194 172 G4double << 195 //'nearly' Gaussian fluctuation if a1>nmaxCont2&&a2>nmaxCont2&&a3>nmaxCont2 173 G4UniversalFluctuation::SampleGlandz(CLHEP::He << 174 const G4M << 175 const G4d << 176 { << 177 G4double a1(0.0), a3(0.0); << 178 G4double loss = 0.0; << 179 G4double e1 = ipotFluct; << 180 << 181 if(tcut > e1) { << 182 a1 = meanLoss*(1.-rate)/e1; << 183 if(a1 < a0) { << 184 const G4double fwnow = 0.1+(fw-0.1)*std: << 185 a1 /= fwnow; << 186 e1 *= fwnow; << 187 } else { << 188 a1 /= fw; << 189 e1 *= fw; << 190 } << 191 } << 192 << 193 const G4double w1 = tcut/e0; << 194 a3 = rate*meanLoss*(tcut - e0)/(e0*tcut*G4Lo << 195 if(a1 <= 0.) { a3 /= rate; } << 196 << 197 //'nearly' Gaussian fluctuation if a1>nmaxCo << 198 G4double emean = 0.; 196 G4double emean = 0.; 199 G4double sig2e = 0.; << 197 G4double sig2e = 0., sige = 0.; 200 << 198 G4double p1 = 0., p2 = 0., p3 = 0.; >> 199 201 // excitation of type 1 200 // excitation of type 1 202 if(a1 > 0.0) { AddExcitation(rndmEngineF, a1 << 201 if(a1 > nmaxCont2) >> 202 { >> 203 emean += a1*e1Fluct; >> 204 sig2e += a1*e1Fluct*e1Fluct; >> 205 } >> 206 else if(a1 > 0.) >> 207 { >> 208 p1 = G4double(G4Poisson(a1)); >> 209 loss += p1*e1Fluct; >> 210 if(p1 > 0.) >> 211 loss += (1.-2.*G4UniformRand())*e1Fluct; >> 212 } 203 213 204 if(sig2e > 0.0) { SampleGauss(rndmEngineF, e << 214 // excitation of type 2 >> 215 if(a2 > nmaxCont2) >> 216 { >> 217 emean += a2*e2Fluct; >> 218 sig2e += a2*e2Fluct*e2Fluct; >> 219 } >> 220 else if(a2 > 0.) >> 221 { >> 222 p2 = G4double(G4Poisson(a2)); >> 223 loss += p2*e2Fluct; >> 224 if(p2 > 0.) >> 225 loss += (1.-2.*G4UniformRand())*e2Fluct; >> 226 } 205 227 206 // ionisation 228 // ionisation 207 if(a3 > 0.) { << 229 G4double lossc = 0.; 208 emean = 0.; << 230 if(a3 > 0.) 209 sig2e = 0.; << 231 { 210 G4double p3 = a3; << 232 p3 = a3; 211 G4double alfa = 1.; 233 G4double alfa = 1.; 212 if(a3 > nmaxCont) { << 234 if(a3 > nmaxCont2) 213 alfa = w1*(nmaxCont+a3)/(w1*nmaxCont+a3) << 235 { 214 const G4double alfa1 = alfa*G4Log(alfa) << 236 alfa = w1*(nmaxCont2+a3)/(w1*nmaxCont2+a3); 215 const G4double namean = a3*w1*(alfa-1.)/ << 237 G4double alfa1 = alfa*log(alfa)/(alfa-1.); 216 emean += namean*e0*alfa1; << 238 G4double namean = a3*w1*(alfa-1.)/((w1-1.)*alfa); 217 sig2e += e0*e0*namean*(alfa-alfa1*alfa1) << 239 emean += namean*e0*alfa1; 218 p3 = a3 - namean; << 240 sig2e += e0*e0*namean*(alfa-alfa1*alfa1); >> 241 p3 = a3-namean; 219 } 242 } 220 << 243 221 const G4double w3 = alfa*e0; << 244 G4double w2 = alfa*e0; 222 if(tcut > w3) { << 245 G4double w = (tmax-w2)/tmax; 223 const G4double w = (tcut-w3)/tcut; << 246 G4double scale = 1.; 224 const G4int nnb = (G4int)G4Poisson(p3); << 247 G4int nb = 0; 225 if(nnb > 0) { << 248 if(p3 < nmaxCont2) 226 if(nnb > sizearray) { << 249 nb = G4Poisson(p3); 227 sizearray = nnb; << 250 else 228 delete [] rndmarray; << 251 { 229 rndmarray = new G4double[nnb]; << 252 nb = G4Poisson(nmaxCont2); 230 } << 253 scale = p3/nmaxCont2; 231 rndmEngineF->flatArray(nnb, rndmarray) << 232 for (G4int k=0; k<nnb; ++k) { loss += << 233 } << 234 } 254 } 235 if(sig2e > 0.0) { SampleGauss(rndmEngineF, << 255 if(nb > 0) >> 256 for (G4int k=0; k<nb; k++) lossc += scale*w2/(1.-w*G4UniformRand()); 236 } 257 } 237 //G4cout << "### loss=" << loss << G4endl; << 258 >> 259 if(emean > 0.) >> 260 { >> 261 sige = sqrt(sig2e); >> 262 loss += max(0.,G4RandGauss::shoot(emean,sige)); >> 263 } >> 264 >> 265 loss += lossc; >> 266 238 return loss; 267 return loss; >> 268 239 } 269 } 240 270 241 //....oooOO0OOooo........oooOO0OOooo........oo 271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 242 272 243 273 244 G4double G4UniversalFluctuation::Dispersion( 274 G4double G4UniversalFluctuation::Dispersion( 245 const G4Material* ma 275 const G4Material* material, 246 const G4DynamicParti 276 const G4DynamicParticle* dp, 247 const G4double tcut, << 277 G4double& tmax, 248 const G4double tmax, << 278 G4double& length) 249 const G4double lengt << 250 { 279 { 251 if(dp->GetDefinition() != particle) { Initia << 280 if(!particle) InitialiseMe(dp->GetDefinition()); 252 const G4double beta = dp->GetBeta(); << 253 return (tmax/(beta*beta) - 0.5*tcut) * CLHEP << 254 * material->GetElectronDensity() * chargeS << 255 } << 256 281 257 //....oooOO0OOooo........oooOO0OOooo........oo << 282 electronDensity = material->GetElectronDensity(); 258 283 259 void << 284 G4double gam = (dp->GetKineticEnergy())/particleMass + 1.0; 260 G4UniversalFluctuation::SetParticleAndCharge(c << 285 G4double beta2 = 1.0 - 1.0/(gam*gam); 261 G << 286 262 { << 287 G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length 263 if(part != particle) { << 288 * electronDensity * chargeSquare; 264 particle = part; << 289 265 particleMass = part->GetPDGMass(); << 290 return siga; 266 << 267 // Derived quantities << 268 m_Inv_particleMass = 1.0 / particleMass; << 269 m_massrate = CLHEP::electron_mass_c2 * m_I << 270 } << 271 chargeSquare = q2; << 272 } 291 } 273 292 274 //....oooOO0OOooo........oooOO0OOooo........oo 293 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 275 294