Geant4  10.01.p03
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 85023 2014-10-23 09:56:39Z 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 "G4LossTableManager.hh"
78 #include "G4Log.hh"
79 #include "G4Exp.hh"
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
84 using namespace std;
85 
87  {0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
89  {0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
91 
93  const G4String& nam)
94  : G4VEmModel(nam),
95  particle(0),
96  sqrte(sqrt(G4Exp(1.))),
97  bh(202.4),
98  bh1(446.),
99  btf(183.),
100  btf1(1429.),
101  fParticleChange(0),
102  lowestKinEnergy(1.0*GeV),
103  minThreshold(0.9*keV)
104 {
107 
108  lowestKinEnergy = 1.*GeV;
109 
110  mass = rmass = cc = coeff = 1.0;
111 
112  if(0.0 == fDN[1]) {
113  for(G4int i=1; i<93; ++i) {
114  G4double dn = 1.54*nist->GetA27(i);
115  fDN[i] = dn;
116  if(1 < i) {
117  fDN[i] /= std::pow(dn, 1./G4double(i));
118  }
119  }
120  }
121 
122  if(p) { SetParticle(p); }
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
128 {}
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131 
133  const G4MaterialCutsCouple*)
134 {
135  return minThreshold;
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139 
141  const G4ParticleDefinition*,
142  G4double cut)
143 {
144  return std::max(lowestKinEnergy,cut);
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
150  const G4DataVector& cuts)
151 {
152  if(p) { SetParticle(p); }
153 
154  // define pointer to G4ParticleChange
156 
157  if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
158  InitialiseElementSelectors(p, cuts);
159  }
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
165  G4VEmModel* masterModel)
166 {
167  if(p == particle && lowestKinEnergy < HighEnergyLimit()) {
169  }
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173 
175  const G4Material* material,
176  const G4ParticleDefinition*,
177  G4double kineticEnergy,
178  G4double cutEnergy)
179 {
180  G4double dedx = 0.0;
181  if (kineticEnergy <= lowestKinEnergy) { return dedx; }
182 
183  G4double tmax = kineticEnergy;
184  G4double cut = std::min(cutEnergy,tmax);
185  if(cut < minThreshold) { cut = minThreshold; }
186 
187  const G4ElementVector* theElementVector = material->GetElementVector();
188  const G4double* theAtomicNumDensityVector =
189  material->GetAtomicNumDensityVector();
190 
191  // loop for elements in the material
192  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
193 
194  G4double loss =
195  ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
196 
197  dedx += loss*theAtomicNumDensityVector[i];
198  }
199  // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
200  if(dedx < 0.) dedx = 0.;
201  return dedx;
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205 
207  G4double tkin, G4double cut)
208 {
209  G4double totalEnergy = mass + tkin;
210  static const G4double ak1 = 0.05;
211  static const G4int k2=5;
212  G4double loss = 0.;
213 
214  G4double vcut = cut/totalEnergy;
215  G4double vmax = tkin/totalEnergy;
216 
217  G4double aaa = 0.;
218  G4double bbb = vcut;
219  if(vcut>vmax) { bbb = vmax; }
220  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2;
221  if(kkk < 1) { kkk = 1; }
222 
223  G4double hhh=(bbb-aaa)/G4double(kkk);
224 
225  G4double aa = aaa;
226  for(G4int l=0; l<kkk; l++)
227  {
228  for(G4int i=0; i<6; i++)
229  {
230  G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
231  loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
232  }
233  aa += hhh;
234  }
235 
236  loss *=hhh*totalEnergy ;
237 
238  return loss;
239 }
240 
241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
242 
244  G4double tkin,
245  G4double Z,
246  G4double cut)
247 {
248  G4double totalEnergy = tkin + mass;
249  static const G4double ak1 = 2.3;
250  static const G4int k2 = 4;
251  G4double cross = 0.;
252 
253  if(cut >= tkin) return cross;
254 
255  G4double vcut = cut/totalEnergy;
256  G4double vmax = tkin/totalEnergy;
257 
258  G4double aaa = G4Log(vcut);
259  G4double bbb = G4Log(vmax);
260  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
261  if(kkk < 1) { kkk = 1; }
262 
263  G4double hhh = (bbb-aaa)/G4double(kkk);
264 
265  G4double aa = aaa;
266 
267  for(G4int l=0; l<kkk; l++)
268  {
269  for(G4int i=0; i<6; i++)
270  {
271  G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
272  cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
273  }
274  aa += hhh;
275  }
276 
277  cross *=hhh;
278 
279  //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
280 
281  return cross;
282 }
283 
284 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
285 
287  G4double tkin,
288  G4double Z,
289  G4double gammaEnergy)
290 // differential cross section
291 {
292  G4double dxsection = 0.;
293 
294  if(gammaEnergy > tkin) { return dxsection; }
295 
296  G4double E = tkin + mass ;
297  G4double v = gammaEnergy/E ;
298  G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
299  G4double rab0 = delta*sqrte ;
300 
301  G4int iz = G4lrint(Z);
302  if(iz < 1) { iz = 1; }
303  else if(iz > 92) { iz = 92; }
304 
305  G4double z13 = 1.0/nist->GetZ13(iz);
306  G4double dnstar = fDN[iz];
307 
308  G4double b,b1;
309 
310  if(1 == iz) {
311  b = bh;
312  b1 = bh1;
313  } else {
314  b = btf;
315  b1 = btf1;
316  }
317 
318  // nucleus contribution logarithm
319  G4double rab1=b*z13;
320  G4double fn=G4Log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
321  (mass+delta*(dnstar*sqrte-2.))) ;
322  if(fn <0.) { fn = 0.; }
323  // electron contribution logarithm
324  G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
325  G4double fe=0.;
326  if(gammaEnergy<epmax1)
327  {
328  G4double rab2=b1*z13*z13 ;
329  fe=G4Log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
330  (electron_mass_c2+rab0*rab2))) ;
331  if(fe<0.) { fe=0.; }
332  }
333 
334  dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
335 
336  return dxsection;
337 }
338 
339 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
340 
342  const G4ParticleDefinition*,
343  G4double kineticEnergy,
344  G4double Z, G4double,
345  G4double cutEnergy,
346  G4double maxEnergy)
347 {
348  G4double cross = 0.0;
349  if (kineticEnergy <= lowestKinEnergy) return cross;
350  G4double tmax = std::min(maxEnergy, kineticEnergy);
351  G4double cut = std::min(cutEnergy, kineticEnergy);
352  if(cut < minThreshold) cut = minThreshold;
353  if (cut >= tmax) return cross;
354 
355  cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
356  if(tmax < kineticEnergy) {
357  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
358  }
359  return cross;
360 }
361 
362 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
363 
365  std::vector<G4DynamicParticle*>* vdp,
366  const G4MaterialCutsCouple* couple,
367  const G4DynamicParticle* dp,
368  G4double minEnergy,
369  G4double maxEnergy)
370 {
371  G4double kineticEnergy = dp->GetKineticEnergy();
372  // check against insufficient energy
373  G4double tmax = std::min(kineticEnergy, maxEnergy);
374  G4double tmin = std::min(kineticEnergy, minEnergy);
375  if(tmin < minThreshold) tmin = minThreshold;
376  if(tmin >= tmax) return;
377 
378  // ===== sampling of energy transfer ======
379 
380  G4ParticleMomentum partDirection = dp->GetMomentumDirection();
381 
382  // select randomly one element constituing the material
383  const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
384  G4double Z = anElement->GetZ();
385 
386  G4double totalEnergy = kineticEnergy + mass;
387  G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
388 
389  G4double func1 = tmin*
390  ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
391 
392  G4double lnepksi, epksi;
393  G4double func2;
394 
395  G4double xmin = G4Log(tmin/MeV);
396  G4double xmax = G4Log(kineticEnergy/tmin);
397 
398  do {
399  lnepksi = xmin + G4UniformRand()*xmax;
400  epksi = MeV*G4Exp(lnepksi);
401  func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
402 
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;
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......
static const double MeV
Definition: G4SIunits.hh:193
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=0, const G4String &nam="MuBrem")
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
int func1(int i)
Definition: XFunc.cc:40
std::vector< G4Element * > G4ElementVector
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
G4ParticleDefinition * theGamma
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:616
int func2(int i)
Definition: XFunc.cc:51
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetZ13(G4double Z)
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
G4double GetA27(G4int Z)
void SetParticle(const G4ParticleDefinition *)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
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 const G4double xgi[6]
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
static const G4double wgi[6]
#define G4UniformRand()
Definition: Randomize.hh:93
static const G4double sqrte
const G4ParticleDefinition * particle
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static const double GeV
Definition: G4SIunits.hh:196
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:784
void SetProposedMomentumDirection(const G4ThreeVector &dir)
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
G4double G4Log(G4double x)
Definition: G4Log.hh:230
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:792
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ParticleChangeForLoss * fParticleChange
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const G4double b1
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
static const double keV
Definition: G4SIunits.hh:195
double G4double
Definition: G4Types.hh:76
G4ThreeVector G4ParticleMomentum
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:526
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4fissionEvent * fe