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

Last change on this file since 1348 was 1340, checked in by garnier, 14 years ago

update ti head

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