Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MollerBhabhaModel.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: G4MollerBhabhaModel.cc 96934 2016-05-18 09:10:41Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4MollerBhabhaModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 03.01.2002
38 //
39 // Modifications:
40 //
41 // 13-11-02 Minor fix - use normalised direction (V.Ivanchenko)
42 // 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (V.Ivanchenko)
43 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
44 // 27-01-03 Make models region aware (V.Ivanchenko)
45 // 13-02-03 Add name (V.Ivanchenko)
46 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
47 // 25-07-05 Add protection in calculation of recoil direction for the case
48 // of complete energy transfer from e+ to e- (V.Ivanchenko)
49 // 06-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
50 // 15-05-06 Fix MinEnergyCut (V.Ivanchenko)
51 //
52 //
53 // Class Description:
54 //
55 // Implementation of energy loss and delta-electron production by e+/e-
56 //
57 // -------------------------------------------------------------------
58 //
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61 
62 #include "G4MollerBhabhaModel.hh"
63 #include "G4PhysicalConstants.hh"
64 #include "G4SystemOfUnits.hh"
65 #include "G4Electron.hh"
66 #include "G4Positron.hh"
67 #include "Randomize.hh"
69 #include "G4Log.hh"
70 #include "G4DeltaAngle.hh"
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73 
74 using namespace std;
75 
77  const G4String& nam)
78  : G4VEmModel(nam),
79  particle(nullptr),
80  isElectron(true),
81  twoln10(2.0*G4Log(10.0)),
82  lowLimit(0.02*keV),
83  isInitialised(false)
84 {
86  if(p) { SetParticle(p); }
87  fParticleChange = nullptr;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 {}
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98  G4double kinEnergy)
99 {
100  G4double tmax = kinEnergy;
101  if(isElectron) { tmax *= 0.5; }
102  return tmax;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 
108  const G4DataVector&)
109 {
110  if(!particle) { SetParticle(p); }
111 
112  if(isInitialised) { return; }
113 
114  isInitialised = true;
118  }
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
123 G4double
126  G4double cutEnergy,
127  G4double maxEnergy)
128 {
129  if(!particle) { SetParticle(p); }
130 
131  G4double cross = 0.0;
132  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
133  tmax = std::min(maxEnergy, tmax);
134  //G4cout << "E= " << kineticEnergy << " cut= " << cutEnergy
135  // << " Emax= " << tmax << G4endl;
136  if(cutEnergy < tmax) {
137 
138  G4double xmin = cutEnergy/kineticEnergy;
139  G4double xmax = tmax/kineticEnergy;
140  G4double tau = kineticEnergy/electron_mass_c2;
141  G4double gam = tau + 1.0;
142  G4double gamma2= gam*gam;
143  G4double beta2 = tau*(tau + 2)/gamma2;
144 
145  //Moller (e-e-) scattering
146  if (isElectron) {
147 
148  G4double gg = (2.0*gam - 1.0)/gamma2;
149  cross = ((xmax - xmin)*(1.0 - gg + 1.0/(xmin*xmax)
150  + 1.0/((1.0-xmin)*(1.0 - xmax)))
151  - gg*G4Log( xmax*(1.0 - xmin)/(xmin*(1.0 - xmax)) ) ) / beta2;
152 
153  //Bhabha (e+e-) scattering
154  } else {
155 
156  G4double y = 1.0/(1.0 + gam);
157  G4double y2 = y*y;
158  G4double y12 = 1.0 - 2.0*y;
159  G4double b1 = 2.0 - y2;
160  G4double b2 = y12*(3.0 + y2);
161  G4double y122= y12*y12;
162  G4double b4 = y122*y12;
163  G4double b3 = b4 + y122;
164 
165  cross = (xmax - xmin)*(1.0/(beta2*xmin*xmax) + b2
166  - 0.5*b3*(xmin + xmax)
167  + b4*(xmin*xmin + xmin*xmax + xmax*xmax)/3.0)
168  - b1*G4Log(xmax/xmin);
169  }
170 
171  cross *= twopi_mc2_rcl2/kineticEnergy;
172  }
173  return cross;
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177 
179  const G4ParticleDefinition* p,
182  G4double cutEnergy,
183  G4double maxEnergy)
184 {
185  return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189 
191  const G4Material* material,
192  const G4ParticleDefinition* p,
193  G4double kinEnergy,
194  G4double cutEnergy,
195  G4double maxEnergy)
196 {
197  G4double eDensity = material->GetElectronDensity();
198  return eDensity*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
199 }
200 
201 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202 
204  const G4Material* material,
205  const G4ParticleDefinition* p,
207  G4double cut)
208 {
209  if(nullptr == particle) { SetParticle(p); }
210  // calculate the dE/dx due to the ionization by Seltzer-Berger formula
211  // checl low-energy limit
212  G4double electronDensity = material->GetElectronDensity();
213 
214  G4double Zeff = material->GetIonisation()->GetZeffective();
215  G4double th = 0.25*sqrt(Zeff)*keV;
216  G4double tkin = std::max(kineticEnergy, th);
217 
218  G4double tau = tkin/electron_mass_c2;
219  G4double gam = tau + 1.0;
220  G4double gamma2= gam*gam;
221  G4double bg2 = tau*(tau + 2);
222  G4double beta2 = bg2/gamma2;
223 
224  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
225  eexc /= electron_mass_c2;
226  G4double eexc2 = eexc*eexc;
227 
229  G4double dedx;
230 
231  // electron
232  if (isElectron) {
233 
234  dedx = G4Log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
235  + G4Log((tau-d)*d) + tau/(tau-d)
236  + (0.5*d*d + (2.0*tau + 1.)*G4Log(1. - d/tau))/gamma2;
237 
238  //positron
239  } else {
240 
241  G4double d2 = d*d*0.5;
242  G4double d3 = d2*d/1.5;
243  G4double d4 = d3*d*0.75;
244  G4double y = 1.0/(1.0 + gam);
245  dedx = G4Log(2.0*(tau + 2.0)/eexc2) + G4Log(tau*d)
246  - beta2*(tau + 2.0*d - y*(3.0*d2
247  + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
248  }
249 
250  //density correction
251  G4double x = G4Log(bg2)/twoln10;
252  dedx -= material->GetIonisation()->DensityCorrection(x);
253 
254  // now you can compute the total ionization loss
255  dedx *= twopi_mc2_rcl2*electronDensity/beta2;
256  if (dedx < 0.0) { dedx = 0.0; }
257 
258  // lowenergy extrapolation
259 
260  if (kineticEnergy < th) {
261  x = kineticEnergy/th;
262  if(x > 0.25) { dedx /= sqrt(x); }
263  else { dedx *= 1.4*sqrt(x)/(0.1 + x); }
264  }
265  return dedx;
266 }
267 
268 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
269 
270 void
271 G4MollerBhabhaModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
272  const G4MaterialCutsCouple* couple,
273  const G4DynamicParticle* dp,
274  G4double cutEnergy,
275  G4double maxEnergy)
276 {
278  //G4cout << "G4MollerBhabhaModel::SampleSecondaries: E= " << kineticEnergy
279  // << " in " << couple->GetMaterial()->GetName() << G4endl;
280  G4double tmax;
281  G4double tmin = cutEnergy;
282  if(isElectron) {
283  tmax = 0.5*kineticEnergy;
284  } else {
285  tmax = kineticEnergy;
286  }
287  if(maxEnergy < tmax) { tmax = maxEnergy; }
288  if(tmin >= tmax) { return; }
289 
290  G4double energy = kineticEnergy + electron_mass_c2;
291  G4double xmin = tmin/kineticEnergy;
292  G4double xmax = tmax/kineticEnergy;
293  G4double gam = energy/electron_mass_c2;
294  G4double gamma2 = gam*gam;
295  G4double beta2 = 1.0 - 1.0/gamma2;
296  G4double x, z, grej;
297  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
298  G4double rndm[2];
299 
300  //Moller (e-e-) scattering
301  if (isElectron) {
302 
303  G4double gg = (2.0*gam - 1.0)/gamma2;
304  G4double y = 1.0 - xmax;
305  grej = 1.0 - gg*xmax + xmax*xmax*(1.0 - gg + (1.0 - gg*y)/(y*y));
306 
307  do {
308  rndmEngine->flatArray(2, rndm);
309  x = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
310  y = 1.0 - x;
311  z = 1.0 - gg*x + x*x*(1.0 - gg + (1.0 - gg*y)/(y*y));
312  /*
313  if(z > grej) {
314  G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
315  << "Majorant " << grej << " < "
316  << z << " for x= " << x
317  << " e-e- scattering"
318  << G4endl;
319  }
320  */
321  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
322  } while(grej * rndm[1] > z);
323 
324  //Bhabha (e+e-) scattering
325  } else {
326 
327  G4double y = 1.0/(1.0 + gam);
328  G4double y2 = y*y;
329  G4double y12 = 1.0 - 2.0*y;
330  G4double b1 = 2.0 - y2;
331  G4double b2 = y12*(3.0 + y2);
332  G4double y122= y12*y12;
333  G4double b4 = y122*y12;
334  G4double b3 = b4 + y122;
335 
336  y = xmax*xmax;
337  grej = 1.0 + (y*y*b4 - xmin*xmin*xmin*b3 + y*b2 - xmin*b1)*beta2;
338  do {
339  rndmEngine->flatArray(2, rndm);
340  x = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
341  y = x*x;
342  z = 1.0 + (y*y*b4 - x*y*b3 + y*b2 - x*b1)*beta2;
343  /*
344  if(z > grej) {
345  G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
346  << "Majorant " << grej << " < "
347  << z << " for x= " << x
348  << " e+e- scattering"
349  << G4endl;
350  }
351  */
352  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
353  } while(grej * rndm[1] > z);
354  }
355 
356  G4double deltaKinEnergy = x * kineticEnergy;
357 
358  G4ThreeVector deltaDirection;
359 
361  const G4Material* mat = couple->GetMaterial();
363 
364  deltaDirection =
365  GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
366 
367  } else {
368 
369  G4double deltaMomentum =
370  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
371  G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
372  (deltaMomentum * dp->GetTotalMomentum());
373  if(cost > 1.0) { cost = 1.0; }
374  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
375 
376  G4double phi = twopi * rndmEngine->flat() ;
377 
378  deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
379  deltaDirection.rotateUz(dp->GetMomentumDirection());
380  }
381 
382  // create G4DynamicParticle object for delta ray
383  G4DynamicParticle* delta =
384  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
385  vdp->push_back(delta);
386 
387  // primary change
388  kineticEnergy -= deltaKinEnergy;
389  G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
390  finalP = finalP.unit();
391 
394 }
395 
396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:313
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
static constexpr double keV
Definition: G4SIunits.hh:216
static constexpr double twopi_mc2_rcl2
const G4ThreeVector & GetMomentumDirection() const
Float_t y
Definition: compare.C:6
void SetParticle(const G4ParticleDefinition *p)
const char * p
Definition: xmltok.h:285
Double_t z
G4double DensityCorrection(G4double x)
virtual void flatArray(const int size, double *vect)=0
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Float_t y2[n_points_geant4]
Definition: compare.C:26
Float_t Z
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:116
double G4double
Definition: G4Types.hh:76
virtual double flat()=0
static const G4double d2
static constexpr double electron_mass_c2
static constexpr double twopi
Definition: G4SIunits.hh:76
double energy
Definition: plottest35.C:25
Float_t d
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:582
Float_t mat
G4double GetMeanExcitationEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
Hep3Vector unit() const
G4ThreeVector GetMomentum() const
int G4int
Definition: G4Types.hh:78
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:672
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:227
G4double GetZeffective() const
G4ParticleDefinition * theElectron
G4double GetKineticEnergy() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4MollerBhabhaModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MollerBhabha")
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:575
const G4ParticleDefinition * particle
const G4Material * GetMaterial() const
G4double GetTotalMomentum() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4bool isElectron(G4int ityp)
G4double GetElectronDensity() const
Definition: G4Material.hh:218
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4ParticleChangeForLoss * fParticleChange