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

Last change on this file since 836 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 12.8 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.23 2007/05/22 17:35:58 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-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
91  if(p) SetParticle(p);
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
96G4MuBetheBlochModel::~G4MuBetheBlochModel()
97{}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
101void G4MuBetheBlochModel::SetParticle(const G4ParticleDefinition* p)
102{
103  if(!particle) {
104    particle = p;
105    mass = particle->GetPDGMass();
106    massSquare = mass*mass;
107    ratio = electron_mass_c2/mass;
108  }
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
113G4double G4MuBetheBlochModel::MinEnergyCut(const G4ParticleDefinition*,
114                                           const G4MaterialCutsCouple* couple)
115{
116  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
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(pParticleChange)
127    fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*>
128                                                             (pParticleChange);
129  else
130    fParticleChange = new G4ParticleChangeForLoss();
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
135G4double G4MuBetheBlochModel::ComputeCrossSectionPerElectron(
136                                           const G4ParticleDefinition* p,
137                                                 G4double kineticEnergy,
138                                                 G4double cutEnergy,
139                                                 G4double maxKinEnergy)
140{
141  G4double cross = 0.0;
142  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
143  G4double maxEnergy = min(tmax,maxKinEnergy);
144  if(cutEnergy < maxEnergy) {
145
146    G4double totEnergy = kineticEnergy + mass;
147    G4double energy2   = totEnergy*totEnergy;
148    G4double beta2     = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
149
150    cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax
151          + 0.5*(maxEnergy - cutEnergy)/energy2;
152
153    // radiative corrections of R. Kokoulin
154    if (maxEnergy > limitKinEnergy) {
155
156      G4double logtmax = log(maxEnergy);
157      G4double logtmin = log(max(cutEnergy,limitKinEnergy));
158      G4double logstep = logtmax - logtmin;
159      G4double dcross  = 0.0;
160
161      for (G4int ll=0; ll<8; ll++)
162      {
163        G4double ep = exp(logtmin + xgi[ll]*logstep);
164        G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2);
165        G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare);
166        dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
167      }
168
169      cross += dcross*logstep*alphaprime;
170    }
171
172    cross *= twopi_mc2_rcl2/beta2;
173
174  }
175
176  //  G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
177  //         << " cross= " << cross << G4endl;
178 
179  return cross;
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183
184G4double G4MuBetheBlochModel::ComputeCrossSectionPerAtom(
185                                           const G4ParticleDefinition* p,
186                                                 G4double kineticEnergy,
187                                                 G4double Z, G4double,
188                                                 G4double cutEnergy,
189                                                 G4double maxEnergy)
190{
191  G4double cross = Z*ComputeCrossSectionPerElectron
192                                         (p,kineticEnergy,cutEnergy,maxEnergy);
193  return cross;
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197
198G4double G4MuBetheBlochModel::CrossSectionPerVolume(
199                                           const G4Material* material,
200                                           const G4ParticleDefinition* p,
201                                                 G4double kineticEnergy,
202                                                 G4double cutEnergy,
203                                                 G4double maxEnergy)
204{
205  G4double eDensity = material->GetElectronDensity();
206  G4double cross = eDensity*ComputeCrossSectionPerElectron
207                                         (p,kineticEnergy,cutEnergy,maxEnergy);
208  return cross;
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
212
213G4double G4MuBetheBlochModel::ComputeDEDXPerVolume(const G4Material* material,
214                                                  const G4ParticleDefinition* p,
215                                                  G4double kineticEnergy,
216                                                  G4double cut)
217{
218  G4double tmax  = MaxSecondaryEnergy(p, kineticEnergy);
219  G4double tau   = kineticEnergy/mass;
220  G4double cutEnergy = min(cut,tmax);
221  G4double gam   = tau + 1.0;
222  G4double bg2   = tau * (tau+2.0);
223  G4double beta2 = bg2/(gam*gam);
224
225  G4double eexc  = material->GetIonisation()->GetMeanExcitationEnergy();
226  G4double eexc2 = eexc*eexc;
227  G4double cden  = material->GetIonisation()->GetCdensity();
228  G4double mden  = material->GetIonisation()->GetMdensity();
229  G4double aden  = material->GetIonisation()->GetAdensity();
230  G4double x0den = material->GetIonisation()->GetX0density();
231  G4double x1den = material->GetIonisation()->GetX1density();
232
233  G4double eDensity = material->GetElectronDensity();
234
235  G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
236                 -(1.0 + cutEnergy/tmax)*beta2;
237
238  G4double totEnergy = kineticEnergy + mass;
239  G4double del = 0.5*cutEnergy/totEnergy;
240  dedx += del*del;
241
242  // density correction
243  G4double x = log(bg2)/twoln10;
244  if ( x >= x0den ) {
245    dedx -= twoln10*x - cden ;
246    if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ;
247  }
248
249  // shell correction
250  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
251
252  // now compute the total ionization loss
253
254  if (dedx < 0.0) dedx = 0.0 ;
255
256  // radiative corrections of R. Kokoulin
257  if (cutEnergy > limitKinEnergy) {
258
259    G4double logtmax = log(cutEnergy);
260    G4double logstep = logtmax - logLimitKinEnergy;
261    G4double dloss = 0.0;
262    G4double ftot2= 0.5/(totEnergy*totEnergy);
263
264    for (G4int ll=0; ll<8; ll++)
265    {
266      G4double ep = exp(logLimitKinEnergy + xgi[ll]*logstep);
267      G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2);
268      G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare);
269      dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
270    }
271    dedx += dloss*logstep*alphaprime;
272  }
273
274  dedx *= twopi_mc2_rcl2*eDensity/beta2;
275
276  //High order corrections
277  dedx += corr->HighOrderCorrections(p,material,kineticEnergy);
278
279  return dedx;
280}
281
282//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
283
284void G4MuBetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
285                                            const G4MaterialCutsCouple*,
286                                            const G4DynamicParticle* dp,
287                                            G4double minKinEnergy,
288                                            G4double maxEnergy)
289{
290  G4double tmax = MaxSecondaryKinEnergy(dp);
291  G4double maxKinEnergy = min(maxEnergy,tmax);
292  if(minKinEnergy >= maxKinEnergy) return;
293
294  G4double kineticEnergy = dp->GetKineticEnergy();
295  G4double totEnergy     = kineticEnergy + mass;
296  G4double etot2         = totEnergy*totEnergy;
297  G4double beta2         = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
298 
299  G4double grej  = 1.;
300  if(tmax > limitKinEnergy) {
301    G4double a0    = log(2.*totEnergy/mass);
302    grej  += alphaprime*a0*a0;
303  }
304
305  G4double deltaKinEnergy, f;
306
307  // sampling follows ...
308  do {
309    G4double q = G4UniformRand();
310    deltaKinEnergy = minKinEnergy*maxKinEnergy
311                    /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
312
313
314    f = 1.0 - beta2*deltaKinEnergy/tmax
315            + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
316
317    if(deltaKinEnergy > limitKinEnergy) {
318      G4double a1 = log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
319      G4double a3 = log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
320      f *= (1. + alphaprime*a1*(a3 - a1));
321    }
322
323    if(f > grej) {
324        G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
325               << "Majorant " << grej << " < "
326               << f << " for edelta= " << deltaKinEnergy
327               << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
328               << G4endl;
329    }
330
331
332  } while( grej*G4UniformRand() > f );
333
334  G4double deltaMomentum =
335           sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
336  G4double totalMomentum = totEnergy*sqrt(beta2);
337  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
338                                   (deltaMomentum * totalMomentum);
339
340  G4double sint = sqrt(1.0 - cost*cost);
341
342  G4double phi = twopi * G4UniformRand() ;
343
344  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
345  G4ThreeVector direction = dp->GetMomentumDirection();
346  deltaDirection.rotateUz(direction);
347
348  // primary change
349  kineticEnergy -= deltaKinEnergy;
350  G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
351  direction = dir.unit();
352  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
353  fParticleChange->SetProposedMomentumDirection(direction);
354
355  // create G4DynamicParticle object for delta ray
356  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
357                                                 deltaDirection,deltaKinEnergy);
358  vdp->push_back(delta);
359}
360
361//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.