source: trunk/source/processes/electromagnetic/standard/src/G4BetheBlochModel.cc@ 1005

Last change on this file since 1005 was 1005, checked in by garnier, 17 years ago

fichier oublies

File size: 15.9 KB
RevLine 
[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//
[1005]26// $Id: G4BetheBlochModel.cc,v 1.25 2009/02/20 12:06:37 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
[819]28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4BetheBlochModel
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 03.01.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// 24-03-05 Add G4EmCorrections (V.Ivanchenko)
[961]47// 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
[819]48// 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49// 12-02-06 move G4LossTableManager::Instance()->EmCorrections()
50// in constructor (mma)
[961]51// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
52// CorrectionsAlongStep needed for ions(V.Ivanchenko)
[819]53//
54// -------------------------------------------------------------------
55//
56
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
61#include "G4BetheBlochModel.hh"
62#include "Randomize.hh"
63#include "G4Electron.hh"
64#include "G4LossTableManager.hh"
65#include "G4EmCorrections.hh"
66#include "G4ParticleChangeForLoss.hh"
[961]67#include "G4NistManager.hh"
[819]68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71using namespace std;
72
73G4BetheBlochModel::G4BetheBlochModel(const G4ParticleDefinition* p,
74 const G4String& nam)
75 : G4VEmModel(nam),
[961]76 particle(0),
77 tlimit(DBL_MAX),
78 twoln10(2.0*log(10.0)),
79 bg2lim(0.0169),
80 taulim(8.4146e-3),
81 isIon(false),
82 isInitialised(false)
[819]83{
[961]84 fParticleChange = 0;
[819]85 if(p) SetParticle(p);
86 theElectron = G4Electron::Electron();
87 corr = G4LossTableManager::Instance()->EmCorrections();
[961]88 nist = G4NistManager::Instance();
89 SetLowEnergyLimit(2.0*MeV);
[819]90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
94G4BetheBlochModel::~G4BetheBlochModel()
95{}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98
99G4double G4BetheBlochModel::MinEnergyCut(const G4ParticleDefinition*,
100 const G4MaterialCutsCouple* couple)
101{
102 return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106
107void G4BetheBlochModel::Initialise(const G4ParticleDefinition* p,
108 const G4DataVector&)
109{
110 if (!particle) SetParticle(p);
111
[961]112 corrFactor = chargeSquare;
[1005]113 // always false before the run
114 SetDeexcitationFlag(false);
[961]115
116 if(!isInitialised) {
117 isInitialised = true;
118
119 if(!fParticleChange) {
120 if (pParticleChange) {
121 fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*>
122 (pParticleChange);
123 } else {
124 fParticleChange = new G4ParticleChangeForLoss();
125 }
126 }
127 }
[819]128}
129
[961]130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131
132void G4BetheBlochModel::SetParticle(const G4ParticleDefinition* p)
133{
134 if(particle != p) {
135 particle = p;
136 G4String pname = particle->GetParticleName();
137 if (particle->GetParticleType() == "nucleus" &&
138 pname != "deuteron" && pname != "triton") {
139 isIon = true;
140 }
141
142 mass = particle->GetPDGMass();
143 spin = particle->GetPDGSpin();
144 G4double q = particle->GetPDGCharge()/eplus;
145 chargeSquare = q*q;
146 ratio = electron_mass_c2/mass;
147 G4double magmom = particle->GetPDGMagneticMoment()
148 *mass/(0.5*eplus*hbar_Planck*c_squared);
149 magMoment2 = magmom*magmom - 1.0;
150 formfact = 0.0;
151 if(particle->GetLeptonNumber() == 0) {
152 G4double x = 0.8426*GeV;
153 if(spin == 0.0 && mass < GeV) {x = 0.736*GeV;}
154 else if(mass > GeV) {
155 x /= nist->GetZ13(mass/proton_mass_c2);
156 // tlimit = 51.2*GeV*A13[iz]*A13[iz];
157 }
158 formfact = 2.0*electron_mass_c2/(x*x);
159 tlimit = 2.0/formfact;
160 }
161 }
162}
163
[819]164//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165
[961]166G4double G4BetheBlochModel::GetChargeSquareRatio(const G4ParticleDefinition* p,
167 const G4Material* mat,
168 G4double kineticEnergy)
[819]169{
[961]170 // this method is called only for ions
171 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
172 corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
173 return corrFactor;
174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177
178G4double G4BetheBlochModel::GetParticleCharge(const G4ParticleDefinition* p,
179 const G4Material* mat,
180 G4double kineticEnergy)
181{
182 // this method is called only for ions
183 return corr->GetParticleCharge(p,mat,kineticEnergy);
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
187
188G4double
189G4BetheBlochModel::ComputeCrossSectionPerElectron(const G4ParticleDefinition* p,
190 G4double kineticEnergy,
191 G4double cutEnergy,
192 G4double maxKinEnergy)
193{
[819]194 G4double cross = 0.0;
195 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
196 G4double maxEnergy = min(tmax,maxKinEnergy);
197 if(cutEnergy < maxEnergy) {
198
199 G4double totEnergy = kineticEnergy + mass;
200 G4double energy2 = totEnergy*totEnergy;
201 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
202
[961]203 cross = 1.0/cutEnergy - 1.0/maxEnergy
204 - beta2*log(maxEnergy/cutEnergy)/tmax;
[819]205
206 // +term for spin=1/2 particle
207 if( 0.5 == spin ) cross += 0.5*(maxEnergy - cutEnergy)/energy2;
208
[961]209 // High order correction different for hadrons and ions
210 // nevetheless they are applied to reduce high energy transfers
211 // if(!isIon)
212 //cross += corr->FiniteSizeCorrectionXS(p,currentMaterial,
213 // kineticEnergy,cutEnergy);
214
[819]215 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
216 }
217
218 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
219 // << " tmax= " << tmax << " cross= " << cross << G4endl;
220
221 return cross;
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
225
226G4double G4BetheBlochModel::ComputeCrossSectionPerAtom(
227 const G4ParticleDefinition* p,
228 G4double kineticEnergy,
229 G4double Z, G4double,
230 G4double cutEnergy,
231 G4double maxEnergy)
232{
233 G4double cross = Z*ComputeCrossSectionPerElectron
234 (p,kineticEnergy,cutEnergy,maxEnergy);
235 return cross;
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
239
240G4double G4BetheBlochModel::CrossSectionPerVolume(
241 const G4Material* material,
242 const G4ParticleDefinition* p,
243 G4double kineticEnergy,
244 G4double cutEnergy,
245 G4double maxEnergy)
246{
[961]247 currentMaterial = material;
[819]248 G4double eDensity = material->GetElectronDensity();
249 G4double cross = eDensity*ComputeCrossSectionPerElectron
250 (p,kineticEnergy,cutEnergy,maxEnergy);
251 return cross;
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255
256G4double G4BetheBlochModel::ComputeDEDXPerVolume(const G4Material* material,
257 const G4ParticleDefinition* p,
258 G4double kineticEnergy,
259 G4double cut)
260{
261 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
262 G4double cutEnergy = min(cut,tmax);
263
264 G4double tau = kineticEnergy/mass;
265 G4double gam = tau + 1.0;
266 G4double bg2 = tau * (tau+2.0);
267 G4double beta2 = bg2/(gam*gam);
268
269 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
270 G4double eexc2 = eexc*eexc;
271 G4double cden = material->GetIonisation()->GetCdensity();
272 G4double mden = material->GetIonisation()->GetMdensity();
273 G4double aden = material->GetIonisation()->GetAdensity();
274 G4double x0den = material->GetIonisation()->GetX0density();
275 G4double x1den = material->GetIonisation()->GetX1density();
276
277 G4double eDensity = material->GetElectronDensity();
278
279 G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
280 - (1.0 + cutEnergy/tmax)*beta2;
281
282 if(0.5 == spin) {
283 G4double del = 0.5*cutEnergy/(kineticEnergy + mass);
284 dedx += del*del;
285 }
286
287 // density correction
288 G4double x = log(bg2)/twoln10;
289 if ( x >= x0den ) {
290 dedx -= twoln10*x - cden ;
291 if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
292 }
293
294 // shell correction
295 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
296
297 // now compute the total ionization loss
298
299 if (dedx < 0.0) dedx = 0.0 ;
300
301 dedx *= twopi_mc2_rcl2*chargeSquare*eDensity/beta2;
302
[961]303 //High order correction different for hadrons and ions
304 if(isIon) {
305 dedx += corr->IonBarkasCorrection(p,material,kineticEnergy);
306 } else {
307 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
308 }
[819]309 return dedx;
310}
311
312//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
313
[961]314void G4BetheBlochModel::CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
315 const G4DynamicParticle* dp,
316 G4double& eloss,
317 G4double&,
318 G4double length)
319{
320 const G4ParticleDefinition* p = dp->GetDefinition();
321 const G4Material* mat = couple->GetMaterial();
322 G4double preKinEnergy = dp->GetKineticEnergy();
323 G4double e = preKinEnergy - eloss*0.5;
324 if(e < 0.0) e = preKinEnergy*0.5;
325
326 if(isIon) {
327 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
328 GetModelOfFluctuations()->SetParticleAndCharge(p, q2);
329 eloss *= q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor;
330 eloss += length*corr->IonHighOrderCorrections(p,couple,e);
331 }
332
333 if(nuclearStopping && preKinEnergy*proton_mass_c2/mass < chargeSquare*100.*MeV) {
334
335 G4double nloss = length*corr->NuclearDEDX(p,mat,e,false);
336
337 // too big energy loss
338 if(eloss + nloss > preKinEnergy) {
339 nloss *= (preKinEnergy/(eloss + nloss));
340 eloss = preKinEnergy;
341 } else {
342 eloss += nloss;
343 }
344 /*
345 G4cout << "G4ionIonisation::CorrectionsAlongStep: e= " << preKinEnergy
346 << " de= " << eloss << " NIEL= " << nloss
347 << " dynQ= " << dp->GetCharge()/eplus << G4endl;
348 */
349 fParticleChange->ProposeNonIonizingEnergyDeposit(nloss);
350 }
351
352}
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
355
[819]356void G4BetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
357 const G4MaterialCutsCouple*,
358 const G4DynamicParticle* dp,
359 G4double minKinEnergy,
360 G4double maxEnergy)
361{
362 G4double kineticEnergy = dp->GetKineticEnergy();
363 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
364
[961]365 G4double maxKinEnergy = std::min(maxEnergy,tmax);
[819]366 if(minKinEnergy >= maxKinEnergy) return;
367
368 G4double totEnergy = kineticEnergy + mass;
369 G4double etot2 = totEnergy*totEnergy;
370 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
371
[961]372 G4double deltaKinEnergy, f;
373 G4double f1 = 0.0;
374 G4double fmax = 1.0;
375 if( 0.5 == spin ) fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2;
[819]376
[961]377 // sampling without nuclear size effect
[819]378 do {
379 G4double q = G4UniformRand();
380 deltaKinEnergy = minKinEnergy*maxKinEnergy
381 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
382
383 f = 1.0 - beta2*deltaKinEnergy/tmax;
[961]384 if( 0.5 == spin ) {
385 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
386 f += f1;
387 }
[819]388
[961]389 } while( fmax*G4UniformRand() > f);
390
391 // projectile formfactor - suppresion of high energy
392 // delta-electron production at high energy
393
394 G4double x = formfact*deltaKinEnergy;
395 if(x > 1.e-6) {
396
397 G4double x1 = 1.0 + x;
398 G4double g = 1.0/(x1*x1);
399 if( 0.5 == spin ) {
400 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
401 g *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
[819]402 }
[961]403 if(g > 1.0) {
404 G4cout << "### G4BetheBlochModel WARNING: g= " << g
405 << dp->GetDefinition()->GetParticleName()
406 << " Ekin(MeV)= " << kineticEnergy
407 << " delEkin(MeV)= " << deltaKinEnergy
408 << G4endl;
409 }
410 if(G4UniformRand() > g) return;
411 }
[819]412
[961]413 // delta-electron is produced
[819]414 G4double totMomentum = totEnergy*sqrt(beta2);
415 G4double deltaMomentum =
416 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
417 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
418 (deltaMomentum * totMomentum);
419 if(cost > 1.0) {
420 G4cout << "### G4BetheBlochModel WARNING: cost= "
421 << cost << " > 1 for "
422 << dp->GetDefinition()->GetParticleName()
423 << " Ekin(MeV)= " << kineticEnergy
424 << " p(MeV/c)= " << totMomentum
425 << " delEkin(MeV)= " << deltaKinEnergy
426 << " delMom(MeV/c)= " << deltaMomentum
427 << " tmin(MeV)= " << minKinEnergy
428 << " tmax(MeV)= " << maxKinEnergy
429 << " dir= " << dp->GetMomentumDirection()
430 << G4endl;
431 cost = 1.0;
432 }
433 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
434
435 G4double phi = twopi * G4UniformRand() ;
436
437
438 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
439 G4ThreeVector direction = dp->GetMomentumDirection();
440 deltaDirection.rotateUz(direction);
441
442 // create G4DynamicParticle object for delta ray
443 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
444 deltaDirection,deltaKinEnergy);
445
446 vdp->push_back(delta);
447
448 // Change kinematics of primary particle
449 kineticEnergy -= deltaKinEnergy;
450 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
451 finalP = finalP.unit();
452
453 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
454 fParticleChange->SetProposedMomentumDirection(finalP);
455}
456
[1005]457//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
458
459G4double G4BetheBlochModel::MaxSecondaryEnergy(const G4ParticleDefinition* pd,
460 G4double kinEnergy)
461{
462 if(isIon) SetParticle(pd);
463 G4double tau = kinEnergy/mass;
464 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
465 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
466 return std::min(tmax,tlimit);
467}
468
[819]469//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.