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

Last change on this file since 1055 was 1055, checked in by garnier, 15 years ago

maj sur la beta de geant 4.9.3

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