Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MuBremsstrahlungModel.cc
이 파일의 문서화 페이지로 가기
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 108424 2018-02-13 11:19:25Z 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 
86  {0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
88  {0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
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 {
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  const G4MaterialCutsCouple*)
128 {
129  return minThreshold;
130 }
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133 
135  const G4ParticleDefinition*,
136  G4double cut)
137 {
138  return std::max(lowestKinEnergy,cut);
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
144  const G4DataVector& cuts)
145 {
146  if(p) { SetParticle(p); }
147 
148  // define pointer to G4ParticleChange
150 
151  if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
152  InitialiseElementSelectors(p, cuts);
153  }
154 }
155 
156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
157 
159  G4VEmModel* masterModel)
160 {
161  if(p == particle && lowestKinEnergy < HighEnergyLimit()) {
163  }
164 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167 
169  const G4Material* material,
170  const G4ParticleDefinition*,
172  G4double cutEnergy)
173 {
174  G4double dedx = 0.0;
175  if (kineticEnergy <= lowestKinEnergy) { return dedx; }
176 
177  G4double tmax = kineticEnergy;
178  G4double cut = std::min(cutEnergy,tmax);
179  if(cut < minThreshold) { cut = minThreshold; }
180 
181  const G4ElementVector* theElementVector = material->GetElementVector();
182  const G4double* theAtomicNumDensityVector =
183  material->GetAtomicNumDensityVector();
184 
185  // loop for elements in the material
186  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
187 
188  G4double loss =
189  ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
190 
191  dedx += loss*theAtomicNumDensityVector[i];
192  }
193  // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
194  if(dedx < 0.) dedx = 0.;
195  return dedx;
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199 
201  G4double tkin, G4double cut)
202 {
203  G4double totalEnergy = mass + tkin;
204  static const G4double ak1 = 0.05;
205  static const G4int k2=5;
206  G4double loss = 0.;
207 
208  G4double vcut = cut/totalEnergy;
209  G4double vmax = tkin/totalEnergy;
210 
211  G4double aaa = 0.;
212  G4double bbb = vcut;
213  if(vcut>vmax) { bbb = vmax; }
214  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2;
215  if(kkk < 1) { kkk = 1; }
216 
217  G4double hhh=(bbb-aaa)/G4double(kkk);
218 
219  G4double aa = aaa;
220  for(G4int l=0; l<kkk; l++)
221  {
222  for(G4int i=0; i<6; i++)
223  {
224  G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
225  loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
226  }
227  aa += hhh;
228  }
229 
230  loss *=hhh*totalEnergy ;
231 
232  return loss;
233 }
234 
235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
236 
238  G4double tkin,
239  G4double Z,
240  G4double cut)
241 {
242  G4double totalEnergy = tkin + mass;
243  static const G4double ak1 = 2.3;
244  static const G4int k2 = 4;
245  G4double cross = 0.;
246 
247  if(cut >= tkin) return cross;
248 
249  G4double vcut = cut/totalEnergy;
250  G4double vmax = tkin/totalEnergy;
251 
252  G4double aaa = G4Log(vcut);
253  G4double bbb = G4Log(vmax);
254  G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
255  if(kkk < 1) { kkk = 1; }
256 
257  G4double hhh = (bbb-aaa)/G4double(kkk);
258 
259  G4double aa = aaa;
260 
261  for(G4int l=0; l<kkk; l++)
262  {
263  for(G4int i=0; i<6; i++)
264  {
265  G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
266  cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
267  }
268  aa += hhh;
269  }
270 
271  cross *=hhh;
272 
273  //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
274 
275  return cross;
276 }
277 
278 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
279 
281  G4double tkin,
282  G4double Z,
283  G4double gammaEnergy)
284 // differential cross section
285 {
286  G4double dxsection = 0.;
287 
288  if(gammaEnergy > tkin) { return dxsection; }
289 
290  G4double E = tkin + mass ;
291  G4double v = gammaEnergy/E ;
292  G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
293  G4double rab0 = delta*sqrte ;
294 
295  G4int iz = G4lrint(Z);
296  if(iz < 1) { iz = 1; }
297  else if(iz > 92) { iz = 92; }
298 
299  G4double z13 = 1.0/nist->GetZ13(iz);
300  G4double dnstar = fDN[iz];
301 
302  G4double b,b1;
303 
304  if(1 == iz) {
305  b = bh;
306  b1 = bh1;
307  } else {
308  b = btf;
309  b1 = btf1;
310  }
311 
312  // nucleus contribution logarithm
313  G4double rab1=b*z13;
314  G4double fn=G4Log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
315  (mass+delta*(dnstar*sqrte-2.))) ;
316  if(fn <0.) { fn = 0.; }
317  // electron contribution logarithm
318  G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
319  G4double fe=0.;
320  if(gammaEnergy<epmax1)
321  {
322  G4double rab2=b1*z13*z13 ;
323  fe=G4Log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
324  (electron_mass_c2+rab0*rab2))) ;
325  if(fe<0.) { fe=0.; }
326  }
327 
328  dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
329 
330  return dxsection;
331 }
332 
333 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
334 
336  const G4ParticleDefinition*,
339  G4double cutEnergy,
340  G4double maxEnergy)
341 {
342  G4double cross = 0.0;
343  if (kineticEnergy <= lowestKinEnergy) return cross;
344  G4double tmax = std::min(maxEnergy, kineticEnergy);
345  G4double cut = std::min(cutEnergy, kineticEnergy);
346  if(cut < minThreshold) cut = minThreshold;
347  if (cut >= tmax) return cross;
348 
349  cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
350  if(tmax < kineticEnergy) {
351  cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
352  }
353  return cross;
354 }
355 
356 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
357 
359  std::vector<G4DynamicParticle*>* vdp,
360  const G4MaterialCutsCouple* couple,
361  const G4DynamicParticle* dp,
362  G4double minEnergy,
363  G4double maxEnergy)
364 {
366  // check against insufficient energy
367  G4double tmax = std::min(kineticEnergy, maxEnergy);
368  G4double tmin = std::min(kineticEnergy, minEnergy);
369  if(tmin < minThreshold) tmin = minThreshold;
370  if(tmin >= tmax) return;
371 
372  // ===== sampling of energy transfer ======
373 
374  G4ParticleMomentum partDirection = dp->GetMomentumDirection();
375 
376  // select randomly one element constituing the material
377  const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
378  G4double Z = anElement->GetZ();
379 
380  G4double totalEnergy = kineticEnergy + mass;
381  G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
382 
383  G4double func1 = tmin*
384  ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
385 
386  G4double lnepksi, epksi;
387  G4double func2;
388 
389  G4double xmin = G4Log(tmin/MeV);
390  G4double xmax = G4Log(kineticEnergy/tmin);
391 
392  do {
393  lnepksi = xmin + G4UniformRand()*xmax;
394  epksi = MeV*G4Exp(lnepksi);
395  func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
396 
397  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
398  } while(func2 < func1*G4UniformRand());
399 
400  G4double gEnergy = epksi;
401 
402  // ===== sample angle =====
403 
404  G4double gam = totalEnergy/mass;
405  G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
406  G4double rmax2= rmax*rmax;
407  G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
408 
409  G4double theta = sqrt(x/(1.0 - x))/gam;
410  G4double sint = sin(theta);
411  G4double phi = twopi * G4UniformRand() ;
412  G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
413 
414  G4ThreeVector gDirection(dirx, diry, dirz);
415  gDirection.rotateUz(partDirection);
416 
417  partDirection *= totalMomentum;
418  partDirection -= gEnergy*gDirection;
419  partDirection = partDirection.unit();
420 
421  // primary change
422  kineticEnergy -= gEnergy;
425 
426  // save secondary
427  G4DynamicParticle* aGamma =
428  new G4DynamicParticle(theGamma,gDirection,gEnergy);
429  vdp->push_back(aGamma);
430 }
431 
432 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:792
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double MeV
Definition: G4SIunits.hh:214
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
static constexpr double keV
Definition: G4SIunits.hh:216
static const G4double xgi[6]
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
int func2(int i)
Definition: XFunc.cc:51
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:191
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:552
void SetProposedMomentumDirection(const G4ThreeVector &dir)
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:784
static const G4double ak1
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Float_t Z
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:116
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
double G4double
Definition: G4Types.hh:76
const G4ParticleDefinition * particle
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
static constexpr double electron_mass_c2
G4fissionEvent * fe
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
int func1(int i)
Definition: XFunc.cc:40
static const G4double sqrte
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
Hep3Vector unit() const
static const G4double wgi[6]
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override
int G4lrint(double ad)
Definition: templates.hh:151
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
G4double GetA27(G4int Z) const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double GetKineticEnergy() const
G4ParticleDefinition * theGamma
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
G4double GetZ13(G4double Z) const
G4double GetZ() const
Definition: G4Element.hh:131
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:217
TH1F * hhh
Definition: AddMC.C:5
void SetParticle(const G4ParticleDefinition *)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
static constexpr double GeV
Definition: G4SIunits.hh:217
G4ParticleChangeForLoss * fParticleChange
size_t GetNumberOfElements() const
Definition: G4Material.hh:187
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4NistManager * Instance()
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)