Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MuBetheBlochModel.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: G4MuBetheBlochModel.cc 108424 2018-02-13 11:19:25Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4MuBetheBlochModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 09.08.2002
38 //
39 // Modifications:
40 //
41 // 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko)
42 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
43 // 27-01-03 Make models region aware (V.Ivanchenko)
44 // 13-02-03 Add name (V.Ivanchenko)
45 // 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I)
46 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
47 // 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko)
48 // 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49 //
50 
51 //
52 // -------------------------------------------------------------------
53 //
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 
58 #include "G4MuBetheBlochModel.hh"
59 #include "G4PhysicalConstants.hh"
60 #include "G4SystemOfUnits.hh"
61 #include "Randomize.hh"
62 #include "G4Electron.hh"
63 #include "G4LossTableManager.hh"
64 #include "G4EmCorrections.hh"
66 #include "G4Log.hh"
67 #include "G4Exp.hh"
68 
69 G4double G4MuBetheBlochModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917,
70  0.7628, 0.8983, 0.9801 };
71 
72 G4double G4MuBetheBlochModel::wgi[]={ 0.0506, 0.1112, 0.1569, 0.1813, 0.1813,
73  0.1569, 0.1112, 0.0506 };
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76 
77 using namespace std;
78 
80  const G4String& nam)
81  : G4VEmModel(nam),
82  particle(nullptr),
83  limitKinEnergy(100.*keV),
84  logLimitKinEnergy(G4Log(limitKinEnergy)),
85  twoln10(2.0*G4Log(10.0)),
86  //bg2lim(0.0169),
87  //taulim(8.4146e-3),
88  alphaprime(fine_structure_const/twopi)
89 {
92  fParticleChange = nullptr;
93 
94  // initial initialisation of memeber should be overwritten
95  // by SetParticle
96  mass = massSquare = ratio = 1.0;
97 
98  if(p) { SetParticle(p); }
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102 
104  const G4MaterialCutsCouple* couple)
105 {
106  return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy();
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
110 
112  G4double kinEnergy)
113 {
114  G4double tau = kinEnergy/mass;
115  G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
116  (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
117  return tmax;
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
121 
123  const G4DataVector&)
124 {
125  if(p) { SetParticle(p); }
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
132  const G4ParticleDefinition* p,
134  G4double cutEnergy,
135  G4double maxKinEnergy)
136 {
137  G4double cross = 0.0;
138  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
139  G4double maxEnergy = min(tmax,maxKinEnergy);
140  if(cutEnergy < maxEnergy) {
141 
142  G4double totEnergy = kineticEnergy + mass;
143  G4double energy2 = totEnergy*totEnergy;
144  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
145 
146  cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*G4Log(maxEnergy/cutEnergy)/tmax
147  + 0.5*(maxEnergy - cutEnergy)/energy2;
148 
149  // radiative corrections of R. Kokoulin
150  if (maxEnergy > limitKinEnergy) {
151 
152  G4double logtmax = G4Log(maxEnergy);
153  G4double logtmin = G4Log(max(cutEnergy,limitKinEnergy));
154  G4double logstep = logtmax - logtmin;
155  G4double dcross = 0.0;
156 
157  for (G4int ll=0; ll<8; ll++)
158  {
159  G4double ep = G4Exp(logtmin + xgi[ll]*logstep);
160  G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
161  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
162  dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
163  }
164 
165  cross += dcross*logstep*alphaprime;
166  }
167 
168  cross *= twopi_mc2_rcl2/beta2;
169 
170  }
171 
172  // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
173  // << " cross= " << cross << G4endl;
174 
175  return cross;
176 }
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179 
181  const G4ParticleDefinition* p,
184  G4double cutEnergy,
185  G4double maxEnergy)
186 {
188  (p,kineticEnergy,cutEnergy,maxEnergy);
189  return cross;
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 
195  const G4Material* material,
196  const G4ParticleDefinition* p,
198  G4double cutEnergy,
199  G4double maxEnergy)
200 {
201  G4double eDensity = material->GetElectronDensity();
203  (p,kineticEnergy,cutEnergy,maxEnergy);
204  return cross;
205 }
206 
207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
208 
210  const G4ParticleDefinition* p,
212  G4double cut)
213 {
214  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
215  G4double tau = kineticEnergy/mass;
216  G4double cutEnergy = min(cut,tmax);
217  G4double gam = tau + 1.0;
218  G4double bg2 = tau * (tau+2.0);
219  G4double beta2 = bg2/(gam*gam);
220 
221  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
222  G4double eexc2 = eexc*eexc;
223 
224  G4double eDensity = material->GetElectronDensity();
225 
226  G4double dedx = G4Log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2)
227  -(1.0 + cutEnergy/tmax)*beta2;
228 
229  G4double totEnergy = kineticEnergy + mass;
230  G4double del = 0.5*cutEnergy/totEnergy;
231  dedx += del*del;
232 
233  // density correction
234  G4double x = G4Log(bg2)/twoln10;
235  dedx -= material->GetIonisation()->DensityCorrection(x);
236 
237  // shell correction
238  dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
239 
240  // now compute the total ionization loss
241 
242  if (dedx < 0.0) dedx = 0.0 ;
243 
244  // radiative corrections of R. Kokoulin
245  if (cutEnergy > limitKinEnergy) {
246 
247  G4double logtmax = G4Log(cutEnergy);
248  G4double logstep = logtmax - logLimitKinEnergy;
249  G4double dloss = 0.0;
250  G4double ftot2= 0.5/(totEnergy*totEnergy);
251 
252  for (G4int ll=0; ll<8; ll++)
253  {
254  G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
255  G4double a1 = G4Log(1.0 + 2.0*ep/electron_mass_c2);
256  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
257  dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
258  }
259  dedx += dloss*logstep*alphaprime;
260  }
261 
262  dedx *= twopi_mc2_rcl2*eDensity/beta2;
263 
264  //High order corrections
265  dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
266 
267  return dedx;
268 }
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
271 
272 void G4MuBetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
273  const G4MaterialCutsCouple*,
274  const G4DynamicParticle* dp,
275  G4double minKinEnergy,
276  G4double maxEnergy)
277 {
278  G4double tmax = MaxSecondaryKinEnergy(dp);
279  G4double maxKinEnergy = min(maxEnergy,tmax);
280  if(minKinEnergy >= maxKinEnergy) { return; }
281 
283  G4double totEnergy = kineticEnergy + mass;
284  G4double etot2 = totEnergy*totEnergy;
285  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
286 
287  G4double grej = 1.;
288  if(tmax > limitKinEnergy) {
289  G4double a0 = G4Log(2.*totEnergy/mass);
290  grej += alphaprime*a0*a0;
291  }
292 
293  G4double deltaKinEnergy, f;
294 
295  // sampling follows ...
296  do {
297  G4double q = G4UniformRand();
298  deltaKinEnergy = minKinEnergy*maxKinEnergy
299  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
300 
301 
302  f = 1.0 - beta2*deltaKinEnergy/tmax
303  + 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
304 
305  if(deltaKinEnergy > limitKinEnergy) {
306  G4double a1 = G4Log(1.0 + 2.0*deltaKinEnergy/electron_mass_c2);
307  G4double a3 = G4Log(4.0*totEnergy*(totEnergy - deltaKinEnergy)/massSquare);
308  f *= (1. + alphaprime*a1*(a3 - a1));
309  }
310 
311  if(f > grej) {
312  G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
313  << "Majorant " << grej << " < "
314  << f << " for edelta= " << deltaKinEnergy
315  << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
316  << G4endl;
317  }
318  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
319  } while( grej*G4UniformRand() > f );
320 
321  G4double deltaMomentum =
322  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
323  G4double totalMomentum = totEnergy*sqrt(beta2);
324  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
325  (deltaMomentum * totalMomentum);
326 
327  G4double sint = sqrt(1.0 - cost*cost);
328 
329  G4double phi = twopi * G4UniformRand() ;
330 
331  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
332  G4ThreeVector direction = dp->GetMomentumDirection();
333  deltaDirection.rotateUz(direction);
334 
335  // primary change
336  kineticEnergy -= deltaKinEnergy;
337  G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
338  direction = dir.unit();
341 
342  // create G4DynamicParticle object for delta ray
344  deltaDirection,deltaKinEnergy);
345  vdp->push_back(delta);
346 }
347 
348 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
G4ParticleDefinition * theElectron
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double keV
Definition: G4SIunits.hh:216
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:492
static constexpr double twopi_mc2_rcl2
#define G4endl
Definition: G4ios.hh:61
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4ThreeVector & GetMomentumDirection() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
const char * p
Definition: xmltok.h:285
G4double DensityCorrection(G4double x)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void SetProposedMomentumDirection(const G4ThreeVector &dir)
static G4double xgi[8]
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Float_t Z
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:116
double G4double
Definition: G4Types.hh:76
G4EmCorrections * EmCorrections()
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double GetMeanExcitationEnergy() const
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetParticle(const G4ParticleDefinition *p)
Hep3Vector unit() const
const G4double a0
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
int G4int
Definition: G4Types.hh:78
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
void SetProposedKineticEnergy(G4double proposedKinEnergy)
TDirectory * dir
G4EmCorrections * corr
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:227
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
static G4LossTableManager * Instance()
const G4Material * GetMaterial() const
static constexpr double fine_structure_const
static G4double wgi[8]
G4double GetElectronDensity() const
Definition: G4Material.hh:218
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4ParticleChangeForLoss * fParticleChange
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4MuBetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBetheBloch")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override