Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MuBremsstrahlungModel.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: G4MuBremsstrahlungModel
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 PostStepDoIt (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 name (V.Ivanchenko)
46 // 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
47 // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
48 // 03-08-05 Angular correlations according to PRM (V.Ivanchenko)
49 // 13-02-06 add ComputeCrossSectionPerAtom (mma)
50 // 21-03-06 Fix problem of initialisation in case when cuts are not defined (VI)
51 // 07-11-07 Improve sampling of final state (A.Bogdanov)
52 // 28-02-08 Use precomputed Z^1/3 and Log(A) (V.Ivanchenko)
53 //
54 
55 //
56 // Class Description:
57 //
58 //
59 // -------------------------------------------------------------------
60 //
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
65 #include "G4PhysicalConstants.hh"
66 #include "G4SystemOfUnits.hh"
67 #include "G4Gamma.hh"
68 #include "G4MuonMinus.hh"
69 #include "G4MuonPlus.hh"
70 #include "Randomize.hh"
71 #include "G4Material.hh"
72 #include "G4Element.hh"
73 #include "G4ElementVector.hh"
74 #include "G4ProductionCutsTable.hh"
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
80 using namespace std;
81 
83  const G4String& nam)
84  : G4VEmModel(nam),
85  particle(0),
86  sqrte(sqrt(exp(1.))),
87  bh(202.4),
88  bh1(446.),
89  btf(183.),
90  btf1(1429.),
91  fParticleChange(0),
92  lowestKinEnergy(1.0*GeV),
93  minThreshold(1.0*keV)
94 {
95  theGamma = G4Gamma::Gamma();
97 
98  mass = rmass = cc = coeff = 1.0;
99 
100  fDN[0] = 0.0;
101  for(G4int i=1; i<93; ++i) {
102  G4double dn = 1.54*nist->GetA27(i);
103  fDN[i] = dn;
104  if(1 < i) {
105  fDN[i] /= std::pow(dn, 1./G4double(i));
106  }
107  }
108 
109  if(p) { SetParticle(p); }
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
115 {
116  size_t n = partialSumSigma.size();
117  if(n > 0) {
118  for(size_t i=0; i<n; i++) {
119  delete partialSumSigma[i];
120  }
121  }
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125 
127  const G4MaterialCutsCouple*)
128 {
129  return minThreshold;
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133 
135  const G4DataVector& cuts)
136 {
137  if(p) { SetParticle(p); }
138 
139  // partial cross section is computed for fixed energy
140  G4double fixedEnergy = 0.5*HighEnergyLimit();
141 
142  const G4ProductionCutsTable* theCoupleTable=
144  if(theCoupleTable) {
145  G4int numOfCouples = theCoupleTable->GetTableSize();
146 
147  G4int nn = partialSumSigma.size();
148  G4int nc = cuts.size();
149 
150  // do we need to perform initialisation?
151  if(nn == numOfCouples) { return; }
152 
153  // clear old data
154  if(nn > 0) {
155  for (G4int ii=0; ii<nn; ii++){
156  G4DataVector* a = partialSumSigma[ii];
157  if ( a ) { delete a; }
158  }
159  partialSumSigma.clear();
160  }
161  // fill new data
162  if (numOfCouples>0) {
163  for (G4int i=0; i<numOfCouples; i++) {
164  G4double cute = DBL_MAX;
165 
166  // protection for usage with extrapolator
167  if(i < nc) { cute = cuts[i]; }
168 
169  const G4MaterialCutsCouple* couple =
170  theCoupleTable->GetMaterialCutsCouple(i);
171  const G4Material* material = couple->GetMaterial();
172  G4DataVector* dv = ComputePartialSumSigma(material,fixedEnergy,cute);
173  partialSumSigma.push_back(dv);
174  }
175  }
176  }
177 
178  // define pointer to G4ParticleChange
179  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
180 }
181 
182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183 
185  const G4Material* material,
186  const G4ParticleDefinition*,
187  G4double kineticEnergy,
188  G4double cutEnergy)
189 {
190  G4double dedx = 0.0;
191  if (kineticEnergy <= lowestKinEnergy) return dedx;
192 
193  G4double tmax = kineticEnergy;
194  G4double cut = std::min(cutEnergy,tmax);
195  if(cut < minThreshold) cut = minThreshold;
196 
197  const G4ElementVector* theElementVector = material->GetElementVector();
198  const G4double* theAtomicNumDensityVector =
199  material->GetAtomicNumDensityVector();
200 
201  // loop for elements in the material
202  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
203 
204  G4double loss =
205  ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
206 
207  dedx += loss*theAtomicNumDensityVector[i];
208  }
209  // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
210  if(dedx < 0.) dedx = 0.;
211  return dedx;
212 }
213 
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215 
217  G4double tkin, G4double cut)
218 {
219  G4double totalEnergy = mass + tkin;
220  G4double ak1 = 0.05;
221  G4int k2=5;
222  G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
223  G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
224  G4double loss = 0.;
225 
226  G4double vcut = cut/totalEnergy;
227  G4double vmax = tkin/totalEnergy;
228 
229  G4double aaa = 0.;
230  G4double bbb = vcut;
231  if(vcut>vmax) bbb=vmax ;
232  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
233  if (kkk < 1) { kkk = 1; }
234  G4double hhh=(bbb-aaa)/float(kkk) ;
235 
236  G4double aa = aaa;
237  for(G4int l=0; l<kkk; l++)
238  {
239  for(G4int i=0; i<6; i++)
240  {
241  G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
242  loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
243  }
244  aa += hhh;
245  }
246 
247  loss *=hhh*totalEnergy ;
248 
249  return loss;
250 }
251 
252 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
253 
255  G4double tkin,
256  G4double Z,
257  G4double cut)
258 {
259  G4double totalEnergy = tkin + mass;
260  G4double ak1 = 2.3;
261  G4int k2 = 4;
262  G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
263  G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
264  G4double cross = 0.;
265 
266  if(cut >= tkin) return cross;
267 
268  G4double vcut = cut/totalEnergy;
269  G4double vmax = tkin/totalEnergy;
270 
271  G4double aaa = log(vcut);
272  G4double bbb = log(vmax);
273  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
274  G4double hhh = (bbb-aaa)/G4double(kkk);
275 
276  G4double aa = aaa;
277 
278  for(G4int l=0; l<kkk; l++)
279  {
280  for(G4int i=0; i<6; i++)
281  {
282  G4double ep = exp(aa + xgi[i]*hhh)*totalEnergy;
283  cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
284  }
285  aa += hhh;
286  }
287 
288  cross *=hhh;
289 
290  //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
291 
292  return cross;
293 }
294 
295 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
296 
298  G4double tkin,
299  G4double Z,
300  G4double gammaEnergy)
301 // differential cross section
302 {
303  G4double dxsection = 0.;
304 
305  if( gammaEnergy > tkin) return dxsection ;
306 
307  G4double E = tkin + mass ;
308  G4double v = gammaEnergy/E ;
309  G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
310  G4double rab0=delta*sqrte ;
311 
312  G4int iz = G4int(Z);
313  if(iz < 1) iz = 1;
314  else if(iz > 92) iz = 92;
315 
316  G4double z13 = 1.0/nist->GetZ13(iz);
317  G4double dnstar = fDN[iz];
318 
319  G4double b,b1;
320 
321  if(1 == iz) {
322  b = bh;
323  b1 = bh1;
324  } else {
325  b = btf;
326  b1 = btf1;
327  }
328 
329  // nucleus contribution logarithm
330  G4double rab1=b*z13;
331  G4double fn=log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
332  (mass+delta*(dnstar*sqrte-2.))) ;
333  if(fn <0.) fn = 0. ;
334  // electron contribution logarithm
335  G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
336  G4double fe=0.;
337  if(gammaEnergy<epmax1)
338  {
339  G4double rab2=b1*z13*z13 ;
340  fe=log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
341  (electron_mass_c2+rab0*rab2))) ;
342  if(fe<0.) fe=0. ;
343  }
344 
345  dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
346 
347  return dxsection;
348 }
349 
350 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
351 
353  const G4ParticleDefinition*,
354  G4double kineticEnergy,
356  G4double cutEnergy,
357  G4double maxEnergy)
358 {
359  G4double cross = 0.0;
360  if (kineticEnergy <= lowestKinEnergy) return cross;
361  G4double tmax = std::min(maxEnergy, kineticEnergy);
362  G4double cut = std::min(cutEnergy, kineticEnergy);
363  if(cut < minThreshold) cut = minThreshold;
364  if (cut >= tmax) return cross;
365 
366  cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
367  if(tmax < kineticEnergy) {
368  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
369  }
370  return cross;
371 }
372 
373 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374 
375 G4DataVector* G4MuBremsstrahlungModel::ComputePartialSumSigma(
376  const G4Material* material,
377  G4double kineticEnergy,
378  G4double cut)
379 
380 // Build the table of cross section per element.
381 // The table is built for material
382 // This table is used to select randomly an element in the material.
383 {
384  G4int nElements = material->GetNumberOfElements();
385  const G4ElementVector* theElementVector = material->GetElementVector();
386  const G4double* theAtomNumDensityVector =
387  material->GetAtomicNumDensityVector();
388 
389  G4DataVector* dv = new G4DataVector();
390 
391  G4double cross = 0.0;
392 
393  for (G4int i=0; i<nElements; i++ ) {
394  cross += theAtomNumDensityVector[i]
395  * ComputeMicroscopicCrossSection(kineticEnergy,
396  (*theElementVector)[i]->GetZ(), cut);
397  dv->push_back(cross);
398  }
399  return dv;
400 }
401 
402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
403 
405  std::vector<G4DynamicParticle*>* vdp,
406  const G4MaterialCutsCouple* couple,
407  const G4DynamicParticle* dp,
408  G4double minEnergy,
409  G4double maxEnergy)
410 {
411  G4double kineticEnergy = dp->GetKineticEnergy();
412  // check against insufficient energy
413  G4double tmax = std::min(kineticEnergy, maxEnergy);
414  G4double tmin = std::min(kineticEnergy, minEnergy);
415  if(tmin < minThreshold) tmin = minThreshold;
416  if(tmin >= tmax) return;
417 
418  // ===== sampling of energy transfer ======
419 
420  G4ParticleMomentum partDirection = dp->GetMomentumDirection();
421 
422  // select randomly one element constituing the material
423  const G4Element* anElement = SelectRandomAtom(couple);
424  G4double Z = anElement->GetZ();
425 
426  G4double totalEnergy = kineticEnergy + mass;
427  G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
428 
429  G4double func1 = tmin*
430  ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
431 
432  G4double lnepksi, epksi;
433  G4double func2;
434 
435  G4double xmin = log(tmin/MeV);
436  G4double xmax = log(kineticEnergy/tmin);
437 
438  do {
439  lnepksi = xmin + G4UniformRand()*xmax;
440  epksi = MeV*exp(lnepksi);
441  func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
442 
443  } while(func2 < func1*G4UniformRand());
444 
445  G4double gEnergy = epksi;
446 
447  // ===== sample angle =====
448 
449  G4double gam = totalEnergy/mass;
450  G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
451  G4double rmax2= rmax*rmax;
452  G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
453 
454  G4double theta = sqrt(x/(1.0 - x))/gam;
455  G4double sint = sin(theta);
456  G4double phi = twopi * G4UniformRand() ;
457  G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
458 
459  G4ThreeVector gDirection(dirx, diry, dirz);
460  gDirection.rotateUz(partDirection);
461 
462  partDirection *= totalMomentum;
463  partDirection -= gEnergy*gDirection;
464  partDirection = partDirection.unit();
465 
466  // primary change
467  kineticEnergy -= gEnergy;
468  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
469  fParticleChange->SetProposedMomentumDirection(partDirection);
470 
471  // save secondary
472  G4DynamicParticle* aGamma =
473  new G4DynamicParticle(theGamma,gDirection,gEnergy);
474  vdp->push_back(aGamma);
475 }
476 
477 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
478 
479 const G4Element* G4MuBremsstrahlungModel::SelectRandomAtom(
480  const G4MaterialCutsCouple* couple) const
481 {
482  // select randomly 1 element within the material
483 
484  const G4Material* material = couple->GetMaterial();
485  G4int nElements = material->GetNumberOfElements();
486  const G4ElementVector* theElementVector = material->GetElementVector();
487  if(1 == nElements) { return (*theElementVector)[0]; }
488  else if(1 > nElements) { return 0; }
489 
490  G4DataVector* dv = partialSumSigma[couple->GetIndex()];
491  G4double rval = G4UniformRand()*((*dv)[nElements-1]);
492  for (G4int i=0; i<nElements; i++) {
493  if (rval <= (*dv)[i]) { return (*theElementVector)[i]; }
494  }
495  return (*theElementVector)[nElements-1];
496 }
497 
498 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......