Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros 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: G4MuBremsstrahlungModel.cc 97392 2016-06-02 10:10:32Z gcosmo $
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 // 31-05-13 Use element selectors instead of local data structure (V.Ivanchenko)
54 //
55 
56 //
57 // Class Description:
58 //
59 //
60 // -------------------------------------------------------------------
61 //
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 
66 #include "G4PhysicalConstants.hh"
67 #include "G4SystemOfUnits.hh"
68 #include "G4Gamma.hh"
69 #include "G4MuonMinus.hh"
70 #include "G4MuonPlus.hh"
71 #include "Randomize.hh"
72 #include "G4Material.hh"
73 #include "G4Element.hh"
74 #include "G4ElementVector.hh"
75 #include "G4ProductionCutsTable.hh"
77 #include "G4Log.hh"
78 #include "G4Exp.hh"
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82 
83 using namespace std;
84 
85 const G4double G4MuBremsstrahlungModel::xgi[] =
86  {0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
87 const G4double G4MuBremsstrahlungModel::wgi[] =
88  {0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
89 G4double G4MuBremsstrahlungModel::fDN[] = {0.0};
90 
92  const G4String& nam)
93  : G4VEmModel(nam),
94  particle(nullptr),
95  sqrte(sqrt(G4Exp(1.))),
96  bh(202.4),
97  bh1(446.),
98  btf(183.),
99  btf1(1429.),
100  fParticleChange(nullptr),
101  lowestKinEnergy(1.0*GeV),
102  minThreshold(0.9*keV)
103 {
104  theGamma = G4Gamma::Gamma();
106 
107  lowestKinEnergy = 1.*GeV;
108 
109  mass = rmass = cc = coeff = 1.0;
110 
111  if(0.0 == fDN[1]) {
112  for(G4int i=1; i<93; ++i) {
113  G4double dn = 1.54*nist->GetA27(i);
114  fDN[i] = dn;
115  if(1 < i) {
116  fDN[i] /= std::pow(dn, 1./G4double(i));
117  }
118  }
119  }
120 
121  if(p) { SetParticle(p); }
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125 
127 {}
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
132  const G4MaterialCutsCouple*)
133 {
134  return minThreshold;
135 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138 
140  const G4ParticleDefinition*,
141  G4double cut)
142 {
143  return std::max(lowestKinEnergy,cut);
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 
149  const G4DataVector& cuts)
150 {
151  if(p) { SetParticle(p); }
152 
153  // define pointer to G4ParticleChange
154  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
155 
156  if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
157  InitialiseElementSelectors(p, cuts);
158  }
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162 
164  G4VEmModel* masterModel)
165 {
166  if(p == particle && lowestKinEnergy < HighEnergyLimit()) {
168  }
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
172 
174  const G4Material* material,
175  const G4ParticleDefinition*,
176  G4double kineticEnergy,
177  G4double cutEnergy)
178 {
179  G4double dedx = 0.0;
180  if (kineticEnergy <= lowestKinEnergy) { return dedx; }
181 
182  G4double tmax = kineticEnergy;
183  G4double cut = std::min(cutEnergy,tmax);
184  if(cut < minThreshold) { cut = minThreshold; }
185 
186  const G4ElementVector* theElementVector = material->GetElementVector();
187  const G4double* theAtomicNumDensityVector =
188  material->GetAtomicNumDensityVector();
189 
190  // loop for elements in the material
191  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
192 
193  G4double loss =
194  ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
195 
196  dedx += loss*theAtomicNumDensityVector[i];
197  }
198  // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
199  if(dedx < 0.) dedx = 0.;
200  return dedx;
201 }
202 
203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204 
206  G4double tkin, G4double cut)
207 {
208  G4double totalEnergy = mass + tkin;
209  static const G4double ak1 = 0.05;
210  static const G4int k2=5;
211  G4double loss = 0.;
212 
213  G4double vcut = cut/totalEnergy;
214  G4double vmax = tkin/totalEnergy;
215 
216  G4double aaa = 0.;
217  G4double bbb = vcut;
218  if(vcut>vmax) { bbb = vmax; }
219  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2;
220  if(kkk < 1) { kkk = 1; }
221 
222  G4double hhh=(bbb-aaa)/G4double(kkk);
223 
224  G4double aa = aaa;
225  for(G4int l=0; l<kkk; l++)
226  {
227  for(G4int i=0; i<6; i++)
228  {
229  G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
230  loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
231  }
232  aa += hhh;
233  }
234 
235  loss *=hhh*totalEnergy ;
236 
237  return loss;
238 }
239 
240 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
241 
243  G4double tkin,
244  G4double Z,
245  G4double cut)
246 {
247  G4double totalEnergy = tkin + mass;
248  static const G4double ak1 = 2.3;
249  static const G4int k2 = 4;
250  G4double cross = 0.;
251 
252  if(cut >= tkin) return cross;
253 
254  G4double vcut = cut/totalEnergy;
255  G4double vmax = tkin/totalEnergy;
256 
257  G4double aaa = G4Log(vcut);
258  G4double bbb = G4Log(vmax);
259  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
260  if(kkk < 1) { kkk = 1; }
261 
262  G4double hhh = (bbb-aaa)/G4double(kkk);
263 
264  G4double aa = aaa;
265 
266  for(G4int l=0; l<kkk; l++)
267  {
268  for(G4int i=0; i<6; i++)
269  {
270  G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
271  cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
272  }
273  aa += hhh;
274  }
275 
276  cross *=hhh;
277 
278  //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
279 
280  return cross;
281 }
282 
283 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
284 
286  G4double tkin,
287  G4double Z,
288  G4double gammaEnergy)
289 // differential cross section
290 {
291  G4double dxsection = 0.;
292 
293  if(gammaEnergy > tkin) { return dxsection; }
294 
295  G4double E = tkin + mass ;
296  G4double v = gammaEnergy/E ;
297  G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
298  G4double rab0 = delta*sqrte ;
299 
300  G4int iz = G4lrint(Z);
301  if(iz < 1) { iz = 1; }
302  else if(iz > 92) { iz = 92; }
303 
304  G4double z13 = 1.0/nist->GetZ13(iz);
305  G4double dnstar = fDN[iz];
306 
307  G4double b,b1;
308 
309  if(1 == iz) {
310  b = bh;
311  b1 = bh1;
312  } else {
313  b = btf;
314  b1 = btf1;
315  }
316 
317  // nucleus contribution logarithm
318  G4double rab1=b*z13;
319  G4double fn=G4Log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
320  (mass+delta*(dnstar*sqrte-2.))) ;
321  if(fn <0.) { fn = 0.; }
322  // electron contribution logarithm
323  G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
324  G4double fe=0.;
325  if(gammaEnergy<epmax1)
326  {
327  G4double rab2=b1*z13*z13 ;
328  fe=G4Log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
329  (electron_mass_c2+rab0*rab2))) ;
330  if(fe<0.) { fe=0.; }
331  }
332 
333  dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
334 
335  return dxsection;
336 }
337 
338 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
339 
341  const G4ParticleDefinition*,
342  G4double kineticEnergy,
344  G4double cutEnergy,
345  G4double maxEnergy)
346 {
347  G4double cross = 0.0;
348  if (kineticEnergy <= lowestKinEnergy) return cross;
349  G4double tmax = std::min(maxEnergy, kineticEnergy);
350  G4double cut = std::min(cutEnergy, kineticEnergy);
351  if(cut < minThreshold) cut = minThreshold;
352  if (cut >= tmax) return cross;
353 
354  cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
355  if(tmax < kineticEnergy) {
356  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
357  }
358  return cross;
359 }
360 
361 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
362 
364  std::vector<G4DynamicParticle*>* vdp,
365  const G4MaterialCutsCouple* couple,
366  const G4DynamicParticle* dp,
367  G4double minEnergy,
368  G4double maxEnergy)
369 {
370  G4double kineticEnergy = dp->GetKineticEnergy();
371  // check against insufficient energy
372  G4double tmax = std::min(kineticEnergy, maxEnergy);
373  G4double tmin = std::min(kineticEnergy, minEnergy);
374  if(tmin < minThreshold) tmin = minThreshold;
375  if(tmin >= tmax) return;
376 
377  // ===== sampling of energy transfer ======
378 
379  G4ParticleMomentum partDirection = dp->GetMomentumDirection();
380 
381  // select randomly one element constituing the material
382  const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
383  G4double Z = anElement->GetZ();
384 
385  G4double totalEnergy = kineticEnergy + mass;
386  G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
387 
388  G4double func1 = tmin*
389  ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
390 
391  G4double lnepksi, epksi;
392  G4double func2;
393 
394  G4double xmin = G4Log(tmin/MeV);
395  G4double xmax = G4Log(kineticEnergy/tmin);
396 
397  do {
398  lnepksi = xmin + G4UniformRand()*xmax;
399  epksi = MeV*G4Exp(lnepksi);
400  func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
401 
402  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
403  } while(func2 < func1*G4UniformRand());
404 
405  G4double gEnergy = epksi;
406 
407  // ===== sample angle =====
408 
409  G4double gam = totalEnergy/mass;
410  G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
411  G4double rmax2= rmax*rmax;
412  G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
413 
414  G4double theta = sqrt(x/(1.0 - x))/gam;
415  G4double sint = sin(theta);
416  G4double phi = twopi * G4UniformRand() ;
417  G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
418 
419  G4ThreeVector gDirection(dirx, diry, dirz);
420  gDirection.rotateUz(partDirection);
421 
422  partDirection *= totalMomentum;
423  partDirection -= gEnergy*gDirection;
424  partDirection = partDirection.unit();
425 
426  // primary change
427  kineticEnergy -= gEnergy;
428  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
429  fParticleChange->SetProposedMomentumDirection(partDirection);
430 
431  // save secondary
432  G4DynamicParticle* aGamma =
433  new G4DynamicParticle(theGamma,gDirection,gEnergy);
434  vdp->push_back(aGamma);
435 }
436 
437 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
const char * p
Definition: xmltok.h:285
G4double GetA27(G4int Z) const
G4double GetZ() const
Definition: G4Element.hh:131
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
void SetParticle(const G4ParticleDefinition *)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
G4double GetZ13(G4double Z) const
#define G4UniformRand()
Definition: Randomize.hh:97
static const G4double sqrte
const G4ParticleDefinition * particle
const G4ThreeVector & GetMomentumDirection() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
static const G4double ak1
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
static constexpr double keV
Definition: G4SIunits.hh:216
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
G4fissionEvent * fe