Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MuPairProductionModel.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: G4MuPairProductionModel.cc 74544 2013-10-14 12:40:29Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4MuPairProductionModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 24.06.2002
38 //
39 // Modifications:
40 //
41 // 04-12-02 Change G4DynamicParticle constructor in PostStep (V.Ivanchenko)
42 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
43 // 24-01-03 Fix for compounds (V.Ivanchenko)
44 // 27-01-03 Make models region aware (V.Ivanchenko)
45 // 13-02-03 Add model (V.Ivanchenko)
46 // 06-06-03 Fix in cross section calculation for high energy (V.Ivanchenko)
47 // 20-10-03 2*xi in ComputeDDMicroscopicCrossSection (R.Kokoulin)
48 // 8 integration points in ComputeDMicroscopicCrossSection
49 // 12-01-04 Take min cut of e- and e+ not its sum (V.Ivanchenko)
50 // 10-02-04 Update parameterisation using R.Kokoulin model (V.Ivanchenko)
51 // 28-04-04 For complex materials repeat calculation of max energy for each
52 // material (V.Ivanchenko)
53 // 01-11-04 Fix bug inside ComputeDMicroscopicCrossSection (R.Kokoulin)
54 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
55 // 03-08-05 Add SetParticle method (V.Ivantchenko)
56 // 23-10-05 Add protection in sampling of e+e- pair energy needed for
57 // low cuts (V.Ivantchenko)
58 // 13-02-06 Add ComputeCrossSectionPerAtom (mma)
59 // 24-04-07 Add protection in SelectRandomAtom method (V.Ivantchenko)
60 // 12-05-06 Updated sampling (use cut) in SelectRandomAtom (A.Bogdanov)
61 // 11-10-07 Add ignoreCut flag (V.Ivanchenko)
62 
63 //
64 // Class Description:
65 //
66 //
67 // -------------------------------------------------------------------
68 //
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 
73 #include "G4PhysicalConstants.hh"
74 #include "G4SystemOfUnits.hh"
75 #include "G4Electron.hh"
76 #include "G4Positron.hh"
77 #include "G4MuonMinus.hh"
78 #include "G4MuonPlus.hh"
79 #include "Randomize.hh"
80 #include "G4Material.hh"
81 #include "G4Element.hh"
82 #include "G4ElementVector.hh"
83 #include "G4ProductionCutsTable.hh"
86 #include "G4Log.hh"
87 #include "G4Exp.hh"
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90 
91 // static members
92 //
93 const G4int G4MuPairProductionModel::zdat[]={1, 4, 13, 29, 92};
94 const G4double G4MuPairProductionModel::adat[]={1.01, 9.01,26.98, 63.55, 238.03};
95 const G4double G4MuPairProductionModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083,
96  0.5917, 0.7628, 0.8983, 0.9801 };
97 const G4double G4MuPairProductionModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813,
98  0.1813, 0.1569, 0.1112, 0.0506 };
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101 
102 using namespace std;
103 
105  const G4String& nam)
106  : G4VEmModel(nam),
107  particle(0),
110  sqrte(sqrt(G4Exp(1.))),
111  currentZ(0),
112  fParticleChange(0),
113  minPairEnergy(4.*electron_mass_c2),
114  lowestKinEnergy(1.0*GeV),
115  nzdat(5),
116  nYBinPerDecade(4),
117  nbiny(1000),
118  nbine(0),
119  ymin(-5.),
120  dy(0.005)
121 {
123 
124  theElectron = G4Electron::Electron();
125  thePositron = G4Positron::Positron();
126 
127  particleMass = lnZ = z13 = z23 = 0;
128 
129  // setup lowest limit dependent on particle mass
130  if(p) {
131  SetParticle(p);
132  G4double limit = p->GetPDGMass()*8;
133  if(limit > lowestKinEnergy) { lowestKinEnergy = limit; }
134  }
135  emin = lowestKinEnergy;
136  emax = 10*TeV;
137 }
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140 
142 {}
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145 
147  const G4ParticleDefinition*,
148  G4double cut)
149 {
150  return std::max(lowestKinEnergy,cut);
151 }
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
154 
156  const G4DataVector& cuts)
157 {
158  SetParticle(p);
159 
160  // define scale of internal table for each thread only once
161  if(0 == nbine) {
162  emax = HighEnergyLimit();
163  emin = std::max(lowestKinEnergy, LowEnergyLimit());
164  nbine = size_t(nYBinPerDecade*std::log10(emax/emin));
165  if(nbine < 3) { nbine = 3; }
166 
167  ymin = G4Log(minPairEnergy/emin);
168  dy = -ymin/G4double(nbiny);
169  }
170 
171  if(IsMaster() && p == particle) {
172 
173  if(!fElementData) {
174  fElementData = new G4ElementData();
175  MakeSamplingTables();
176  }
177  InitialiseElementSelectors(p, cuts);
178  }
179 
180  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184 
186  G4VEmModel* masterModel)
187 {
188  if(p == particle) {
190  fElementData = masterModel->GetElementData();
191  }
192 }
193 
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195 
197  const G4Material* material,
198  const G4ParticleDefinition*,
199  G4double kineticEnergy,
200  G4double cutEnergy)
201 {
202  G4double dedx = 0.0;
203  if (cutEnergy <= minPairEnergy || kineticEnergy <= lowestKinEnergy)
204  { return dedx; }
205 
206  const G4ElementVector* theElementVector = material->GetElementVector();
207  const G4double* theAtomicNumDensityVector =
208  material->GetAtomicNumDensityVector();
209 
210  // loop for elements in the material
211  for (size_t i=0; i<material->GetNumberOfElements(); ++i) {
212  G4double Z = (*theElementVector)[i]->GetZ();
213  G4double tmax = MaxSecondaryEnergyForElement(kineticEnergy, Z);
214  G4double loss = ComputMuPairLoss(Z, kineticEnergy, cutEnergy, tmax);
215  dedx += loss*theAtomicNumDensityVector[i];
216  }
217  if (dedx < 0.) { dedx = 0.; }
218  return dedx;
219 }
220 
221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222 
224  G4double tkin,
225  G4double cutEnergy,
226  G4double tmax)
227 {
228  G4double loss = 0.0;
229 
230  G4double cut = std::min(cutEnergy,tmax);
231  if(cut <= minPairEnergy) { return loss; }
232 
233  // calculate the rectricted loss
234  // numerical integration in log(PairEnergy)
235  G4double ak1=6.9;
236  G4double ak2=1.0;
237  G4double aaa = G4Log(minPairEnergy);
238  G4double bbb = G4Log(cut);
239 
240  G4int kkk = (G4int)((bbb-aaa)/ak1+ak2);
241  if (kkk > 8) { kkk = 8; }
242  else if (kkk < 1) { kkk = 1; }
243 
244  G4double hhh = (bbb-aaa)/(G4double)kkk;
245  G4double x = aaa;
246 
247  for (G4int l=0 ; l<kkk; l++)
248  {
249 
250  for (G4int ll=0; ll<8; ll++)
251  {
252  G4double ep = G4Exp(x+xgi[ll]*hhh);
253  loss += wgi[ll]*ep*ep*ComputeDMicroscopicCrossSection(tkin, Z, ep);
254  }
255  x += hhh;
256  }
257  loss *= hhh;
258  if (loss < 0.) loss = 0.;
259  return loss;
260 }
261 
262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
263 
265  G4double tkin,
266  G4double Z,
267  G4double cutEnergy)
268 {
269  G4double cross = 0.;
270  G4double tmax = MaxSecondaryEnergyForElement(tkin, Z);
271  G4double cut = std::max(cutEnergy, minPairEnergy);
272  if (tmax <= cut) { return cross; }
273 
274  G4double ak1=6.9 ;
275  G4double ak2=1.0 ;
276  G4double aaa = G4Log(cut);
277  G4double bbb = G4Log(tmax);
278  G4int kkk = (G4int)((bbb-aaa)/ak1 + ak2);
279  if(kkk > 8) { kkk = 8; }
280  else if (kkk < 1) { kkk = 1; }
281 
282  G4double hhh = (bbb-aaa)/G4double(kkk);
283  G4double x = aaa;
284 
285  for(G4int l=0; l<kkk; ++l)
286  {
287  for(G4int i=0; i<8; ++i)
288  {
289  G4double ep = G4Exp(x + xgi[i]*hhh);
290  cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
291  }
292  x += hhh;
293  }
294 
295  cross *= hhh;
296  if(cross < 0.0) { cross = 0.0; }
297  return cross;
298 }
299 
300 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
301 
303  G4double tkin,
304  G4double Z,
305  G4double pairEnergy)
306  // Calculates the differential (D) microscopic cross section
307  // using the cross section formula of R.P. Kokoulin (18/01/98)
308  // Code modified by R.P. Kokoulin, V.N. Ivanchenko (27/01/04)
309 {
310  static const G4double bbbtf= 183. ;
311  static const G4double bbbh = 202.4 ;
312  static const G4double g1tf = 1.95e-5 ;
313  static const G4double g2tf = 5.3e-5 ;
314  static const G4double g1h = 4.4e-5 ;
315  static const G4double g2h = 4.8e-5 ;
316 
317  G4double totalEnergy = tkin + particleMass;
318  G4double residEnergy = totalEnergy - pairEnergy;
319  G4double massratio = particleMass/electron_mass_c2 ;
320  G4double massratio2 = massratio*massratio ;
321  G4double cross = 0.;
322 
323  G4double c3 = 0.75*sqrte*particleMass;
324  if (residEnergy <= c3*z13) { return cross; }
325 
326  G4double c7 = 4.*CLHEP::electron_mass_c2;
327  G4double c8 = 6.*particleMass*particleMass;
328  G4double alf = c7/pairEnergy;
329  G4double a3 = 1. - alf;
330  if (a3 <= 0.) { return cross; }
331 
332  // zeta calculation
333  G4double bbb,g1,g2;
334  if( Z < 1.5 ) { bbb = bbbh ; g1 = g1h ; g2 = g2h ; }
335  else { bbb = bbbtf; g1 = g1tf; g2 = g2tf; }
336 
337  G4double zeta = 0;
338  G4double zeta1 =
339  0.073*G4Log(totalEnergy/(particleMass+g1*z23*totalEnergy))-0.26;
340  if ( zeta1 > 0.)
341  {
342  G4double zeta2 =
343  0.058*G4Log(totalEnergy/(particleMass+g2*z13*totalEnergy))-0.14;
344  zeta = zeta1/zeta2 ;
345  }
346 
347  G4double z2 = Z*(Z+zeta);
348  G4double screen0 = 2.*electron_mass_c2*sqrte*bbb/(z13*pairEnergy);
349  G4double a0 = totalEnergy*residEnergy;
350  G4double a1 = pairEnergy*pairEnergy/a0;
351  G4double bet = 0.5*a1;
352  G4double xi0 = 0.25*massratio2*a1;
353  G4double del = c8/a0;
354 
355  G4double rta3 = sqrt(a3);
356  G4double tmnexp = alf/(1. + rta3) + del*rta3;
357  if(tmnexp >= 1.0) { return cross; }
358 
359  G4double tmn = G4Log(tmnexp);
360  G4double sum = 0.;
361 
362  // Gaussian integration in ln(1-ro) ( with 8 points)
363  for (G4int i=0; i<8; ++i)
364  {
365  G4double a4 = G4Exp(tmn*xgi[i]); // a4 = (1.-asymmetry)
366  G4double a5 = a4*(2.-a4) ;
367  G4double a6 = 1.-a5 ;
368  G4double a7 = 1.+a6 ;
369  G4double a9 = 3.+a6 ;
370  G4double xi = xi0*a5 ;
371  G4double xii = 1./xi ;
372  G4double xi1 = 1.+xi ;
373  G4double screen = screen0*xi1/a5 ;
374  G4double yeu = 5.-a6+4.*bet*a7 ;
375  G4double yed = 2.*(1.+3.*bet)*G4Log(3.+xii)-a6-a1*(2.-a6) ;
376  G4double ye1 = 1.+yeu/yed ;
377  G4double ale = G4Log(bbb/z13*sqrt(xi1*ye1)/(1.+screen*ye1)) ;
378  G4double cre = 0.5*G4Log(1.+2.25*z23*xi1*ye1/massratio2) ;
379  G4double be;
380 
381  if (xi <= 1.e3) {
382  be = ((2.+a6)*(1.+bet)+xi*a9)*G4Log(1.+xii)+(a5-bet)/xi1-a9;
383  } else {
384  be = (3.-a6+a1*a7)/(2.*xi);
385  }
386  G4double fe = (ale-cre)*be;
387  if ( fe < 0.) fe = 0. ;
388 
389  G4double ymu = 4.+a6 +3.*bet*a7 ;
390  G4double ymd = a7*(1.5+a1)*G4Log(3.+xi)+1.-1.5*a6 ;
391  G4double ym1 = 1.+ymu/ymd ;
392  G4double alm_crm = G4Log(bbb*massratio/(1.5*z23*(1.+screen*ym1)));
393  G4double a10,bm;
394  if ( xi >= 1.e-3)
395  {
396  a10 = (1.+a1)*a5 ;
397  bm = (a7*(1.+1.5*bet)-a10*xii)*G4Log(xi1)+xi*(a5-bet)/xi1+a10;
398  } else {
399  bm = (5.-a6+bet*a9)*(xi/2.);
400  }
401 
402  G4double fm = alm_crm*bm;
403  if ( fm < 0.) { fm = 0.; }
404 
405  sum += wgi[i]*a4*(fe+fm/massratio2);
406  }
407 
408  cross = -tmn*sum*factorForCross*z2*residEnergy/(totalEnergy*pairEnergy);
409  if(cross < 0.0) { cross = 0.0; }
410  return cross;
411 }
412 
413 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
414 
416  const G4ParticleDefinition*,
417  G4double kineticEnergy,
418  G4double Z, G4double,
419  G4double cutEnergy,
420  G4double maxEnergy)
421 {
422  G4double cross = 0.0;
423  if (kineticEnergy <= lowestKinEnergy) { return cross; }
424 
425  G4double maxPairEnergy = MaxSecondaryEnergyForElement(kineticEnergy, Z);
426  G4double tmax = std::min(maxEnergy, maxPairEnergy);
427  G4double cut = std::max(cutEnergy, minPairEnergy);
428  if (cut >= tmax) { return cross; }
429 
430  cross = ComputeMicroscopicCrossSection(kineticEnergy, Z, cut);
431  if(tmax < kineticEnergy) {
432  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
433  }
434  return cross;
435 }
436 
437 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
438 
439 void G4MuPairProductionModel::MakeSamplingTables()
440 {
441  G4double factore = G4Exp(G4Log(emax/emin)/G4double(nbine));
442 
443  for (G4int iz=0; iz<nzdat; ++iz) {
444 
445  G4double Z = zdat[iz];
446  G4Physics2DVector* pv = new G4Physics2DVector(nbiny+1,nbine+1);
447  G4double kinEnergy = emin;
448 
449  for (size_t it=0; it<=nbine; ++it) {
450 
451  pv->PutY(it, G4Log(kinEnergy/MeV));
452  G4double maxPairEnergy = MaxSecondaryEnergyForElement(kinEnergy, Z);
453  /*
454  G4cout << "it= " << it << " E= " << kinEnergy
455  << " " << particle->GetParticleName()
456  << " maxE= " << maxPairEnergy << " minE= " << minPairEnergy
457  << " ymin= " << ymin << G4endl;
458  */
459  G4double coef = G4Log(minPairEnergy/kinEnergy)/ymin;
460  G4double ymax = G4Log(maxPairEnergy/kinEnergy)/coef;
461  G4double fac = (ymax - ymin)/dy;
462  size_t imax = (size_t)fac;
463  fac -= (G4double)imax;
464 
465  G4double xSec = 0.0;
466  G4double x = ymin;
467  /*
468  G4cout << "Z= " << currentZ << " z13= " << z13
469  << " mE= " << maxPairEnergy << " ymin= " << ymin
470  << " dy= " << dy << " c= " << coef << G4endl;
471  */
472  // start from zero
473  pv->PutValue(0, it, 0.0);
474  if(0 == it) { pv->PutX(nbiny, 0.0); }
475 
476  for (size_t i=0; i<nbiny; ++i) {
477 
478  if(0 == it) { pv->PutX(i, x); }
479 
480  if(i < imax) {
481  G4double ep = kinEnergy*G4Exp(coef*(x + dy*0.5));
482 
483  // not multiplied by interval, because table
484  // will be used only for sampling
485  //G4cout << "i= " << i << " x= " << x << "E= " << kinEnergy
486  // << " Egamma= " << ep << G4endl;
487  xSec += ep*ComputeDMicroscopicCrossSection(kinEnergy, Z, ep);
488 
489  // last bin before the kinematic limit
490  } else if(i == imax) {
491  G4double ep = kinEnergy*G4Exp(coef*(x + fac*dy*0.5));
492  xSec += ep*fac*ComputeDMicroscopicCrossSection(kinEnergy, Z, ep);
493  }
494  pv->PutValue(i + 1, it, xSec);
495  x += dy;
496  }
497  kinEnergy *= factore;
498 
499  // to avoid precision lost
500  if(it+1 == nbine) { kinEnergy = emax; }
501  }
503  }
504 }
505 
506 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
507 
509  std::vector<G4DynamicParticle*>* vdp,
510  const G4MaterialCutsCouple* couple,
511  const G4DynamicParticle* aDynamicParticle,
512  G4double tmin,
513  G4double tmax)
514 {
515  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
516  //G4cout << "------- G4MuPairProductionModel::SampleSecondaries E(MeV)= "
517  // << kineticEnergy << " "
518  // << aDynamicParticle->GetDefinition()->GetParticleName() << G4endl;
519  G4double totalEnergy = kineticEnergy + particleMass;
520  G4double totalMomentum =
521  sqrt(kineticEnergy*(kineticEnergy + 2.0*particleMass));
522 
523  G4ThreeVector partDirection = aDynamicParticle->GetMomentumDirection();
524 
525  // select randomly one element constituing the material
526  const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
527 
528  // define interval of energy transfer
529  G4double maxPairEnergy = MaxSecondaryEnergyForElement(kineticEnergy,
530  anElement->GetZ());
531  G4double maxEnergy = std::min(tmax, maxPairEnergy);
532  G4double minEnergy = std::max(tmin, minPairEnergy);
533 
534  if(minEnergy >= maxEnergy) { return; }
535  //G4cout << "emin= " << minEnergy << " emax= " << maxEnergy
536  // << " minPair= " << minPairEnergy << " maxpair= " << maxPairEnergy
537  // << " ymin= " << ymin << " dy= " << dy << G4endl;
538 
539  G4double coeff = G4Log(minPairEnergy/kineticEnergy)/ymin;
540 
541  // compute limits
542  G4double yymin = G4Log(minEnergy/kineticEnergy)/coeff;
543  G4double yymax = G4Log(maxEnergy/kineticEnergy)/coeff;
544 
545  //G4cout << "yymin= " << yymin << " yymax= " << yymax << G4endl;
546 
547  // units should not be used, bacause table was built without
548  G4double logTkin = G4Log(kineticEnergy/MeV);
549 
550  // sample e-e+ energy, pair energy first
551 
552  // select sample table via Z
553  G4int iz1(0), iz2(0);
554  for(G4int iz=0; iz<nzdat; ++iz) {
555  if(currentZ == zdat[iz]) {
556  iz1 = iz2 = currentZ;
557  break;
558  } else if(currentZ < zdat[iz]) {
559  iz2 = zdat[iz];
560  if(iz > 0) { iz1 = zdat[iz-1]; }
561  else { iz1 = iz2; }
562  break;
563  }
564  }
565  if(0 == iz1) { iz1 = iz2 = zdat[nzdat-1]; }
566 
567  G4double PairEnergy = 0.0;
568  G4int count = 0;
569  //G4cout << "start loop Z1= " << iz1 << " Z2= " << iz2 << G4endl;
570  do {
571  ++count;
572  // sampling using only one random number
573  G4double rand = G4UniformRand();
574 
575  G4double x = FindScaledEnergy(iz1, rand, logTkin, yymin, yymax);
576  if(iz1 != iz2) {
577  G4double x2 = FindScaledEnergy(iz2, rand, logTkin, yymin, yymax);
578  G4double lz1= nist->GetLOGZ(iz1);
579  G4double lz2= nist->GetLOGZ(iz2);
580  //G4cout << count << ". x= " << x << " x2= " << x2
581  // << " Z1= " << iz1 << " Z2= " << iz2 << G4endl;
582  x += (x2 - x)*(lnZ - lz1)/(lz2 - lz1);
583  }
584  //G4cout << "x= " << x << " coeff= " << coeff << G4endl;
585  PairEnergy = kineticEnergy*G4Exp(x*coeff);
586 
587  } while((PairEnergy < minEnergy || PairEnergy > maxEnergy) && 10 > count);
588 
589  //G4cout << "## PairEnergy(GeV)= " << PairEnergy/GeV
590  // << " Etot(GeV)= " << totalEnergy/GeV << G4endl;
591 
592  // sample r=(E+-E-)/PairEnergy ( uniformly .....)
593  G4double rmax =
594  (1.-6.*particleMass*particleMass/(totalEnergy*(totalEnergy-PairEnergy)))
595  *sqrt(1.-minPairEnergy/PairEnergy);
596  G4double r = rmax * (-1.+2.*G4UniformRand()) ;
597 
598  // compute energies from PairEnergy,r
599  G4double ElectronEnergy = (1.-r)*PairEnergy*0.5;
600  G4double PositronEnergy = PairEnergy - ElectronEnergy;
601 
602  // The angle of the emitted virtual photon is sampled
603  // according to the muon bremsstrahlung model
604 
605  G4double gam = totalEnergy/particleMass;
606  G4double gmax = gam*std::min(1.0, totalEnergy/PairEnergy - 1.0);
607  G4double gmax2= gmax*gmax;
608  G4double x = G4UniformRand()*gmax2/(1.0 + gmax2);
609 
610  G4double theta = sqrt(x/(1.0 - x))/gam;
611  G4double sint = sin(theta);
612  G4double phi = twopi * G4UniformRand() ;
613  G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
614 
615  G4ThreeVector gDirection(dirx, diry, dirz);
616  gDirection.rotateUz(partDirection);
617 
618  // the angles of e- and e+ assumed to be the same as virtual gamma
619 
620  // create G4DynamicParticle object for the particle1
621  G4DynamicParticle* aParticle1 =
622  new G4DynamicParticle(theElectron, gDirection,
623  ElectronEnergy - electron_mass_c2);
624 
625  // create G4DynamicParticle object for the particle2
626  G4DynamicParticle* aParticle2 =
627  new G4DynamicParticle(thePositron, gDirection,
628  PositronEnergy - electron_mass_c2);
629 
630  // primary change
631  kineticEnergy -= (ElectronEnergy + PositronEnergy);
632  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
633 
634  partDirection *= totalMomentum;
635  partDirection -= (aParticle1->GetMomentum() + aParticle2->GetMomentum());
636  partDirection = partDirection.unit();
637  fParticleChange->SetProposedMomentumDirection(partDirection);
638 
639  // add secondary
640  vdp->push_back(aParticle1);
641  vdp->push_back(aParticle2);
642  //G4cout << "-- G4MuPairProductionModel::SampleSecondaries done" << G4endl;
643 }
644 
645 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
646 
647 void G4MuPairProductionModel::DataCorrupted(G4int Z, G4double logTkin)
648 {
650  ed << "G4ElementData is not properly initialized Z= " << Z
651  << " Ekin(MeV)= " << G4Exp(logTkin)
652  << " IsMasterThread= " << IsMaster()
653  << " Model " << GetName();
654  G4Exception("G4MuPairProductionModel::()","em0033",FatalException,
655  ed,"");
656 }
657 
658 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
659 
660 
G4double MaxSecondaryEnergyForElement(G4double kineticEnergy, G4double Z)
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static const G4double xgi[8]
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
std::vector< G4Element * > G4ElementVector
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
G4MuPairProductionModel(const G4ParticleDefinition *p=0, const G4String &nam="muPairProd")
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
const G4ParticleDefinition * particle
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
const char * p
Definition: xmltok.h:285
G4double GetZ() const
Definition: G4Element.hh:131
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
void SetParticle(const G4ParticleDefinition *)
G4ElementData * GetElementData()
Definition: G4VEmModel.hh:777
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double ComputMuPairLoss(G4double Z, G4double tkin, G4double cut, G4double tmax)
void PutValue(size_t idx, size_t idy, G4double value)
string material
Definition: eplot.py:19
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
#define G4UniformRand()
Definition: Randomize.hh:87
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
#define fm
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
G4double GetPDGMass() const
void PutX(size_t idx, G4double value)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double pairEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const G4double wgi[8]
G4double GetLOGZ(G4int Z)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static G4Electron * Electron()
Definition: G4Electron.cc:94
const G4String & GetName() const
Definition: G4VEmModel.hh:753
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
void PutY(size_t idy, G4double value)
G4ElementData * fElementData
Definition: G4VEmModel.hh:397
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
G4ThreeVector GetMomentum() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4fissionEvent * fe