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