Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuBetheBlochModel.cc
Go to the documentation of this file.
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$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4MuBetheBlochModel
34//
35// Author: Vladimir Ivanchenko on base of Laszlo Urban code
36//
37// Creation date: 09.08.2002
38//
39// Modifications:
40//
41// 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
42// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
43// 27-01-03 Make models region aware (V.Ivanchenko)
44// 13-02-03 Add name (V.Ivanchenko)
45// 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I)
46// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
47// 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko)
48// 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49//
50
51//
52// -------------------------------------------------------------------
53//
54
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
61#include "G4SystemOfUnits.hh"
62#include "Randomize.hh"
63#include "G4Electron.hh"
64#include "G4LossTableManager.hh"
65#include "G4EmCorrections.hh"
67
68G4double G4MuBetheBlochModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917,
69 0.7628, 0.8983, 0.9801 };
70
71G4double G4MuBetheBlochModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813,
72 0.1569, 0.1112, 0.0506 };
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75
76using namespace std;
77
79 const G4String& nam)
80 : G4VEmModel(nam),
81 particle(0),
82 limitKinEnergy(100.*keV),
83 logLimitKinEnergy(log(limitKinEnergy)),
84 twoln10(2.0*log(10.0)),
85 bg2lim(0.0169),
86 taulim(8.4146e-3),
87 alphaprime(fine_structure_const/twopi)
88{
89 theElectron = G4Electron::Electron();
91 fParticleChange = 0;
92
93 // initial initialisation of memeber should be overwritten
94 // by SetParticle
95 mass = massSquare = ratio = 1.0;
96
97 if(p) { SetParticle(p); }
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
103{}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106
108 const G4MaterialCutsCouple* couple)
109{
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
114
116 G4double kinEnergy)
117{
118 G4double tau = kinEnergy/mass;
119 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
120 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
121 return tmax;
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
127 const G4DataVector&)
128{
129 if(p) { SetParticle(p); }
130 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
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
185 const G4ParticleDefinition* p,
186 G4double kineticEnergy,
188 G4double cutEnergy,
189 G4double maxEnergy)
190{
192 (p,kineticEnergy,cutEnergy,maxEnergy);
193 return cross;
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197
199 const G4Material* material,
200 const G4ParticleDefinition* p,
201 G4double kineticEnergy,
202 G4double cutEnergy,
203 G4double maxEnergy)
204{
205 G4double eDensity = material->GetElectronDensity();
207 (p,kineticEnergy,cutEnergy,maxEnergy);
208 return cross;
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
212
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 dedx -= material->GetIonisation()->DensityCorrection(x);
249
250 // shell correction
251 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
252
253 // now compute the total ionization loss
254
255 if (dedx < 0.0) dedx = 0.0 ;
256
257 // radiative corrections of R. Kokoulin
258 if (cutEnergy > limitKinEnergy) {
259
260 G4double logtmax = log(cutEnergy);
261 G4double logstep = logtmax - logLimitKinEnergy;
262 G4double dloss = 0.0;
263 G4double ftot2= 0.5/(totEnergy*totEnergy);
264
265 for (G4int ll=0; ll<8; ll++)
266 {
267 G4double ep = exp(logLimitKinEnergy + xgi[ll]*logstep);
268 G4double a1 = log(1.0 + 2.0*ep/electron_mass_c2);
269 G4double a3 = log(4.0*totEnergy*(totEnergy - ep)/massSquare);
270 dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
271 }
272 dedx += dloss*logstep*alphaprime;
273 }
274
275 dedx *= twopi_mc2_rcl2*eDensity/beta2;
276
277 //High order corrections
278 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
279
280 return dedx;
281}
282
283//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
284
285void G4MuBetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
287 const G4DynamicParticle* dp,
288 G4double minKinEnergy,
289 G4double maxEnergy)
290{
292 G4double maxKinEnergy = min(maxEnergy,tmax);
293 if(minKinEnergy >= maxKinEnergy) { return; }
294
295 G4double kineticEnergy = dp->GetKineticEnergy();
296 G4double totEnergy = kineticEnergy + mass;
297 G4double etot2 = totEnergy*totEnergy;
298 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
299
300 G4double grej = 1.;
301 if(tmax > limitKinEnergy) {
302 G4double a0 = log(2.*totEnergy/mass);
303 grej += alphaprime*a0*a0;
304 }
305
306 G4double deltaKinEnergy, f;
307
308 // sampling follows ...
309 do {
311 deltaKinEnergy = minKinEnergy*maxKinEnergy
312 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
313
314
315 f = 1.0 - beta2*deltaKinEnergy/tmax
316 + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
317
318 if(deltaKinEnergy > limitKinEnergy) {
319 G4double a1 = log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
320 G4double a3 = log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
321 f *= (1. + alphaprime*a1*(a3 - a1));
322 }
323
324 if(f > grej) {
325 G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
326 << "Majorant " << grej << " < "
327 << f << " for edelta= " << deltaKinEnergy
328 << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
329 << G4endl;
330 }
331
332
333 } while( grej*G4UniformRand() > f );
334
335 G4double deltaMomentum =
336 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
337 G4double totalMomentum = totEnergy*sqrt(beta2);
338 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
339 (deltaMomentum * totalMomentum);
340
341 G4double sint = sqrt(1.0 - cost*cost);
342
343 G4double phi = twopi * G4UniformRand() ;
344
345 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
346 G4ThreeVector direction = dp->GetMomentumDirection();
347 deltaDirection.rotateUz(direction);
348
349 // primary change
350 kineticEnergy -= deltaKinEnergy;
351 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
352 direction = dir.unit();
353 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
354 fParticleChange->SetProposedMomentumDirection(direction);
355
356 // create G4DynamicParticle object for delta ray
357 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
358 deltaDirection,deltaKinEnergy);
359 vdp->push_back(delta);
360}
361
362//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetMeanExcitationEnergy() const
G4double DensityCorrection(G4double x)
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
G4double GetElectronDensity() const
Definition: G4Material.hh:216
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
G4MuBetheBlochModel(const G4ParticleDefinition *p=0, const G4String &nam="MuBetheBloch")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:399
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95