source: trunk/source/processes/electromagnetic/muons/src/G4MuPairProductionModel.cc@ 1317

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 21.7 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: G4MuPairProductionModel.cc,v 1.45 2010/06/01 15:21:59 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4MuPairProductionModel
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 24.06.2002
39//
40// Modifications:
41//
42// 04-12-02 Change G4DynamicParticle constructor in PostStep (V.Ivanchenko)
43// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
44// 24-01-03 Fix for compounds (V.Ivanchenko)
45// 27-01-03 Make models region aware (V.Ivanchenko)
46// 13-02-03 Add model (V.Ivanchenko)
47// 06-06-03 Fix in cross section calculation for high energy (V.Ivanchenko)
48// 20-10-03 2*xi in ComputeDDMicroscopicCrossSection (R.Kokoulin)
49// 8 integration points in ComputeDMicroscopicCrossSection
50// 12-01-04 Take min cut of e- and e+ not its sum (V.Ivanchenko)
51// 10-02-04 Update parameterisation using R.Kokoulin model (V.Ivanchenko)
52// 28-04-04 For complex materials repeat calculation of max energy for each
53// material (V.Ivanchenko)
54// 01-11-04 Fix bug inside ComputeDMicroscopicCrossSection (R.Kokoulin)
55// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
56// 03-08-05 Add SetParticle method (V.Ivantchenko)
57// 23-10-05 Add protection in sampling of e+e- pair energy needed for
58// low cuts (V.Ivantchenko)
59// 13-02-06 Add ComputeCrossSectionPerAtom (mma)
60// 24-04-07 Add protection in SelectRandomAtom method (V.Ivantchenko)
61// 12-05-06 Updated sampling (use cut) in SelectRandomAtom (A.Bogdanov)
62// 11-10-07 Add ignoreCut flag (V.Ivanchenko)
63
64//
65// Class Description:
66//
67//
68// -------------------------------------------------------------------
69//
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72
73#include "G4MuPairProductionModel.hh"
74#include "G4Electron.hh"
75#include "G4Positron.hh"
76#include "G4MuonMinus.hh"
77#include "G4MuonPlus.hh"
78#include "Randomize.hh"
79#include "G4Material.hh"
80#include "G4Element.hh"
81#include "G4ElementVector.hh"
82#include "G4ProductionCutsTable.hh"
83#include "G4ParticleChangeForLoss.hh"
84#include "G4ParticleChangeForGamma.hh"
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
88// static members
89//
90G4double G4MuPairProductionModel::zdat[]={1., 4., 13., 29., 92.};
91G4double G4MuPairProductionModel::adat[]={1.01, 9.01, 26.98, 63.55, 238.03};
92G4double G4MuPairProductionModel::tdat[]={1.e3, 1.e4, 1.e5, 1.e6, 1.e7, 1.e8,
93 1.e9, 1.e10};
94G4double G4MuPairProductionModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083,
95 0.5917, 0.7628, 0.8983, 0.9801 };
96G4double G4MuPairProductionModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813,
97 0.1813, 0.1569, 0.1112, 0.0506 };
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
101using namespace std;
102
103G4MuPairProductionModel::G4MuPairProductionModel(const G4ParticleDefinition* p,
104 const G4String& nam)
105 : G4VEmModel(nam),
106 particle(0),
107 factorForCross(4.*fine_structure_const*fine_structure_const
108 *classic_electr_radius*classic_electr_radius/(3.*pi)),
109 sqrte(sqrt(exp(1.))),
110 currentZ(0),
111 fParticleChange(0),
112 minPairEnergy(4.*electron_mass_c2),
113 lowestKinEnergy(1.*GeV),
114 nzdat(5),
115 ntdat(8),
116 nbiny(1000),
117 nmaxElements(0),
118 ymin(-5.),
119 ymax(0.),
120 dy((ymax-ymin)/nbiny),
121 samplingTablesAreFilled(false)
122{
123 SetLowEnergyLimit(minPairEnergy);
124 nist = G4NistManager::Instance();
125
126 theElectron = G4Electron::Electron();
127 thePositron = G4Positron::Positron();
128
129 if(p) SetParticle(p);
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133
134G4MuPairProductionModel::~G4MuPairProductionModel()
135{}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138
139G4double G4MuPairProductionModel::MinEnergyCut(const G4ParticleDefinition*,
140 const G4MaterialCutsCouple* )
141{
142 return minPairEnergy;
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146
147G4double G4MuPairProductionModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
148 G4double kineticEnergy)
149{
150 G4double maxPairEnergy = kineticEnergy + particleMass*(1.0 - 0.75*sqrte*z13);
151 return maxPairEnergy;
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155
156void G4MuPairProductionModel::Initialise(const G4ParticleDefinition* p,
157 const G4DataVector&)
158{
159 if (!samplingTablesAreFilled) {
160 if(p) SetParticle(p);
161 MakeSamplingTables();
162 }
163 if(!fParticleChange) fParticleChange = GetParticleChangeForLoss();
164}
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167
168G4double G4MuPairProductionModel::ComputeDEDXPerVolume(
169 const G4Material* material,
170 const G4ParticleDefinition*,
171 G4double kineticEnergy,
172 G4double cutEnergy)
173{
174 G4double dedx = 0.0;
175 if (cutEnergy <= minPairEnergy || kineticEnergy <= lowestKinEnergy)
176 return dedx;
177
178 const G4ElementVector* theElementVector = material->GetElementVector();
179 const G4double* theAtomicNumDensityVector =
180 material->GetAtomicNumDensityVector();
181
182 // loop for elements in the material
183 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
184 G4double Z = (*theElementVector)[i]->GetZ();
185 SetCurrentElement(Z);
186 G4double tmax = MaxSecondaryEnergy(particle, kineticEnergy);
187 G4double loss = ComputMuPairLoss(Z, kineticEnergy, cutEnergy, tmax);
188 dedx += loss*theAtomicNumDensityVector[i];
189 }
190 if (dedx < 0.) dedx = 0.;
191 return dedx;
192}
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195
196G4double G4MuPairProductionModel::ComputMuPairLoss(G4double Z,
197 G4double tkin,
198 G4double cutEnergy,
199 G4double tmax)
200{
201 SetCurrentElement(Z);
202 G4double loss = 0.0;
203
204 G4double cut = std::min(cutEnergy,tmax);
205 if(cut <= minPairEnergy) return loss;
206
207 // calculate the rectricted loss
208 // numerical integration in log(PairEnergy)
209 G4double ak1=6.9;
210 G4double ak2=1.0;
211 G4double aaa = log(minPairEnergy);
212 G4double bbb = log(cut);
213 G4int kkk = (G4int)((bbb-aaa)/ak1+ak2);
214 if (kkk > 8) kkk = 8;
215 G4double hhh = (bbb-aaa)/(G4double)kkk;
216 G4double x = aaa;
217
218 for (G4int l=0 ; l<kkk; l++)
219 {
220
221 for (G4int ll=0; ll<8; ll++)
222 {
223 G4double ep = exp(x+xgi[ll]*hhh);
224 loss += wgi[ll]*ep*ep*ComputeDMicroscopicCrossSection(tkin, Z, ep);
225 }
226 x += hhh;
227 }
228 loss *= hhh;
229 if (loss < 0.) loss = 0.;
230 return loss;
231}
232
233//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
234
235G4double G4MuPairProductionModel::ComputeMicroscopicCrossSection(
236 G4double tkin,
237 G4double Z,
238 G4double cut)
239{
240 G4double cross = 0.;
241 SetCurrentElement(Z);
242 G4double tmax = MaxSecondaryEnergy(particle, tkin);
243 if (tmax <= cut) return cross;
244
245 G4double ak1=6.9 ;
246 G4double ak2=1.0 ;
247 G4double aaa = log(cut);
248 G4double bbb = log(tmax);
249 G4int kkk = (G4int)((bbb-aaa)/ak1 + ak2);
250 if(kkk > 8) kkk = 8;
251 G4double hhh = (bbb-aaa)/float(kkk);
252 G4double x = aaa;
253
254 for(G4int l=0; l<kkk; l++)
255 {
256 for(G4int i=0; i<8; i++)
257 {
258 G4double ep = exp(x + xgi[i]*hhh);
259 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
260 }
261 x += hhh;
262 }
263
264 cross *=hhh;
265 if(cross < 0.0) cross = 0.0;
266 return cross;
267}
268
269G4double G4MuPairProductionModel::ComputeDMicroscopicCrossSection(
270 G4double tkin,
271 G4double Z,
272 G4double pairEnergy)
273 // Calculates the differential (D) microscopic cross section
274 // using the cross section formula of R.P. Kokoulin (18/01/98)
275 // Code modified by R.P. Kokoulin, V.N. Ivanchenko (27/01/04)
276{
277 G4double bbbtf= 183. ;
278 G4double bbbh = 202.4 ;
279 G4double g1tf = 1.95e-5 ;
280 G4double g2tf = 5.3e-5 ;
281 G4double g1h = 4.4e-5 ;
282 G4double g2h = 4.8e-5 ;
283
284 G4double totalEnergy = tkin + particleMass;
285 G4double residEnergy = totalEnergy - pairEnergy;
286 G4double massratio = particleMass/electron_mass_c2 ;
287 G4double massratio2 = massratio*massratio ;
288 G4double cross = 0.;
289
290 SetCurrentElement(Z);
291
292 G4double c3 = 0.75*sqrte*particleMass;
293 if (residEnergy <= c3*z13) return cross;
294
295 G4double c7 = 4.*electron_mass_c2;
296 G4double c8 = 6.*particleMass*particleMass;
297 G4double alf = c7/pairEnergy;
298 G4double a3 = 1. - alf;
299 if (a3 <= 0.) return cross;
300
301 // zeta calculation
302 G4double bbb,g1,g2;
303 if( Z < 1.5 ) { bbb = bbbh ; g1 = g1h ; g2 = g2h ; }
304 else { bbb = bbbtf; g1 = g1tf; g2 = g2tf; }
305
306 G4double zeta = 0;
307 G4double zeta1 = 0.073*log(totalEnergy/(particleMass+g1*z23*totalEnergy))-0.26;
308 if ( zeta1 > 0.)
309 {
310 G4double zeta2 = 0.058*log(totalEnergy/(particleMass+g2*z13*totalEnergy))-0.14;
311 zeta = zeta1/zeta2 ;
312 }
313
314 G4double z2 = Z*(Z+zeta);
315 G4double screen0 = 2.*electron_mass_c2*sqrte*bbb/(z13*pairEnergy);
316 G4double a0 = totalEnergy*residEnergy;
317 G4double a1 = pairEnergy*pairEnergy/a0;
318 G4double bet = 0.5*a1;
319 G4double xi0 = 0.25*massratio2*a1;
320 G4double del = c8/a0;
321
322 G4double rta3 = sqrt(a3);
323 G4double tmnexp = alf/(1. + rta3) + del*rta3;
324 if(tmnexp >= 1.0) return cross;
325
326 G4double tmn = log(tmnexp);
327 G4double sum = 0.;
328
329 // Gaussian integration in ln(1-ro) ( with 8 points)
330 for (G4int i=0; i<8; i++)
331 {
332 G4double a4 = exp(tmn*xgi[i]); // a4 = (1.-asymmetry)
333 G4double a5 = a4*(2.-a4) ;
334 G4double a6 = 1.-a5 ;
335 G4double a7 = 1.+a6 ;
336 G4double a9 = 3.+a6 ;
337 G4double xi = xi0*a5 ;
338 G4double xii = 1./xi ;
339 G4double xi1 = 1.+xi ;
340 G4double screen = screen0*xi1/a5 ;
341 G4double yeu = 5.-a6+4.*bet*a7 ;
342 G4double yed = 2.*(1.+3.*bet)*log(3.+xii)-a6-a1*(2.-a6) ;
343 G4double ye1 = 1.+yeu/yed ;
344 G4double ale=log(bbb/z13*sqrt(xi1*ye1)/(1.+screen*ye1)) ;
345 G4double cre = 0.5*log(1.+2.25*z23*xi1*ye1/massratio2) ;
346 G4double be;
347
348 if (xi <= 1.e3) be = ((2.+a6)*(1.+bet)+xi*a9)*log(1.+xii)+(a5-bet)/xi1-a9;
349 else be = (3.-a6+a1*a7)/(2.*xi);
350
351 G4double fe = (ale-cre)*be;
352 if ( fe < 0.) fe = 0. ;
353
354 G4double ymu = 4.+a6 +3.*bet*a7 ;
355 G4double ymd = a7*(1.5+a1)*log(3.+xi)+1.-1.5*a6 ;
356 G4double ym1 = 1.+ymu/ymd ;
357 G4double alm_crm = log(bbb*massratio/(1.5*z23*(1.+screen*ym1)));
358 G4double a10,bm;
359 if ( xi >= 1.e-3)
360 {
361 a10 = (1.+a1)*a5 ;
362 bm = (a7*(1.+1.5*bet)-a10*xii)*log(xi1)+xi*(a5-bet)/xi1+a10;
363 } else {
364 bm = (5.-a6+bet*a9)*(xi/2.);
365 }
366
367 G4double fm = alm_crm*bm;
368 if ( fm < 0.) fm = 0. ;
369
370 sum += wgi[i]*a4*(fe+fm/massratio2);
371 }
372
373 cross = -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy);
374
375 return cross;
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
379
380G4double G4MuPairProductionModel::ComputeCrossSectionPerAtom(
381 const G4ParticleDefinition*,
382 G4double kineticEnergy,
383 G4double Z, G4double,
384 G4double cutEnergy,
385 G4double maxEnergy)
386{
387 G4double cross = 0.0;
388 if (kineticEnergy <= lowestKinEnergy) return cross;
389
390 SetCurrentElement(Z);
391
392 G4double maxPairEnergy = MaxSecondaryEnergy(particle,kineticEnergy);
393 G4double tmax = std::min(maxEnergy, maxPairEnergy);
394 G4double cut = std::max(cutEnergy, minPairEnergy);
395 if (cut >= tmax) return cross;
396
397 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
398 if(tmax < kineticEnergy) {
399 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
400 }
401 return cross;
402}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
405
406void G4MuPairProductionModel::MakeSamplingTables()
407{
408 for (G4int iz=0; iz<nzdat; iz++)
409 {
410 G4double Z = zdat[iz];
411 SetCurrentElement(Z);
412
413 for (G4int it=0; it<ntdat; it++) {
414
415 G4double kineticEnergy = tdat[it];
416 G4double maxPairEnergy = MaxSecondaryEnergy(particle,kineticEnergy);
417 // G4cout << "Z= " << currentZ << " z13= " << z13
418 //<< " mE= " << maxPairEnergy << G4endl;
419 G4double CrossSection = 0.0 ;
420
421 if(maxPairEnergy > minPairEnergy) {
422
423 G4double y = ymin - 0.5*dy ;
424 G4double yy = ymin - dy ;
425 G4double x = exp(y);
426 G4double fac = exp(dy);
427 G4double dx = exp(yy)*(fac - 1.0);
428
429 G4double c = log(maxPairEnergy/minPairEnergy);
430
431 for (G4int i=0 ; i<nbiny; i++) {
432 y += dy ;
433 if(c > 0.0) {
434 x *= fac;
435 dx*= fac;
436 G4double ep = minPairEnergy*exp(c*x) ;
437 CrossSection +=
438 ep*dx*ComputeDMicroscopicCrossSection(kineticEnergy, Z, ep);
439 }
440 ya[i] = y;
441 proba[iz][it][i] = CrossSection;
442 }
443
444 } else {
445 for (G4int i=0 ; i<nbiny; i++) {
446 proba[iz][it][i] = CrossSection;
447 }
448 }
449
450 ya[nbiny]=ymax;
451 proba[iz][it][nbiny] = CrossSection;
452
453 }
454 }
455 samplingTablesAreFilled = true;
456}
457
458//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
459
460void
461G4MuPairProductionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
462 const G4MaterialCutsCouple* couple,
463 const G4DynamicParticle* aDynamicParticle,
464 G4double tmin,
465 G4double tmax)
466{
467 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
468 G4double totalEnergy = kineticEnergy + particleMass;
469 G4double totalMomentum =
470 sqrt(kineticEnergy*(kineticEnergy + 2.0*particleMass));
471
472 G4ThreeVector partDirection = aDynamicParticle->GetMomentumDirection();
473
474 G4int it;
475 for(it=1; it<ntdat; ++it) {if(kineticEnergy <= tdat[it]) break;}
476 if(it == ntdat) { --it; }
477 G4double dt = log(kineticEnergy/tdat[it-1])/log(tdat[it]/tdat[it-1]);
478
479 // select randomly one element constituing the material
480 const G4Element* anElement =
481 SelectRandomAtom(kineticEnergy, dt, it, couple, tmin);
482 SetCurrentElement(anElement->GetZ());
483
484 // define interval of enegry transfer
485 G4double maxPairEnergy = MaxSecondaryEnergy(particle,kineticEnergy);
486 G4double maxEnergy = std::min(tmax, maxPairEnergy);
487 G4double minEnergy = std::max(tmin, minPairEnergy);
488
489 if(minEnergy >= maxEnergy) { return; }
490 //G4cout << "emin= " << minEnergy << " emax= " << maxEnergy
491 // << " minPair= " << minPairEnergy << " maxpair= " << maxPairEnergy
492 // << " ymin= " << ymin << " dy= " << dy << G4endl;
493
494 // select bins
495 G4int iymin = 0;
496 G4int iymax = nbiny-1;
497 if( minEnergy > minPairEnergy)
498 {
499 G4double xc = log(minEnergy/minPairEnergy)/log(maxPairEnergy/minPairEnergy);
500 iymin = (G4int)((log(xc) - ymin)/dy);
501 if(iymin >= nbiny) iymin = nbiny-1;
502 else if(iymin < 0) iymin = 0;
503 xc = log(maxEnergy/minPairEnergy)/log(maxPairEnergy/minPairEnergy);
504 iymax = (G4int)((log(xc) - ymin)/dy) + 1;
505 if(iymax >= nbiny) iymax = nbiny-1;
506 else if(iymax < 0) iymax = 0;
507 }
508
509 // sample e-e+ energy, pair energy first
510 G4int iz, iy;
511
512 for(iz=1; iz<nzdat; ++iz) { if(currentZ <= zdat[iz]) { break; } }
513 if(iz == nzdat) { --iz; }
514
515 G4double dz = log(currentZ/zdat[iz-1])/log(zdat[iz]/zdat[iz-1]);
516
517 G4double pmin = InterpolatedIntegralCrossSection(dt,dz,iz,it,iymin,currentZ);
518 G4double pmax = InterpolatedIntegralCrossSection(dt,dz,iz,it,iymax,currentZ);
519
520 G4double p = pmin+G4UniformRand()*(pmax - pmin);
521
522 // interpolate sampling vector;
523 G4double p1 = pmin;
524 G4double p2 = pmin;
525 for(iy=iymin+1; iy<=iymax; ++iy) {
526 p1 = p2;
527 p2 = InterpolatedIntegralCrossSection(dt, dz, iz, it, iy, currentZ);
528 if(p <= p2) break;
529 }
530 // G4cout << "iy= " << iy << " iymin= " << iymin << " iymax= "
531 // << iymax << " Z= " << currentZ << G4endl;
532 G4double y = ya[iy-1] + dy*(p - p1)/(p2 - p1);
533
534 G4double PairEnergy = minPairEnergy*exp(exp(y)
535 *log(maxPairEnergy/minPairEnergy));
536
537 if(PairEnergy < minEnergy) { PairEnergy = minEnergy; }
538 if(PairEnergy > maxEnergy) { PairEnergy = maxEnergy; }
539
540 // sample r=(E+-E-)/PairEnergy ( uniformly .....)
541 G4double rmax =
542 (1.-6.*particleMass*particleMass/(totalEnergy*(totalEnergy-PairEnergy)))
543 *sqrt(1.-minPairEnergy/PairEnergy);
544 G4double r = rmax * (-1.+2.*G4UniformRand()) ;
545
546 // compute energies from PairEnergy,r
547 G4double ElectronEnergy = (1.-r)*PairEnergy*0.5;
548 G4double PositronEnergy = PairEnergy - ElectronEnergy;
549
550 // The angle of the emitted virtual photon is sampled
551 // according to the muon bremsstrahlung model
552
553 G4double gam = totalEnergy/particleMass;
554 G4double gmax = gam*std::min(1.0, totalEnergy/PairEnergy - 1.0);
555 G4double gmax2= gmax*gmax;
556 G4double x = G4UniformRand()*gmax2/(1.0 + gmax2);
557
558 G4double theta = sqrt(x/(1.0 - x))/gam;
559 G4double sint = sin(theta);
560 G4double phi = twopi * G4UniformRand() ;
561 G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
562
563 G4ThreeVector gDirection(dirx, diry, dirz);
564 gDirection.rotateUz(partDirection);
565
566 // the angles of e- and e+ assumed to be the same as virtual gamma
567
568 // create G4DynamicParticle object for the particle1
569 G4DynamicParticle* aParticle1 =
570 new G4DynamicParticle(theElectron, gDirection,
571 ElectronEnergy - electron_mass_c2);
572
573 // create G4DynamicParticle object for the particle2
574 G4DynamicParticle* aParticle2 =
575 new G4DynamicParticle(thePositron, gDirection,
576 PositronEnergy - electron_mass_c2);
577
578 // primary change
579 kineticEnergy -= (ElectronEnergy + PositronEnergy);
580 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
581
582 partDirection *= totalMomentum;
583 partDirection -= (aParticle1->GetMomentum() + aParticle2->GetMomentum());
584 partDirection = partDirection.unit();
585 fParticleChange->SetProposedMomentumDirection(partDirection);
586
587 // add secondary
588 vdp->push_back(aParticle1);
589 vdp->push_back(aParticle2);
590}
591
592//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
593
594const G4Element* G4MuPairProductionModel::SelectRandomAtom(
595 G4double kinEnergy, G4double dt, G4int it,
596 const G4MaterialCutsCouple* couple, G4double tmin)
597{
598 // select randomly 1 element within the material
599
600 const G4Material* material = couple->GetMaterial();
601 size_t nElements = material->GetNumberOfElements();
602 const G4ElementVector* theElementVector = material->GetElementVector();
603 if (nElements == 1) return (*theElementVector)[0];
604
605 if(nElements > nmaxElements) {
606 nmaxElements = nElements;
607 partialSum.resize(nmaxElements);
608 }
609
610 const G4double* theAtomNumDensityVector=material->GetAtomicNumDensityVector();
611
612 G4double sum = 0.0;
613 G4double dl;
614
615 size_t i;
616 for (i=0; i<nElements; i++) {
617 G4double Z = ((*theElementVector)[i])->GetZ();
618 SetCurrentElement(Z);
619 G4double maxPairEnergy = MaxSecondaryEnergy(particle,kinEnergy);
620 G4double minEnergy = std::max(tmin, minPairEnergy);
621 dl = 0.0;
622 if(minEnergy < maxPairEnergy) {
623
624 G4int iz;
625 for(iz=1; iz<nzdat; iz++) {if(Z <= zdat[iz]) break;}
626 if(iz == nzdat) iz--;
627 G4double dz = log(Z/zdat[iz-1])/log(zdat[iz]/zdat[iz-1]);
628
629 G4double sigcut;
630 if(minEnergy <= minPairEnergy)
631 sigcut = 0.;
632 else
633 {
634 G4double xc = log(minEnergy/minPairEnergy)/log(maxPairEnergy/minPairEnergy);
635 G4int iy = (G4int)((log(xc) - ymin)/dy);
636 if(iy < 0) iy = 0;
637 if(iy >= nbiny) iy = nbiny-1;
638 sigcut = InterpolatedIntegralCrossSection(dt,dz,iz,it,iy, Z);
639 }
640
641 G4double sigtot = InterpolatedIntegralCrossSection(dt,dz,iz,it,nbiny,Z);
642 dl = (sigtot - sigcut)*theAtomNumDensityVector[i];
643 }
644 // protection
645 if(dl < 0.0) dl = 0.0;
646 sum += dl;
647 partialSum[i] = sum;
648 }
649
650 G4double rval = G4UniformRand()*sum;
651 for (i=0; i<nElements; i++) {
652 if(rval<=partialSum[i]) return (*theElementVector)[i];
653 }
654
655 return (*theElementVector)[nElements - 1];
656
657}
658
659//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
660
661
Note: See TracBrowser for help on using the repository browser.