Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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#include "G4DeltaAngle.hh"
68
69G4double G4MuBetheBlochModel::xgi[]={ 0.0199, 0.1017, 0.2372, 0.4083, 0.5917,
70 0.7628, 0.8983, 0.9801 };
71
72G4double 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
78 const G4String& nam)
79 : G4VEmModel(nam),
80 limitRadCorrection(250.*CLHEP::MeV),
81 limitKinEnergy(100.*CLHEP::keV),
82 logLimitKinEnergy(G4Log(limitKinEnergy)),
83 twoln10(2.0*G4Log(10.0)),
84 alphaprime(CLHEP::fine_structure_const/CLHEP::twopi)
85{
86 theElectron = G4Electron::Electron();
88 if(nullptr != p) { SetParticle(p); }
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
100
102 G4double kinEnergy)
103{
104 G4double tau = kinEnergy/mass;
105 G4double tmax = 2.0*CLHEP::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
112void G4MuBetheBlochModel::SetParticle(const G4ParticleDefinition* p)
113{
114 if(nullptr == particle) {
115 particle = p;
116 mass = particle->GetPDGMass();
117 massSquare = mass*mass;
118 ratio = CLHEP::electron_mass_c2/mass;
119 }
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123
125 const G4DataVector&)
126{
127 SetParticle(p);
128 if(nullptr == fParticleChange) {
129 fParticleChange = GetParticleChangeForLoss();
130 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
132 }
133 }
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
139 const G4ParticleDefinition* p,
140 G4double kineticEnergy,
141 G4double cutEnergy,
142 G4double maxKinEnergy)
143{
144 G4double cross = 0.0;
145 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
146 G4double maxEnergy = std::min(tmax, maxKinEnergy);
147 if(cutEnergy < maxEnergy) {
148
149 G4double totEnergy = kineticEnergy + mass;
150 G4double energy2 = totEnergy*totEnergy;
151 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
152
153 cross = 1.0/cutEnergy - 1.0/maxEnergy -
154 beta2*G4Log(maxEnergy/cutEnergy)/tmax +
155 0.5*(maxEnergy - cutEnergy)/energy2;
156
157 // radiative corrections of R. Kokoulin
158 if (maxEnergy > limitKinEnergy && kineticEnergy > limitRadCorrection) {
159
160 G4double logtmax = G4Log(maxEnergy);
161 G4double logtmin = G4Log(std::max(cutEnergy,limitKinEnergy));
162 G4double logstep = logtmax - logtmin;
163 G4double dcross = 0.0;
164
165 for (G4int ll=0; ll<8; ++ll) {
166 G4double ep = G4Exp(logtmin + xgi[ll]*logstep);
167 G4double a1 = G4Log(1.0 + 2.0*ep/CLHEP::electron_mass_c2);
168 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
169 dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
170 }
171 cross += dcross*logstep*alphaprime;
172 }
173 cross *= CLHEP::twopi_mc2_rcl2/beta2;
174 }
175 // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
176 // << " cross= " << cross << G4endl;
177 return cross;
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181
183 const G4ParticleDefinition* p,
184 G4double kineticEnergy,
186 G4double cutEnergy,
187 G4double maxEnergy)
188{
190 (p,kineticEnergy,cutEnergy,maxEnergy);
191 return cross;
192}
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195
197 const G4Material* material,
198 const G4ParticleDefinition* p,
199 G4double kineticEnergy,
200 G4double cutEnergy,
201 G4double maxEnergy)
202{
203 G4double eDensity = material->GetElectronDensity();
205 (p,kineticEnergy,cutEnergy,maxEnergy);
206 return cross;
207}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
212 const G4ParticleDefinition* p,
213 G4double kineticEnergy,
214 G4double cut)
215{
216 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
217 G4double tau = kineticEnergy/mass;
218 G4double cutEnergy = std::min(cut, tmax);
219 G4double gam = tau + 1.0;
220 G4double bg2 = tau * (tau+2.0);
221 G4double beta2 = bg2/(gam*gam);
222
223 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
224 G4double eexc2 = eexc*eexc;
225
226 G4double eDensity = material->GetElectronDensity();
227
228 G4double dedx = G4Log(2.0*CLHEP::electron_mass_c2*bg2*cutEnergy/eexc2)
229 -(1.0 + cutEnergy/tmax)*beta2;
230
231 G4double totEnergy = kineticEnergy + mass;
232 G4double del = 0.5*cutEnergy/totEnergy;
233 dedx += del*del;
234
235 // density correction
236 G4double x = G4Log(bg2)/twoln10;
237 dedx -= material->GetIonisation()->DensityCorrection(x);
238
239 // shell and high order corrections
240 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
241
242 // radiative corrections of R. Kokoulin
243 if (cutEnergy > limitKinEnergy && kineticEnergy > limitRadCorrection) {
244
245 G4double logtmax = G4Log(cutEnergy);
246 G4double logstep = logtmax - logLimitKinEnergy;
247 G4double dloss = 0.0;
248 G4double ftot2= 0.5/(totEnergy*totEnergy);
249
250 for (G4int ll=0; ll<8; ++ll) {
251 G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
252 G4double a1 = G4Log(1.0 + 2.0*ep/CLHEP::electron_mass_c2);
253 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
254 dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
255 }
256 dedx += dloss*logstep*alphaprime;
257 }
258 dedx *= CLHEP::twopi_mc2_rcl2*eDensity/beta2;
259
260 //High order corrections
261 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
262 dedx = std::max(dedx, 0.);
263 return dedx;
264}
265
266//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
267
269 std::vector<G4DynamicParticle*>* vdp,
270 const G4MaterialCutsCouple* couple,
271 const G4DynamicParticle* dp,
272 G4double minKinEnergy,
273 G4double maxEnergy)
274{
275 G4double kineticEnergy = dp->GetKineticEnergy();
276 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kineticEnergy);
277 G4double maxKinEnergy = std::min(maxEnergy, tmax);
278 if(minKinEnergy >= maxKinEnergy) { return; }
279
280 G4double totEnergy = kineticEnergy + mass;
281 G4double etot2 = totEnergy*totEnergy;
282 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
283
284 G4double grej = 1.;
285 G4bool radC = (tmax > limitKinEnergy && kineticEnergy > limitRadCorrection);
286 if(radC) {
287 G4double a0 = G4Log(2.*totEnergy/mass);
288 grej += alphaprime*a0*a0;
289 }
290
291 G4double tkin, f;
292
293 // sampling follows ...
294 do {
296 tkin = minKinEnergy*maxKinEnergy/(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
297 f = 1.0 - beta2*tkin/tmax + 0.5*tkin*tkin/etot2;
298
299 if(radC && tkin > limitKinEnergy) {
300 G4double a1 = G4Log(1.0 + 2.0*tkin/CLHEP::electron_mass_c2);
301 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - tkin)/massSquare);
302 f *= (1. + alphaprime*a1*(a3 - a1));
303 }
304
305 if(f > grej) {
306 G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
307 << "Majorant " << grej << " < "
308 << f << " for edelta= " << tkin
309 << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
310 << G4endl;
311 }
312 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
313 } while( grej*G4UniformRand() > f );
314
315 G4ThreeVector deltaDirection;
316
318 const G4Material* mat = couple->GetMaterial();
319 deltaDirection = GetAngularDistribution()->SampleDirection(dp, tkin,
320 SelectRandomAtomNumber(mat), mat);
321 } else {
322
323 G4double deltaMom = std::sqrt(tkin * (tkin + 2.0*CLHEP::electron_mass_c2));
324 G4double totalMom = totEnergy*std::sqrt(beta2);
325 G4double cost = tkin * (totEnergy + CLHEP::electron_mass_c2) /
326 (deltaMom * totalMom);
327 cost = std::min(cost, 1.0);
328 const G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
329 const G4double phi = twopi*G4UniformRand();
330
331 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
332 deltaDirection.rotateUz(dp->GetMomentumDirection());
333 }
334 // create G4DynamicParticle object for delta ray
335 auto delta = new G4DynamicParticle(theElectron, deltaDirection, tkin);
336 vdp->push_back(delta);
337
338 // primary change
339 kineticEnergy -= tkin;
340 G4ThreeVector dir = dp->GetMomentum() - delta->GetMomentum();
341 dir = dir.unit();
342 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
343 fParticleChange->SetProposedMomentumDirection(dir);
344}
345
346//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
const G4double a0
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4Electron * Electron()
Definition G4Electron.cc:91
G4double HighOrderCorrections(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy, const G4double cutEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)
G4double DensityCorrection(G4double x) const
G4double GetMeanExcitationEnergy() const
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
G4double GetElectronDensity() const
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
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4MuBetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBetheBloch")
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)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4VEmAngularDistribution * GetAngularDistribution()
G4int SelectRandomAtomNumber(const G4Material *) const
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()