Geant4 11.1.1
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
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
77 const G4String& nam)
78 : G4VEmModel(nam),
79 limitKinEnergy(100.*CLHEP::keV),
80 logLimitKinEnergy(G4Log(limitKinEnergy)),
81 twoln10(2.0*G4Log(10.0)),
82 alphaprime(CLHEP::fine_structure_const/CLHEP::twopi)
83{
84 theElectron = G4Electron::Electron();
86 if(nullptr != p) { SetParticle(p); }
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
92 const G4MaterialCutsCouple* couple)
93{
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
98
100 G4double kinEnergy)
101{
102 G4double tau = kinEnergy/mass;
103 G4double tmax = 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
104 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
105 return tmax;
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109
110void G4MuBetheBlochModel::SetParticle(const G4ParticleDefinition* p)
111{
112 if(nullptr == particle) {
113 particle = p;
114 mass = particle->GetPDGMass();
115 massSquare = mass*mass;
116 ratio = CLHEP::electron_mass_c2/mass;
117 }
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
121
123 const G4DataVector&)
124{
125 SetParticle(p);
126 if(nullptr == fParticleChange) {
127 fParticleChange = GetParticleChangeForLoss();
128 }
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132
134 const G4ParticleDefinition* p,
135 G4double kineticEnergy,
136 G4double cutEnergy,
137 G4double maxKinEnergy)
138{
139 G4double cross = 0.0;
140 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
141 G4double maxEnergy = std::min(tmax, maxKinEnergy);
142 if(cutEnergy < maxEnergy) {
143
144 G4double totEnergy = kineticEnergy + mass;
145 G4double energy2 = totEnergy*totEnergy;
146 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
147
148 cross = 1.0/cutEnergy - 1.0/maxEnergy -
149 beta2*G4Log(maxEnergy/cutEnergy)/tmax +
150 0.5*(maxEnergy - cutEnergy)/energy2;
151
152 // radiative corrections of R. Kokoulin
153 if (maxEnergy > limitKinEnergy) {
154
155 G4double logtmax = G4Log(maxEnergy);
156 G4double logtmin = G4Log(std::max(cutEnergy,limitKinEnergy));
157 G4double logstep = logtmax - logtmin;
158 G4double dcross = 0.0;
159
160 for (G4int ll=0; ll<8; ++ll) {
161 G4double ep = G4Exp(logtmin + xgi[ll]*logstep);
162 G4double a1 = G4Log(1.0 + 2.0*ep/CLHEP::electron_mass_c2);
163 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
164 dcross += wgi[ll]*(1.0/ep - beta2/tmax + 0.5*ep/energy2)*a1*(a3 - a1);
165 }
166 cross += dcross*logstep*alphaprime;
167 }
168 cross *= CLHEP::twopi_mc2_rcl2/beta2;
169 }
170 // G4cout << "tmin= " << cutEnergy << " tmax= " << tmax
171 // << " cross= " << cross << G4endl;
172 return cross;
173}
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
176
178 const G4ParticleDefinition* p,
179 G4double kineticEnergy,
181 G4double cutEnergy,
182 G4double maxEnergy)
183{
185 (p,kineticEnergy,cutEnergy,maxEnergy);
186 return cross;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
192 const G4Material* material,
193 const G4ParticleDefinition* p,
194 G4double kineticEnergy,
195 G4double cutEnergy,
196 G4double maxEnergy)
197{
198 G4double eDensity = material->GetElectronDensity();
200 (p,kineticEnergy,cutEnergy,maxEnergy);
201 return cross;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
207 const G4ParticleDefinition* p,
208 G4double kineticEnergy,
209 G4double cut)
210{
211 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
212 G4double tau = kineticEnergy/mass;
213 G4double cutEnergy = std::min(cut, tmax);
214 G4double gam = tau + 1.0;
215 G4double bg2 = tau * (tau+2.0);
216 G4double beta2 = bg2/(gam*gam);
217
218 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
219 G4double eexc2 = eexc*eexc;
220
221 G4double eDensity = material->GetElectronDensity();
222
223 G4double dedx = G4Log(2.0*CLHEP::electron_mass_c2*bg2*cutEnergy/eexc2)
224 -(1.0 + cutEnergy/tmax)*beta2;
225
226 G4double totEnergy = kineticEnergy + mass;
227 G4double del = 0.5*cutEnergy/totEnergy;
228 dedx += del*del;
229
230 // density correction
231 G4double x = G4Log(bg2)/twoln10;
232 dedx -= material->GetIonisation()->DensityCorrection(x);
233
234 // shell correction
235 dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy);
236 dedx = std::max(dedx, 0.0);
237
238 // radiative corrections of R. Kokoulin
239 if (cutEnergy > limitKinEnergy) {
240
241 G4double logtmax = G4Log(cutEnergy);
242 G4double logstep = logtmax - logLimitKinEnergy;
243 G4double dloss = 0.0;
244 G4double ftot2= 0.5/(totEnergy*totEnergy);
245
246 for (G4int ll=0; ll<8; ++ll) {
247 G4double ep = G4Exp(logLimitKinEnergy + xgi[ll]*logstep);
248 G4double a1 = G4Log(1.0 + 2.0*ep/CLHEP::electron_mass_c2);
249 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - ep)/massSquare);
250 dloss += wgi[ll]*(1.0 - beta2*ep/tmax + ep*ep*ftot2)*a1*(a3 - a1);
251 }
252 dedx += dloss*logstep*alphaprime;
253 }
254
255 dedx *= CLHEP::twopi_mc2_rcl2*eDensity/beta2;
256
257 //High order corrections
258 dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy);
259 dedx = std::max(dedx, 0.);
260 return dedx;
261}
262
263//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
264
266 std::vector<G4DynamicParticle*>* vdp,
268 const G4DynamicParticle* dp,
269 G4double minKinEnergy,
270 G4double maxEnergy)
271{
272 G4double kineticEnergy = dp->GetKineticEnergy();
273 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kineticEnergy);
274 G4double maxKinEnergy = std::min(maxEnergy, tmax);
275 if(minKinEnergy >= maxKinEnergy) { return; }
276
277 G4double totEnergy = kineticEnergy + mass;
278 G4double etot2 = totEnergy*totEnergy;
279 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
280
281 G4double grej = 1.;
282 if(tmax > limitKinEnergy) {
283 G4double a0 = G4Log(2.*totEnergy/mass);
284 grej += alphaprime*a0*a0;
285 }
286
287 G4double tkin, f;
288
289 // sampling follows ...
290 do {
292 tkin = minKinEnergy*maxKinEnergy/(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
293 f = 1.0 - beta2*tkin/tmax + 0.5*tkin*tkin/etot2;
294
295 if(tkin > limitKinEnergy) {
296 G4double a1 = G4Log(1.0 + 2.0*tkin/CLHEP::electron_mass_c2);
297 G4double a3 = G4Log(4.0*totEnergy*(totEnergy - tkin)/massSquare);
298 f *= (1. + alphaprime*a1*(a3 - a1));
299 }
300
301 if(f > grej) {
302 G4cout << "G4MuBetheBlochModel::SampleSecondary Warning! "
303 << "Majorant " << grej << " < "
304 << f << " for edelta= " << tkin
305 << " tmin= " << minKinEnergy << " max= " << maxKinEnergy
306 << G4endl;
307 }
308 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
309 } while( grej*G4UniformRand() > f );
310
311 G4double deltaMomentum =
312 std::sqrt(tkin * (tkin + 2.0*CLHEP::electron_mass_c2));
313 G4double totalMomentum = totEnergy*std::sqrt(beta2);
314 G4double cost = tkin * (totEnergy + CLHEP::electron_mass_c2) /
315 (deltaMomentum * totalMomentum);
316
317 G4double sint = std::sqrt(1.0 - cost*cost);
318 G4double phi = CLHEP::twopi * G4UniformRand();
319 G4ThreeVector deltaDirection(sint*std::cos(phi), sint*std::sin(phi), cost);
320 G4ThreeVector direction = dp->GetMomentumDirection();
321 deltaDirection.rotateUz(direction);
322
323 // primary change
324 kineticEnergy -= tkin;
325 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
326 direction = dir.unit();
327 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
328 fParticleChange->SetProposedMomentumDirection(direction);
329
330 // create G4DynamicParticle object for delta ray
331 G4DynamicParticle* delta =
332 new G4DynamicParticle(theElectron, deltaDirection, tkin);
333 vdp->push_back(delta);
334}
335
336//....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
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
G4ParticleDefinition * GetDefinition() 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 DensityCorrection(G4double x) const
G4double GetMeanExcitationEnergy() const
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
G4double GetElectronDensity() const
Definition: G4Material.hh:212
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)
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
Definition: DoubConv.h:17