source: trunk/source/processes/electromagnetic/muons/src/G4MuBetheBlochModel.cc @ 961

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

update processes

File size: 12.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: G4MuBetheBlochModel.cc,v 1.25 2009/02/20 14:48:16 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4MuBetheBlochModel
35//
36// Author:        Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 09.08.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// 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I)
47// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
48// 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko)
49// 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
50//
51
52//
53// -------------------------------------------------------------------
54//
55
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
60#include "G4MuBetheBlochModel.hh"
61#include "Randomize.hh"
62#include "G4Electron.hh"
63#include "G4LossTableManager.hh"
64#include "G4EmCorrections.hh"
65#include "G4ParticleChangeForLoss.hh"
66
67G4double G4MuBetheBlochModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917,
68                                      0.7628, 0.8983, 0.9801 };
69                                     
70G4double G4MuBetheBlochModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813,
71                                      0.1569, 0.1112, 0.0506 };
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74
75using namespace std;
76
77G4MuBetheBlochModel::G4MuBetheBlochModel(const G4ParticleDefinition* p,
78                                         const G4String& nam)
79  : G4VEmModel(nam),
80  particle(0),
81  limitKinEnergy(100.*keV),
82  logLimitKinEnergy(log(limitKinEnergy)),
83  twoln10(2.0*log(10.0)),
84  bg2lim(0.0169),
85  taulim(8.4146e-3),
86  alphaprime(fine_structure_const/twopi)
87{
88  theElectron = G4Electron::Electron();
89  corr = G4LossTableManager::Instance()->EmCorrections();
90  fParticleChange = 0;
91
92  if(p) SetParticle(p);
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96
97G4MuBetheBlochModel::~G4MuBetheBlochModel()
98{}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
102G4double G4MuBetheBlochModel::MinEnergyCut(const G4ParticleDefinition*,
103                                           const G4MaterialCutsCouple* couple)
104{
105  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
109
110G4double G4MuBetheBlochModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
111                                                 G4double kinEnergy) 
112{
113  G4double tau  = kinEnergy/mass;
114  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
115                  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
116  return tmax;
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
121void G4MuBetheBlochModel::Initialise(const G4ParticleDefinition* p,
122                                     const G4DataVector&)
123{
124  if(p) SetParticle(p);
125
126  if(!fParticleChange) {
127    if(pParticleChange) {
128      fParticleChange = 
129        reinterpret_cast<G4ParticleChangeForLoss*>(pParticleChange);
130    } else {
131      fParticleChange = new G4ParticleChangeForLoss();
132    }
133  }
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
138G4double G4MuBetheBlochModel::ComputeCrossSectionPerElectron(
139                                           const G4ParticleDefinition* p,
140                                                 G4double kineticEnergy,
141                                                 G4double cutEnergy,
142                                                 G4double maxKinEnergy)
143{
144  G4double cross = 0.0;
145  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
146  G4double maxEnergy = min(tmax,maxKinEnergy);
147  if(cutEnergy < maxEnergy) {
148
149    G4double totEnergy = kineticEnergy + mass;
150    G4double energy2   = totEnergy*totEnergy;
151    G4double beta2     = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
152
153    cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax
154          + 0.5*(maxEnergy - cutEnergy)/energy2;
155
156    // radiative corrections of R. Kokoulin
157    if (maxEnergy > limitKinEnergy) {
158
159      G4double logtmax = log(maxEnergy);
160      G4double logtmin = log(max(cutEnergy,limitKinEnergy));
161      G4double logstep = logtmax - logtmin;
162      G4double dcross  = 0.0;
163
164      for (G4int ll=0; ll<8; ll++)
165      {
166        G4double ep = exp(logtmin + xgi[ll]*logstep);
167        G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2);
168        G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare);
169        dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
170      }
171
172      cross += dcross*logstep*alphaprime;
173    }
174
175    cross *= twopi_mc2_rcl2/beta2;
176
177  }
178
179  //  G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
180  //         << " cross= " << cross << G4endl;
181 
182  return cross;
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186
187G4double G4MuBetheBlochModel::ComputeCrossSectionPerAtom(
188                                           const G4ParticleDefinition* p,
189                                                 G4double kineticEnergy,
190                                                 G4double Z, G4double,
191                                                 G4double cutEnergy,
192                                                 G4double maxEnergy)
193{
194  G4double cross = Z*ComputeCrossSectionPerElectron
195                                         (p,kineticEnergy,cutEnergy,maxEnergy);
196  return cross;
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
200
201G4double G4MuBetheBlochModel::CrossSectionPerVolume(
202                                           const G4Material* material,
203                                           const G4ParticleDefinition* p,
204                                                 G4double kineticEnergy,
205                                                 G4double cutEnergy,
206                                                 G4double maxEnergy)
207{
208  G4double eDensity = material->GetElectronDensity();
209  G4double cross = eDensity*ComputeCrossSectionPerElectron
210                                         (p,kineticEnergy,cutEnergy,maxEnergy);
211  return cross;
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
216G4double G4MuBetheBlochModel::ComputeDEDXPerVolume(const G4Material* material,
217                                                  const G4ParticleDefinition* p,
218                                                  G4double kineticEnergy,
219                                                  G4double cut)
220{
221  G4double tmax  = MaxSecondaryEnergy(p, kineticEnergy);
222  G4double tau   = kineticEnergy/mass;
223  G4double cutEnergy = min(cut,tmax);
224  G4double gam   = tau + 1.0;
225  G4double bg2   = tau * (tau+2.0);
226  G4double beta2 = bg2/(gam*gam);
227
228  G4double eexc  = material->GetIonisation()->GetMeanExcitationEnergy();
229  G4double eexc2 = eexc*eexc;
230  G4double cden  = material->GetIonisation()->GetCdensity();
231  G4double mden  = material->GetIonisation()->GetMdensity();
232  G4double aden  = material->GetIonisation()->GetAdensity();
233  G4double x0den = material->GetIonisation()->GetX0density();
234  G4double x1den = material->GetIonisation()->GetX1density();
235
236  G4double eDensity = material->GetElectronDensity();
237
238  G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
239                 -(1.0 + cutEnergy/tmax)*beta2;
240
241  G4double totEnergy = kineticEnergy + mass;
242  G4double del = 0.5*cutEnergy/totEnergy;
243  dedx += del*del;
244
245  // density correction
246  G4double x = log(bg2)/twoln10;
247  if ( x >= x0den ) {
248    dedx -= twoln10*x - cden ;
249    if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
250  }
251
252  // shell correction
253  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
254
255  // now compute the total ionization loss
256
257  if (dedx < 0.0) dedx = 0.0 ;
258
259  // radiative corrections of R. Kokoulin
260  if (cutEnergy > limitKinEnergy) {
261
262    G4double logtmax = log(cutEnergy);
263    G4double logstep = logtmax - logLimitKinEnergy;
264    G4double dloss = 0.0;
265    G4double ftot2= 0.5/(totEnergy*totEnergy);
266
267    for (G4int ll=0; ll<8; ll++)
268    {
269      G4double ep = exp(logLimitKinEnergy + xgi[ll]*logstep);
270      G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2);
271      G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare);
272      dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
273    }
274    dedx += dloss*logstep*alphaprime;
275  }
276
277  dedx *= twopi_mc2_rcl2*eDensity/beta2;
278
279  //High order corrections
280  dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
281
282  return dedx;
283}
284
285//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
286
287void G4MuBetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
288                                            const G4MaterialCutsCouple*,
289                                            const G4DynamicParticle* dp,
290                                            G4double minKinEnergy,
291                                            G4double maxEnergy)
292{
293  G4double tmax = MaxSecondaryKinEnergy(dp);
294  G4double maxKinEnergy = min(maxEnergy,tmax);
295  if(minKinEnergy >= maxKinEnergy) return;
296
297  G4double kineticEnergy = dp->GetKineticEnergy();
298  G4double totEnergy     = kineticEnergy + mass;
299  G4double etot2         = totEnergy*totEnergy;
300  G4double beta2         = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
301 
302  G4double grej  = 1.;
303  if(tmax > limitKinEnergy) {
304    G4double a0    = log(2.*totEnergy/mass);
305    grej  += alphaprime*a0*a0;
306  }
307
308  G4double deltaKinEnergy, f;
309
310  // sampling follows ...
311  do {
312    G4double q = G4UniformRand();
313    deltaKinEnergy = minKinEnergy*maxKinEnergy
314                    /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
315
316
317    f = 1.0 - beta2*deltaKinEnergy/tmax
318            + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
319
320    if(deltaKinEnergy > limitKinEnergy) {
321      G4double a1 = log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
322      G4double a3 = log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
323      f *= (1. + alphaprime*a1*(a3 - a1));
324    }
325
326    if(f > grej) {
327        G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
328               << "Majorant " << grej << " < "
329               << f << " for edelta= " << deltaKinEnergy
330               << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
331               << G4endl;
332    }
333
334
335  } while( grej*G4UniformRand() > f );
336
337  G4double deltaMomentum =
338           sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
339  G4double totalMomentum = totEnergy*sqrt(beta2);
340  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
341                                   (deltaMomentum * totalMomentum);
342
343  G4double sint = sqrt(1.0 - cost*cost);
344
345  G4double phi = twopi * G4UniformRand() ;
346
347  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
348  G4ThreeVector direction = dp->GetMomentumDirection();
349  deltaDirection.rotateUz(direction);
350
351  // primary change
352  kineticEnergy -= deltaKinEnergy;
353  G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
354  direction = dir.unit();
355  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
356  fParticleChange->SetProposedMomentumDirection(direction);
357
358  // create G4DynamicParticle object for delta ray
359  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
360                                                 deltaDirection,deltaKinEnergy);
361  vdp->push_back(delta);
362}
363
364//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.