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

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

update

File size: 21.4 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.39 2008/07/22 16:11:34 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
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
139void G4MuPairProductionModel::Initialise(const G4ParticleDefinition* p,
140                                         const G4DataVector&)
141{ 
142  if (!samplingTablesAreFilled) {
143    if(p) SetParticle(p);
144    MakeSamplingTables();
145  }
146  if(!fParticleChange) {
147    if(pParticleChange) 
148      fParticleChange = 
149        reinterpret_cast<G4ParticleChangeForLoss*>(pParticleChange);
150    else 
151      fParticleChange = new G4ParticleChangeForLoss();
152  }
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
156
157G4double G4MuPairProductionModel::ComputeDEDXPerVolume(
158                                              const G4Material* material,
159                                              const G4ParticleDefinition*,
160                                                    G4double kineticEnergy,
161                                                    G4double cutEnergy)
162{
163  G4double dedx = 0.0;
164  if (cutEnergy <= minPairEnergy || kineticEnergy <= lowestKinEnergy)
165    return dedx;
166
167  const G4ElementVector* theElementVector = material->GetElementVector();
168  const G4double* theAtomicNumDensityVector =
169                                   material->GetAtomicNumDensityVector();
170
171  //  loop for elements in the material
172  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
173     G4double Z = (*theElementVector)[i]->GetZ();
174     SetCurrentElement(Z);
175     G4double tmax = MaxSecondaryEnergy(particle, kineticEnergy);
176     G4double loss = ComputMuPairLoss(Z, kineticEnergy, cutEnergy, tmax);
177     dedx += loss*theAtomicNumDensityVector[i];
178  }
179  if (dedx < 0.) dedx = 0.;
180  return dedx;
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184
185G4double G4MuPairProductionModel::ComputMuPairLoss(G4double Z, 
186                                                   G4double tkin,
187                                                   G4double cutEnergy, 
188                                                   G4double tmax)
189{
190  SetCurrentElement(Z);
191  G4double loss = 0.0;
192
193  G4double cut = std::min(cutEnergy,tmax);
194  if(cut <= minPairEnergy) return loss;
195
196  // calculate the rectricted loss
197  // numerical integration in log(PairEnergy)
198  G4double ak1=6.9;
199  G4double ak2=1.0;
200  G4double aaa = log(minPairEnergy);
201  G4double bbb = log(cut);
202  G4int    kkk = (G4int)((bbb-aaa)/ak1+ak2);
203  if (kkk > 8) kkk = 8;
204  G4double hhh = (bbb-aaa)/(G4double)kkk;
205  G4double x = aaa;
206
207  for (G4int l=0 ; l<kkk; l++)
208  {
209
210    for (G4int ll=0; ll<8; ll++)
211    {
212      G4double ep = exp(x+xgi[ll]*hhh);
213      loss += wgi[ll]*ep*ep*ComputeDMicroscopicCrossSection(tkin, Z, ep);
214    }
215    x += hhh;
216  }
217  loss *= hhh;
218  if (loss < 0.) loss = 0.;
219  return loss;
220}
221
222//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223
224G4double G4MuPairProductionModel::ComputeMicroscopicCrossSection(
225                                           G4double tkin,
226                                           G4double Z,
227                                           G4double cut)
228{
229  G4double cross = 0.;
230  SetCurrentElement(Z);
231  G4double tmax = MaxSecondaryEnergy(particle, tkin);
232  if (tmax <= cut) return cross;
233
234  G4double ak1=6.9 ;
235  G4double ak2=1.0 ;
236  G4double aaa = log(cut);
237  G4double bbb = log(tmax);
238  G4int kkk = (G4int)((bbb-aaa)/ak1 + ak2);
239  if(kkk > 8) kkk = 8;
240  G4double hhh = (bbb-aaa)/float(kkk);
241  G4double x = aaa;
242
243  for(G4int l=0; l<kkk; l++)
244  {
245    for(G4int i=0; i<8; i++)
246    {
247      G4double ep = exp(x + xgi[i]*hhh);
248      cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
249    }
250    x += hhh;
251  }
252
253  cross *=hhh;
254  if(cross < 0.0) cross = 0.0;
255  return cross;
256}
257
258G4double G4MuPairProductionModel::ComputeDMicroscopicCrossSection(
259                                           G4double tkin,
260                                           G4double Z,
261                                           G4double pairEnergy)
262 // Calculates the  differential (D) microscopic cross section
263 // using the cross section formula of R.P. Kokoulin (18/01/98)
264 // Code modified by R.P. Kokoulin, V.N. Ivanchenko (27/01/04)
265{
266  G4double bbbtf= 183. ;
267  G4double bbbh = 202.4 ;
268  G4double g1tf = 1.95e-5 ;
269  G4double g2tf = 5.3e-5 ;
270  G4double g1h  = 4.4e-5 ;
271  G4double g2h  = 4.8e-5 ;
272
273  G4double totalEnergy  = tkin + particleMass;
274  G4double residEnergy  = totalEnergy - pairEnergy;
275  G4double massratio    = particleMass/electron_mass_c2 ;
276  G4double massratio2   = massratio*massratio ;
277  G4double cross = 0.;
278
279  SetCurrentElement(Z);
280
281  G4double c3 = 0.75*sqrte*particleMass;
282  if (residEnergy <= c3*z13) return cross;
283
284  G4double c7 = 4.*electron_mass_c2;
285  G4double c8 = 6.*particleMass*particleMass;
286  G4double alf = c7/pairEnergy;
287  G4double a3 = 1. - alf;
288  if (a3 <= 0.) return cross;
289
290  // zeta calculation
291  G4double bbb,g1,g2;
292  if( Z < 1.5 ) { bbb = bbbh ; g1 = g1h ; g2 = g2h ; }
293  else          { bbb = bbbtf; g1 = g1tf; g2 = g2tf; }
294
295  G4double zeta = 0;
296  G4double zeta1 = 0.073*log(totalEnergy/(particleMass+g1*z23*totalEnergy))-0.26;
297  if ( zeta1 > 0.)
298  {
299    G4double zeta2 = 0.058*log(totalEnergy/(particleMass+g2*z13*totalEnergy))-0.14;
300    zeta  = zeta1/zeta2 ;
301  }
302
303  G4double z2 = Z*(Z+zeta);
304  G4double screen0 = 2.*electron_mass_c2*sqrte*bbb/(z13*pairEnergy);
305  G4double a0 = totalEnergy*residEnergy;
306  G4double a1 = pairEnergy*pairEnergy/a0;
307  G4double bet = 0.5*a1;
308  G4double xi0 = 0.25*massratio2*a1;
309  G4double del = c8/a0;
310
311  G4double rta3 = sqrt(a3);
312  G4double tmnexp = alf/(1. + rta3) + del*rta3;
313  if(tmnexp >= 1.0) return cross;
314
315  G4double tmn = log(tmnexp);
316  G4double sum = 0.;
317
318  // Gaussian integration in ln(1-ro) ( with 8 points)
319  for (G4int i=0; i<8; i++)
320  {
321    G4double a4 = exp(tmn*xgi[i]);     // a4 = (1.-asymmetry)
322    G4double a5 = a4*(2.-a4) ;
323    G4double a6 = 1.-a5 ;
324    G4double a7 = 1.+a6 ;
325    G4double a9 = 3.+a6 ;
326    G4double xi = xi0*a5 ;
327    G4double xii = 1./xi ;
328    G4double xi1 = 1.+xi ;
329    G4double screen = screen0*xi1/a5 ;
330    G4double yeu = 5.-a6+4.*bet*a7 ;
331    G4double yed = 2.*(1.+3.*bet)*log(3.+xii)-a6-a1*(2.-a6) ;
332    G4double ye1 = 1.+yeu/yed ;
333    G4double ale=log(bbb/z13*sqrt(xi1*ye1)/(1.+screen*ye1)) ;
334    G4double cre = 0.5*log(1.+2.25*z23*xi1*ye1/massratio2) ;
335    G4double be;
336
337    if (xi <= 1.e3) be = ((2.+a6)*(1.+bet)+xi*a9)*log(1.+xii)+(a5-bet)/xi1-a9;
338    else            be = (3.-a6+a1*a7)/(2.*xi);
339
340    G4double fe = (ale-cre)*be;
341    if ( fe < 0.) fe = 0. ;
342
343    G4double ymu = 4.+a6 +3.*bet*a7 ;
344    G4double ymd = a7*(1.5+a1)*log(3.+xi)+1.-1.5*a6 ;
345    G4double ym1 = 1.+ymu/ymd ;
346    G4double alm_crm = log(bbb*massratio/(1.5*z23*(1.+screen*ym1)));
347    G4double a10,bm;
348    if ( xi >= 1.e-3)
349    {
350      a10 = (1.+a1)*a5 ;
351      bm  = (a7*(1.+1.5*bet)-a10*xii)*log(xi1)+xi*(a5-bet)/xi1+a10;
352    } else {
353      bm = (5.-a6+bet*a9)*(xi/2.);
354    }
355
356    G4double fm = alm_crm*bm;
357    if ( fm < 0.) fm = 0. ;
358
359    sum += wgi[i]*a4*(fe+fm/massratio2);
360  }
361
362  cross = -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy);
363
364  return cross;
365}
366
367//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
368
369G4double G4MuPairProductionModel::ComputeCrossSectionPerAtom(
370                                           const G4ParticleDefinition*,
371                                                 G4double kineticEnergy,
372                                                 G4double Z, G4double,
373                                                 G4double cutEnergy,
374                                                 G4double maxEnergy)
375{
376  G4double cross = 0.0;
377  if (kineticEnergy <= lowestKinEnergy) return cross;
378
379  SetCurrentElement(Z);
380  G4double tmax = std::min(maxEnergy, kineticEnergy);
381  G4double cut  = std::min(cutEnergy, kineticEnergy);
382  if(cut < minPairEnergy) cut = minPairEnergy;
383  if (cut >= tmax) return cross;
384
385  cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
386  if(tmax < kineticEnergy) {
387    cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
388  }
389  return cross;
390}
391
392//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
393
394void G4MuPairProductionModel::MakeSamplingTables()
395{
396  for (G4int iz=0; iz<nzdat; iz++)
397  {
398    G4double Z = zdat[iz];
399    SetCurrentElement(Z);
400
401    for (G4int it=0; it<ntdat; it++) {
402
403      G4double kineticEnergy = tdat[it];
404      G4double maxPairEnergy = MaxSecondaryEnergy(particle,kineticEnergy);
405      // G4cout << "Z= " << currentZ << " z13= " << z13
406      //<< " mE= " << maxPairEnergy << G4endl;
407      G4double CrossSection = 0.0 ;
408
409      if(maxPairEnergy > minPairEnergy) {
410
411        G4double y = ymin - 0.5*dy ;
412        G4double yy = ymin - dy ;
413        G4double x = exp(y);
414        G4double fac = exp(dy);
415        G4double dx = exp(yy)*(fac - 1.0);
416
417        G4double c = log(maxPairEnergy/minPairEnergy);
418
419        for (G4int i=0 ; i<nbiny; i++) {
420          y += dy ;
421          if(c > 0.0) {
422            x *= fac;
423            dx*= fac;
424            G4double ep = minPairEnergy*exp(c*x) ;
425            CrossSection += 
426              ep*dx*ComputeDMicroscopicCrossSection(kineticEnergy, Z, ep);
427          }
428          ya[i] = y;
429          proba[iz][it][i] = CrossSection;
430        }
431       
432      } else {
433        for (G4int i=0 ; i<nbiny; i++) {
434          proba[iz][it][i] = CrossSection;
435        }
436      }
437
438      ya[nbiny]=ymax;
439      proba[iz][it][nbiny] = CrossSection;
440
441    }
442  }
443  samplingTablesAreFilled = true;
444}
445
446//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
447
448void G4MuPairProductionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp, 
449                                                const G4MaterialCutsCouple* couple,
450                                                const G4DynamicParticle* aDynamicParticle,
451                                                G4double tmin,
452                                                G4double tmax)
453{
454  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
455  G4double totalEnergy   = kineticEnergy + particleMass ;
456  G4ParticleMomentum ParticleDirection = 
457    aDynamicParticle->GetMomentumDirection();
458
459  G4int it;
460  for(it=1; it<ntdat; it++) {if(kineticEnergy <= tdat[it]) break;}
461  if(it == ntdat) it--;
462  G4double dt = log(kineticEnergy/tdat[it-1])/log(tdat[it]/tdat[it-1]);
463
464  // select randomly one element constituing the material
465  const G4Element* anElement = SelectRandomAtom(kineticEnergy, dt, it, couple, tmin);
466  SetCurrentElement(anElement->GetZ());
467
468  // define interval of enegry transfer
469  G4double maxPairEnergy = MaxSecondaryEnergy(particle,kineticEnergy);
470  G4double maxEnergy     = std::min(tmax, maxPairEnergy);
471  G4double minEnergy     = std::max(tmin, minPairEnergy);
472
473  if(minEnergy >= maxEnergy) return;
474  //G4cout << "emin= " << minEnergy << " emax= " << maxEnergy
475  //     << " minPair= " << minPairEnergy << " maxpair= " << maxPairEnergy
476  //       << " ymin= " << ymin << " dy= " << dy << G4endl;
477
478  // select bins
479  G4int iymin = 0;
480  G4int iymax = nbiny-1;
481  if( minEnergy > minPairEnergy)
482  {
483    G4double xc = log(minEnergy/minPairEnergy)/log(maxPairEnergy/minPairEnergy);
484    iymin = (G4int)((log(xc) - ymin)/dy);
485    if(iymin >= nbiny) iymin = nbiny-1;
486    else if(iymin < 0) iymin = 0;
487    xc = log(maxEnergy/minPairEnergy)/log(maxPairEnergy/minPairEnergy);
488    iymax = (G4int)((log(xc) - ymin)/dy) + 1;
489    if(iymax >= nbiny) iymax = nbiny-1;
490    else if(iymax < 0) iymax = 0;
491  }
492
493  // sample e-e+ energy, pair energy first
494  G4int iz, iy;
495
496  for(iz=1; iz<nzdat; iz++) {if(currentZ <= zdat[iz]) break;}
497  if(iz == nzdat) iz--;
498
499  G4double dz = log(currentZ/zdat[iz-1])/log(zdat[iz]/zdat[iz-1]);
500
501  G4double pmin = InterpolatedIntegralCrossSection(dt,dz,iz,it,iymin,currentZ);
502  G4double pmax = InterpolatedIntegralCrossSection(dt,dz,iz,it,iymax,currentZ);
503
504  G4double p = pmin+G4UniformRand()*(pmax - pmin);
505
506  // interpolate sampling vector;
507  G4double p1 = pmin;
508  G4double p2 = pmin;
509  for(iy=iymin+1; iy<=iymax; iy++) {
510    p1 = p2;
511    p2 = InterpolatedIntegralCrossSection(dt, dz, iz, it, iy, currentZ);
512    if(p <= p2) break;
513  }
514  // G4cout << "iy= " << iy << " iymin= " << iymin << " iymax= "
515  //        << iymax << " Z= " << currentZ << G4endl;
516  G4double y = ya[iy-1] + dy*(p - p1)/(p2 - p1);
517
518  G4double PairEnergy = minPairEnergy*exp(exp(y)
519                       *log(maxPairEnergy/minPairEnergy));
520                       
521  if(PairEnergy < minEnergy) PairEnergy = minEnergy;
522  if(PairEnergy > maxEnergy) PairEnergy = maxEnergy;
523
524  // sample r=(E+-E-)/PairEnergy  ( uniformly .....)
525  G4double rmax =
526    (1.-6.*particleMass*particleMass/(totalEnergy*(totalEnergy-PairEnergy)))
527                                       *sqrt(1.-minPairEnergy/PairEnergy);
528  G4double r = rmax * (-1.+2.*G4UniformRand()) ;
529
530  // compute energies from PairEnergy,r
531  G4double ElectronEnergy = (1.-r)*PairEnergy*0.5;
532  G4double PositronEnergy = PairEnergy - ElectronEnergy;
533
534  //  angles of the emitted particles ( Z - axis along the parent particle)
535  //      (mean theta for the moment)
536
537  //
538  // scattered electron (positron) angles. ( Z - axis along the parent photon)
539  //
540  //  universal distribution suggested by L. Urban
541  // (Geant3 manual (1993) Phys211),
542  //  derived from Tsai distribution (Rev Mod Phys 49,421(1977))
543  //  G4cout << "Ee= " << ElectronEnergy << " Ep= " << PositronEnergy << G4endl;
544  G4double u;
545  const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
546
547  if (9./(9.+d) >G4UniformRand()) u= - log(G4UniformRand()*G4UniformRand())/a1;
548  else                            u= - log(G4UniformRand()*G4UniformRand())/a2;
549
550  G4double TetEl = u*electron_mass_c2/ElectronEnergy;
551  G4double TetPo = u*electron_mass_c2/PositronEnergy;
552  G4double Phi  = twopi * G4UniformRand();
553  G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
554  G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
555
556  G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
557  ElectDirection.rotateUz(ParticleDirection);
558
559  // create G4DynamicParticle object for the particle1
560  G4DynamicParticle* aParticle1= new G4DynamicParticle(theElectron,
561                                                       ElectDirection,
562                                             ElectronEnergy - electron_mass_c2);
563
564  G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
565  PositDirection.rotateUz(ParticleDirection);
566
567  // create G4DynamicParticle object for the particle2
568  G4DynamicParticle* aParticle2 = 
569    new G4DynamicParticle(thePositron,
570                          PositDirection,
571                          PositronEnergy - electron_mass_c2);
572
573  // primary change
574  kineticEnergy -= (ElectronEnergy + PositronEnergy);
575  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
576
577  vdp->push_back(aParticle1);
578  vdp->push_back(aParticle2);
579}
580
581//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
582
583const G4Element* G4MuPairProductionModel::SelectRandomAtom(
584                 G4double kinEnergy, G4double dt, G4int it,
585           const G4MaterialCutsCouple* couple, G4double tmin)
586{
587  // select randomly 1 element within the material
588
589  const G4Material* material = couple->GetMaterial();
590  size_t nElements = material->GetNumberOfElements();
591  const G4ElementVector* theElementVector = material->GetElementVector();
592  if (nElements == 1) return (*theElementVector)[0];
593
594  if(nElements > nmaxElements) {
595    nmaxElements = nElements;
596    partialSum.resize(nmaxElements);
597  }
598
599  const G4double* theAtomNumDensityVector=material->GetAtomicNumDensityVector();
600
601  G4double sum = 0.0;
602
603  size_t i;
604  for (i=0; i<nElements; i++) {
605    G4double Z = ((*theElementVector)[i])->GetZ();
606    SetCurrentElement(Z);
607    G4double maxPairEnergy = MaxSecondaryEnergy(particle,kinEnergy);
608    G4double minEnergy     = std::max(tmin, minPairEnergy);
609
610    G4int iz;
611    for(iz=1; iz<nzdat; iz++) {if(Z <= zdat[iz]) break;}
612    if(iz == nzdat) iz--;
613    G4double dz = log(Z/zdat[iz-1])/log(zdat[iz]/zdat[iz-1]);
614
615    G4double sigcut;
616    if(minEnergy <= minPairEnergy)
617      sigcut = 0.;
618    else
619    {
620      G4double xc = log(minEnergy/minPairEnergy)/log(maxPairEnergy/minPairEnergy);
621      G4int iy = (G4int)((log(xc) - ymin)/dy);
622      if(iy < 0) iy = 0;
623      if(iy >= nbiny) iy = nbiny-1;
624      sigcut = InterpolatedIntegralCrossSection(dt,dz,iz,it,iy,   Z);
625    }
626
627    G4double sigtot = InterpolatedIntegralCrossSection(dt,dz,iz,it,nbiny,Z);
628    G4double dl = (sigtot - sigcut)*theAtomNumDensityVector[i];
629
630    // protection
631    if(dl < 0.0) dl = 0.0;
632    sum += dl;
633    partialSum[i] = sum;
634  }
635
636  G4double rval = G4UniformRand()*sum;
637  for (i=0; i<nElements; i++) {
638    if(rval<=partialSum[i]) return (*theElementVector)[i];
639  }
640
641  return (*theElementVector)[nElements - 1];
642
643}
644
645//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
646
647
Note: See TracBrowser for help on using the repository browser.