Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MollerBhabhaModel.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 file
31//
32//
33// File name: G4MollerBhabhaModel
34//
35// Author: Vladimir Ivanchenko on base of Laszlo Urban code
36//
37// Creation date: 03.01.2002
38//
39// Modifications:
40//
41// 13-11-02 Minor fix - use normalised direction (V.Ivanchenko)
42// 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (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// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
47// 25-07-05 Add protection in calculation of recoil direction for the case
48// of complete energy transfer from e+ to e- (V.Ivanchenko)
49// 06-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
50// 15-05-06 Fix MinEnergyCut (V.Ivanchenko)
51//
52//
53// Class Description:
54//
55// Implementation of energy loss and delta-electron production by e+/e-
56//
57// -------------------------------------------------------------------
58//
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
64#include "G4SystemOfUnits.hh"
65#include "G4Electron.hh"
66#include "G4Positron.hh"
67#include "Randomize.hh"
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
72using namespace std;
73
75 const G4String& nam)
76 : G4VEmModel(nam),
77 particle(0),
78 isElectron(true),
79 twoln10(2.0*log(10.0)),
80 lowLimit(0.02*keV),
81 isInitialised(false)
82{
84 if(p) { SetParticle(p); }
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
91{}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94
96 G4double kinEnergy)
97{
98 G4double tmax = kinEnergy;
99 if(isElectron) { tmax *= 0.5; }
100 return tmax;
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104
106 const G4DataVector&)
107{
108 if(!particle) { SetParticle(p); }
109
110 if(isInitialised) { return; }
111
112 isInitialised = true;
114}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117
120 G4double kineticEnergy,
121 G4double cutEnergy,
122 G4double maxEnergy)
123{
124 if(!particle) { SetParticle(p); }
125
126 G4double cross = 0.0;
127 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
128 tmax = std::min(maxEnergy, tmax);
129
130 if(cutEnergy < tmax) {
131
132 G4double xmin = cutEnergy/kineticEnergy;
133 G4double xmax = tmax/kineticEnergy;
134 G4double tau = kineticEnergy/electron_mass_c2;
135 G4double gam = tau + 1.0;
136 G4double gamma2= gam*gam;
137 G4double beta2 = tau*(tau + 2)/gamma2;
138
139 //Moller (e-e-) scattering
140 if (isElectron) {
141
142 G4double gg = (2.0*gam - 1.0)/gamma2;
143 cross = ((xmax - xmin)*(1.0 - gg + 1.0/(xmin*xmax)
144 + 1.0/((1.0-xmin)*(1.0 - xmax)))
145 - gg*log( xmax*(1.0 - xmin)/(xmin*(1.0 - xmax)) ) ) / beta2;
146
147 //Bhabha (e+e-) scattering
148 } else {
149
150 G4double y = 1.0/(1.0 + gam);
151 G4double y2 = y*y;
152 G4double y12 = 1.0 - 2.0*y;
153 G4double b1 = 2.0 - y2;
154 G4double b2 = y12*(3.0 + y2);
155 G4double y122= y12*y12;
156 G4double b4 = y122*y12;
157 G4double b3 = b4 + y122;
158
159 cross = (xmax - xmin)*(1.0/(beta2*xmin*xmax) + b2
160 - 0.5*b3*(xmin + xmax)
161 + b4*(xmin*xmin + xmin*xmax + xmax*xmax)/3.0)
162 - b1*log(xmax/xmin);
163 }
164
165 cross *= twopi_mc2_rcl2/kineticEnergy;
166 }
167 return cross;
168}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171
173 const G4ParticleDefinition* p,
174 G4double kineticEnergy,
176 G4double cutEnergy,
177 G4double maxEnergy)
178{
179 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183
185 const G4Material* material,
186 const G4ParticleDefinition* p,
187 G4double kinEnergy,
188 G4double cutEnergy,
189 G4double maxEnergy)
190{
191 G4double eDensity = material->GetElectronDensity();
192 return eDensity*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
193 /*
194 G4double Zeff = eDensity/material->GetTotNbOfAtomsPerVolume();
195 G4double th = 0.25*sqrt(Zeff)*keV;
196 G4double cut;
197 if(isElectron) { cut = std::max(th*0.5, cutEnergy); }
198 else { cut = std::max(th, cutEnergy); }
199 G4double res = 0.0;
200 // below this threshold no bremsstrahlung
201 if (kinEnergy > th) {
202 res = eDensity*ComputeCrossSectionPerElectron(p,kinEnergy,cut,maxEnergy);
203 }
204 return res;
205 */
206}
207
208//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
209
211 const G4Material* material,
212 const G4ParticleDefinition* p,
213 G4double kineticEnergy,
214 G4double cut)
215{
216 if(!particle) { SetParticle(p); }
217 // calculate the dE/dx due to the ionization by Seltzer-Berger formula
218 // checl low-energy limit
219 G4double electronDensity = material->GetElectronDensity();
220
221 G4double Zeff = electronDensity/material->GetTotNbOfAtomsPerVolume();
222 G4double th = 0.25*sqrt(Zeff)*keV;
223 // G4double cut;
224 // if(isElectron) { cut = std::max(th*0.5, cutEnergy); }
225 // else { cut = std::max(th, cutEnergy); }
226 G4double tkin = kineticEnergy;
227 if (kineticEnergy < th) { tkin = th; }
228
229 G4double tau = tkin/electron_mass_c2;
230 G4double gam = tau + 1.0;
231 G4double gamma2= gam*gam;
232 G4double bg2 = tau*(tau + 2);
233 G4double beta2 = bg2/gamma2;
234
235 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
236 eexc /= electron_mass_c2;
237 G4double eexc2 = eexc*eexc;
238
239 G4double d = std::min(cut, MaxSecondaryEnergy(p, tkin))/electron_mass_c2;
240 G4double dedx;
241
242 // electron
243 if (isElectron) {
244
245 dedx = log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
246 + log((tau-d)*d) + tau/(tau-d)
247 + (0.5*d*d + (2.0*tau + 1.)*log(1. - d/tau))/gamma2;
248
249 //positron
250 } else {
251
252 G4double d2 = d*d*0.5;
253 G4double d3 = d2*d/1.5;
254 G4double d4 = d3*d*0.75;
255 G4double y = 1.0/(1.0 + gam);
256 dedx = log(2.0*(tau + 2.0)/eexc2) + log(tau*d)
257 - beta2*(tau + 2.0*d - y*(3.0*d2
258 + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
259 }
260
261 //density correction
262 G4double x = log(bg2)/twoln10;
263 dedx -= material->GetIonisation()->DensityCorrection(x);
264
265 // now you can compute the total ionization loss
266 dedx *= twopi_mc2_rcl2*electronDensity/beta2;
267 if (dedx < 0.0) { dedx = 0.0; }
268
269 // lowenergy extrapolation
270
271 if (kineticEnergy < th) {
272 x = kineticEnergy/th;
273 if(x > 0.25) { dedx /= sqrt(x); }
274 else { dedx *= 1.4*sqrt(x)/(0.1 + x); }
275 }
276 return dedx;
277}
278
279//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280
281void G4MollerBhabhaModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
283 const G4DynamicParticle* dp,
284 G4double cutEnergy,
285 G4double maxEnergy)
286{
287 G4double kineticEnergy = dp->GetKineticEnergy();
288 //const G4Material* mat = couple->GetMaterial();
289 //G4double Zeff = mat->GetElectronDensity()/mat->GetTotNbOfAtomsPerVolume();
290 // G4double th = 0.25*sqrt(Zeff)*keV;
291 G4double tmax;
292 G4double tmin = cutEnergy;
293 if(isElectron) {
294 tmax = 0.5*kineticEnergy;
295 } else {
296 tmax = kineticEnergy;
297 }
298 if(maxEnergy < tmax) { tmax = maxEnergy; }
299 if(tmin >= tmax) { return; }
300
301 G4double energy = kineticEnergy + electron_mass_c2;
302 G4double totalMomentum = sqrt(kineticEnergy*(energy + electron_mass_c2));
303 G4double xmin = tmin/kineticEnergy;
304 G4double xmax = tmax/kineticEnergy;
305 G4double gam = energy/electron_mass_c2;
306 G4double gamma2 = gam*gam;
307 G4double beta2 = 1.0 - 1.0/gamma2;
308 G4double x, z, q, grej;
309
310 G4ThreeVector direction = dp->GetMomentumDirection();
311
312 //Moller (e-e-) scattering
313 if (isElectron) {
314
315 G4double gg = (2.0*gam - 1.0)/gamma2;
316 G4double y = 1.0 - xmax;
317 grej = 1.0 - gg*xmax + xmax*xmax*(1.0 - gg + (1.0 - gg*y)/(y*y));
318
319 do {
320 q = G4UniformRand();
321 x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
322 y = 1.0 - x;
323 z = 1.0 - gg*x + x*x*(1.0 - gg + (1.0 - gg*y)/(y*y));
324 /*
325 if(z > grej) {
326 G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
327 << "Majorant " << grej << " < "
328 << z << " for x= " << x
329 << " e-e- scattering"
330 << G4endl;
331 }
332 */
333 } while(grej * G4UniformRand() > z);
334
335 //Bhabha (e+e-) scattering
336 } else {
337
338 G4double y = 1.0/(1.0 + gam);
339 G4double y2 = y*y;
340 G4double y12 = 1.0 - 2.0*y;
341 G4double b1 = 2.0 - y2;
342 G4double b2 = y12*(3.0 + y2);
343 G4double y122= y12*y12;
344 G4double b4 = y122*y12;
345 G4double b3 = b4 + y122;
346
347 y = xmax*xmax;
348 grej = 1.0 + (y*y*b4 - xmin*xmin*xmin*b3 + y*b2 - xmin*b1)*beta2;
349 do {
350 q = G4UniformRand();
351 x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
352 y = x*x;
353 z = 1.0 + (y*y*b4 - x*y*b3 + y*b2 - x*b1)*beta2;
354 /*
355 if(z > grej) {
356 G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
357 << "Majorant " << grej << " < "
358 << z << " for x= " << x
359 << " e+e- scattering"
360 << G4endl;
361 }
362 */
363 } while(grej * G4UniformRand() > z);
364 }
365
366 G4double deltaKinEnergy = x * kineticEnergy;
367
368 G4double deltaMomentum =
369 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
370 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
371 (deltaMomentum * totalMomentum);
372 G4double sint = (1.0 - cost)*(1. + cost);
373 if(sint > 0.0) { sint = sqrt(sint); }
374 else { sint = 0.0; }
375
376 G4double phi = twopi * G4UniformRand() ;
377
378 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
379 deltaDirection.rotateUz(direction);
380
381 // primary change
382 kineticEnergy -= deltaKinEnergy;
384
385 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
386 direction = dir.unit();
388
389 // create G4DynamicParticle object for delta ray
391 deltaDirection,deltaKinEnergy);
392 vdp->push_back(delta);
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
double G4double
Definition: G4Types.hh:64
#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 GetMeanExcitationEnergy() const
G4double DensityCorrection(G4double x)
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:208
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
G4double GetElectronDensity() const
Definition: G4Material.hh:216
const G4ParticleDefinition * particle
void SetParticle(const G4ParticleDefinition *p)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
G4ParticleChangeForLoss * fParticleChange
G4ParticleDefinition * theElectron
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4MollerBhabhaModel(const G4ParticleDefinition *p=0, const G4String &nam="MollerBhabha")
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95