Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuBremsstrahlungModel.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 file
30//
31//
32// File name: G4MuBremsstrahlungModel
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 24.06.2002
37//
38// Modifications:
39//
40// 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (V.Ivanchenko)
41// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
42// 24-01-03 Fix for compounds (V.Ivanchenko)
43// 27-01-03 Make models region aware (V.Ivanchenko)
44// 13-02-03 Add name (V.Ivanchenko)
45// 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
46// 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
47// 03-08-05 Angular correlations according to PRM (V.Ivanchenko)
48// 13-02-06 add ComputeCrossSectionPerAtom (mma)
49// 21-03-06 Fix problem of initialisation in case when cuts are not defined (VI)
50// 07-11-07 Improve sampling of final state (A.Bogdanov)
51// 28-02-08 Use precomputed Z^1/3 and Log(A) (V.Ivanchenko)
52// 31-05-13 Use element selectors instead of local data structure (V.Ivanchenko)
53//
54// -------------------------------------------------------------------
55//
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
61#include "G4SystemOfUnits.hh"
62#include "G4Gamma.hh"
63#include "G4MuonMinus.hh"
64#include "G4MuonPlus.hh"
65#include "Randomize.hh"
66#include "G4Material.hh"
67#include "G4Element.hh"
68#include "G4ElementVector.hh"
70#include "G4ModifiedMephi.hh"
72#include "G4Log.hh"
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
78 {0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
80 {0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
82
84 const G4String& nam)
85 : G4VEmModel(nam),
86 sqrte(std::sqrt(G4Exp(1.))),
87 lowestKinEnergy(0.1*CLHEP::GeV),
88 minThreshold(0.9*CLHEP::keV)
89{
92
94
95 if (nullptr != p) { SetParticle(p); }
96 if (0.0 == fDN[1]) {
97 for (G4int i=1; i<93; ++i) {
98 G4double dn = 1.54*nist->GetA27(i);
99 fDN[i] = dn;
100 if(1 < i) {
101 fDN[i] /= std::pow(dn, 1./G4double(i));
102 }
103 }
104 }
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
119 G4double cut)
120{
121 return std::max(lowestKinEnergy, cut);
122}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
127{
128 if(nullptr == particle) {
129 particle = p;
131 rmass = mass/CLHEP::electron_mass_c2 ;
132 cc = CLHEP::classic_electr_radius/rmass ;
133 coeff = 16.*CLHEP::fine_structure_const*cc*cc/3. ;
134 }
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138
140 const G4DataVector& cuts)
141{
142 SetParticle(p);
143 if(nullptr == fParticleChange) {
145 }
146 if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
148 }
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
152
154 G4VEmModel* masterModel)
155{
158 }
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162
164 const G4Material* material,
166 G4double kineticEnergy,
167 G4double cutEnergy)
168{
169 G4double dedx = 0.0;
170 if (kineticEnergy <= lowestKinEnergy) { return dedx; }
171
172 G4double cut = std::max(cutEnergy, minThreshold);
173 cut = std::min(cut, kineticEnergy);
174
175 const G4ElementVector* theElementVector = material->GetElementVector();
176 const G4double* theAtomicNumDensityVector =
177 material->GetAtomicNumDensityVector();
178
179 // loop for elements in the material
180 for (size_t i=0; i<material->GetNumberOfElements(); ++i) {
181 G4double loss =
182 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
183 dedx += loss*theAtomicNumDensityVector[i];
184 }
185 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
186 dedx = std::max(dedx, 0.);
187 return dedx;
188}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191
193 G4double tkin, G4double cut)
194{
195 G4double totalEnergy = mass + tkin;
196 static const G4double ak1 = 0.05;
197 static const G4int k2 = 5;
198 G4double loss = 0.;
199
200 G4double vcut = cut/totalEnergy;
201 G4int kkk = (G4int)(vcut/ak1) + k2;
202 if (kkk > 8) { kkk = 8; }
203 else if (kkk < 1) { kkk = 1; }
204 G4double hhh = vcut/(G4double)(kkk);
205
206 G4double aa = 0.;
207 for(G4int l=0; l<kkk; ++l) {
208 for(G4int i=0; i<6; ++i) {
209 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
210 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
211 }
212 aa += hhh;
213 }
214
215 loss *= hhh*totalEnergy;
216 return loss;
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220
222 G4double tkin,
223 G4double Z,
224 G4double cut)
225{
226 G4double totalEnergy = tkin + mass;
227 static const G4double ak1 = 2.3;
228 static const G4int k2 = 4;
229 G4double cross = 0.;
230
231 if(cut >= tkin) return cross;
232
233 G4double vcut = cut/totalEnergy;
234 G4double vmax = tkin/totalEnergy;
235
236 G4double aaa = G4Log(vcut);
237 G4double bbb = G4Log(vmax);
238 G4int kkk = (G4int)((bbb-aaa)/ak1) + k2 ;
239 if(kkk > 8) { kkk = 8; }
240 else if (kkk < 1) { kkk = 1; }
241 G4double hhh = (bbb-aaa)/(G4double)(kkk);
242 G4double aa = aaa;
243
244 for(G4int l=0; l<kkk; ++l) {
245 for(G4int i=0; i<6; ++i) {
246 G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
247 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
248 }
249 aa += hhh;
250 }
251
252 cross *= hhh;
253 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
254 return cross;
255}
256
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
258
260 G4double tkin,
261 G4double Z,
262 G4double gammaEnergy)
263// differential cross section
264{
265 G4double dxsection = 0.;
266 if(gammaEnergy > tkin) { return dxsection; }
267
268 G4double E = tkin + mass ;
269 G4double v = gammaEnergy/E ;
270 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
271 G4double rab0 = delta*sqrte ;
272
273 G4int iz = G4lrint(Z);
274 if(iz < 1) { iz = 1; }
275 else if(iz > 92) { iz = 92; }
276
277 G4double z13 = 1.0/nist->GetZ13(iz);
278 G4double dnstar = fDN[iz];
279
280 G4double b,b1;
281 if(1 == iz) {
282 b = bh;
283 b1 = bh1;
284 } else {
285 b = btf;
286 b1 = btf1;
287 }
288
289 // nucleus contribution logarithm
290 G4double rab1 = b*z13;
291 G4double fn = G4Log(rab1/(dnstar*(CLHEP::electron_mass_c2+rab0*rab1))*
292 (mass + delta*(dnstar*sqrte-2.)));
293 fn = std::max(fn, 0.);
294 // electron contribution logarithm
295 G4double epmax1 = E/(1.+0.5*mass*rmass/E);
296 G4double fe = 0.;
297 if(gammaEnergy < epmax1) {
298 G4double rab2 = b1*z13*z13;
299 fe = G4Log(rab2*mass/((1.+delta*rmass/(CLHEP::electron_mass_c2*sqrte))*
300 (CLHEP::electron_mass_c2+rab0*rab2)));
301 fe = std::max(fe, 0.);
302 }
303
304 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
305 dxsection = std::max(dxsection, 0.0);
306 return dxsection;
307}
308
309//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
310
313 G4double kineticEnergy,
315 G4double cutEnergy,
316 G4double maxEnergy)
317{
318 G4double cross = 0.0;
319 if (kineticEnergy <= lowestKinEnergy) return cross;
320 G4double tmax = std::min(maxEnergy, kineticEnergy);
321 G4double cut = std::min(cutEnergy, kineticEnergy);
322 if (cut < minThreshold) cut = minThreshold;
323 if (cut >= tmax) return cross;
324
325 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
326 if(tmax < kineticEnergy) {
327 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
328 }
329 return cross;
330}
331
332//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333
335 std::vector<G4DynamicParticle*>* vdp,
336 const G4MaterialCutsCouple* couple,
337 const G4DynamicParticle* dp,
338 G4double minEnergy,
339 G4double maxEnergy)
340{
341 G4double kineticEnergy = dp->GetKineticEnergy();
342 // check against insufficient energy
343 G4double tmax = std::min(kineticEnergy, maxEnergy);
344 G4double tmin = std::min(kineticEnergy, minEnergy);
345 tmin = std::max(tmin, minThreshold);
346 if(tmin >= tmax) return;
347
348 // ===== sampling of energy transfer ======
349
350 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
351
352 // select randomly one element constituing the material
353 const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
354 G4double Z = anElement->GetZ();
355 G4double func1 = tmin*
356 ComputeDMicroscopicCrossSection(kineticEnergy, Z, tmin);
357
358 G4double gEnergy;
359 G4double func2;
360
361 G4double xmin = G4Log(tmin/minThreshold);
362 G4double xmax = G4Log(tmax/tmin);
363
364 do {
365 gEnergy = minThreshold*G4Exp(xmin + G4UniformRand()*xmax);
366 func2 = gEnergy*ComputeDMicroscopicCrossSection(kineticEnergy, Z, gEnergy);
367
368 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
369 } while(func2 < func1*G4UniformRand());
370
371 // angles of the emitted gamma using general interface
372 G4ThreeVector gamDir =
373 GetAngularDistribution()->SampleDirection(dp, gEnergy, Z,
374 couple->GetMaterial());
375 // create G4DynamicParticle object for the Gamma
376 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma, gamDir, gEnergy);
377 vdp->push_back(gamma);
378
379 // compute post-interaction kinematics of primary e-/e+ based on
380 // energy-momentum conservation
381 const G4double totMomentum =
382 std::sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
383 G4ThreeVector dir =
384 (totMomentum*dp->GetMomentumDirection() - gEnergy*gamDir).unit();
385 const G4double finalE = kineticEnergy - gEnergy;
386
387 // if secondary gamma energy is higher than threshold(very high by default)
388 // then stop tracking the primary particle and create new secondary e-/e+
389 // instead of the primary one
390 if (gEnergy > SecondaryThreshold()) {
393 G4DynamicParticle* newdp = new G4DynamicParticle(particle, dir, finalE);
394 vdp->push_back(newdp);
395 } else {
396 // continue tracking the primary e-/e+ otherwise
399 }
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const G4Element * > G4ElementVector
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4fissionEvent * fe
G4double G4Log(G4double x)
Definition G4Log.hh:227
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition G4Element.hh:119
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
const G4double * GetAtomicNumDensityVector() const
std::size_t GetNumberOfElements() const
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
void SetParticle(const G4ParticleDefinition *)
static const G4double xgi[6]
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
static const G4double wgi[6]
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4ParticleDefinition * theGamma
G4ParticleChangeForLoss * fParticleChange
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override
const G4ParticleDefinition * particle
G4double GetA27(G4int Z) const
G4double GetZ13(G4double Z) const
static G4NistManager * Instance()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
G4VEmAngularDistribution * GetAngularDistribution()
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4bool IsMaster() const
G4double HighEnergyLimit() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetAngularDistribution(G4VEmAngularDistribution *)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
G4double SecondaryThreshold() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()
void ProposeTrackStatus(G4TrackStatus status)
int G4lrint(double ad)
Definition templates.hh:134