Geant4 10.7.0
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// Class Description:
57//
58//
59// -------------------------------------------------------------------
60//
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63
66#include "G4SystemOfUnits.hh"
67#include "G4Gamma.hh"
68#include "G4MuonMinus.hh"
69#include "G4MuonPlus.hh"
70#include "Randomize.hh"
71#include "G4Material.hh"
72#include "G4Element.hh"
73#include "G4ElementVector.hh"
75#include "G4ModifiedMephi.hh"
77#include "G4Log.hh"
78#include "G4Exp.hh"
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82
83using namespace std;
84
86 {0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
88 {0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
90
92 const G4String& nam)
93 : G4VEmModel(nam),
94 particle(nullptr),
95 sqrte(sqrt(G4Exp(1.))),
96 bh(202.4),
97 bh1(446.),
98 btf(183.),
99 btf1(1429.),
100 fParticleChange(nullptr),
101 lowestKinEnergy(1.0*GeV),
102 minThreshold(0.9*keV)
103{
106
107 lowestKinEnergy = 1.*GeV;
108
109 mass = rmass = cc = coeff = 1.0;
110
111 if(0.0 == fDN[1]) {
112 for(G4int i=1; i<93; ++i) {
113 G4double dn = 1.54*nist->GetA27(i);
114 fDN[i] = dn;
115 if(1 < i) {
116 fDN[i] /= std::pow(dn, 1./G4double(i));
117 }
118 }
119 }
121
122 if(p) { SetParticle(p); }
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126
129{
130 return minThreshold;
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
137 G4double cut)
138{
139 return std::max(lowestKinEnergy,cut);
140}
141
142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143
145 const G4DataVector& cuts)
146{
147 if(p) { SetParticle(p); }
148
149 // define pointer to G4ParticleChange
151
152 if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
154 }
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158
160 G4VEmModel* masterModel)
161{
164 }
165}
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
168
170 const G4Material* material,
172 G4double kineticEnergy,
173 G4double cutEnergy)
174{
175 G4double dedx = 0.0;
176 if (kineticEnergy <= lowestKinEnergy) { return dedx; }
177
178 G4double tmax = kineticEnergy;
179 G4double cut = std::min(cutEnergy,tmax);
180 if(cut < minThreshold) { cut = minThreshold; }
181
182 const G4ElementVector* theElementVector = material->GetElementVector();
183 const G4double* theAtomicNumDensityVector =
184 material->GetAtomicNumDensityVector();
185
186 // loop for elements in the material
187 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
188
189 G4double loss =
190 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
191
192 dedx += loss*theAtomicNumDensityVector[i];
193 }
194 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
195 if(dedx < 0.) dedx = 0.;
196 return dedx;
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
200
202 G4double tkin, G4double cut)
203{
204 G4double totalEnergy = mass + tkin;
205 static const G4double ak1 = 0.05;
206 static const G4int k2=5;
207 G4double loss = 0.;
208
209 G4double vcut = cut/totalEnergy;
210 G4double vmax = tkin/totalEnergy;
211
212 G4double aaa = 0.;
213 G4double bbb = vcut;
214 if(vcut>vmax) { bbb = vmax; }
215 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2;
216 if(kkk < 1) { kkk = 1; }
217
218 G4double hhh=(bbb-aaa)/G4double(kkk);
219
220 G4double aa = aaa;
221 for(G4int l=0; l<kkk; l++)
222 {
223 for(G4int i=0; i<6; i++)
224 {
225 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
226 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
227 }
228 aa += hhh;
229 }
230
231 loss *=hhh*totalEnergy ;
232
233 return loss;
234}
235
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
237
239 G4double tkin,
240 G4double Z,
241 G4double cut)
242{
243 G4double totalEnergy = tkin + mass;
244 static const G4double ak1 = 2.3;
245 static const G4int k2 = 4;
246 G4double cross = 0.;
247
248 if(cut >= tkin) return cross;
249
250 G4double vcut = cut/totalEnergy;
251 G4double vmax = tkin/totalEnergy;
252
253 G4double aaa = G4Log(vcut);
254 G4double bbb = G4Log(vmax);
255 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
256 if(kkk < 1) { kkk = 1; }
257
258 G4double hhh = (bbb-aaa)/G4double(kkk);
259
260 G4double aa = aaa;
261
262 for(G4int l=0; l<kkk; l++)
263 {
264 for(G4int i=0; i<6; i++)
265 {
266 G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
267 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
268 }
269 aa += hhh;
270 }
271
272 cross *=hhh;
273
274 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
275
276 return cross;
277}
278
279//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280
282 G4double tkin,
283 G4double Z,
284 G4double gammaEnergy)
285// differential cross section
286{
287 G4double dxsection = 0.;
288
289 if(gammaEnergy > tkin) { return dxsection; }
290
291 G4double E = tkin + mass ;
292 G4double v = gammaEnergy/E ;
293 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
294 G4double rab0 = delta*sqrte ;
295
296 G4int iz = G4lrint(Z);
297 if(iz < 1) { iz = 1; }
298 else if(iz > 92) { iz = 92; }
299
300 G4double z13 = 1.0/nist->GetZ13(iz);
301 G4double dnstar = fDN[iz];
302
303 G4double b,b1;
304
305 if(1 == iz) {
306 b = bh;
307 b1 = bh1;
308 } else {
309 b = btf;
310 b1 = btf1;
311 }
312
313 // nucleus contribution logarithm
314 G4double rab1=b*z13;
315 G4double fn=G4Log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
316 (mass+delta*(dnstar*sqrte-2.))) ;
317 if(fn <0.) { fn = 0.; }
318 // electron contribution logarithm
319 G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
320 G4double fe=0.;
321 if(gammaEnergy<epmax1)
322 {
323 G4double rab2=b1*z13*z13 ;
324 fe=G4Log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
325 (electron_mass_c2+rab0*rab2))) ;
326 if(fe<0.) { fe=0.; }
327 }
328
329 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
330
331 return dxsection;
332}
333
334//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
335
338 G4double kineticEnergy,
340 G4double cutEnergy,
341 G4double maxEnergy)
342{
343 G4double cross = 0.0;
344 if (kineticEnergy <= lowestKinEnergy) return cross;
345 G4double tmax = std::min(maxEnergy, kineticEnergy);
346 G4double cut = std::min(cutEnergy, kineticEnergy);
347 if(cut < minThreshold) cut = minThreshold;
348 if (cut >= tmax) return cross;
349
350 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
351 if(tmax < kineticEnergy) {
352 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
353 }
354 return cross;
355}
356
357//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
358
360 std::vector<G4DynamicParticle*>* vdp,
361 const G4MaterialCutsCouple* couple,
362 const G4DynamicParticle* dp,
363 G4double minEnergy,
364 G4double maxEnergy)
365{
366 G4double kineticEnergy = dp->GetKineticEnergy();
367 // check against insufficient energy
368 G4double tmax = std::min(kineticEnergy, maxEnergy);
369 G4double tmin = std::min(kineticEnergy, minEnergy);
370 if(tmin < minThreshold) tmin = minThreshold;
371 if(tmin >= tmax) return;
372
373 // ===== sampling of energy transfer ======
374
375 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
376
377 // select randomly one element constituing the material
378 const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
379 G4double Z = anElement->GetZ();
380
381 G4double func1 = tmin*
382 ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
383
384 G4double lnepksi, epksi;
385 G4double func2;
386
387 G4double xmin = G4Log(tmin/CLHEP::MeV);
388 G4double xmax = G4Log(kineticEnergy/tmin);
389
390 do {
391 lnepksi = xmin + G4UniformRand()*xmax;
392 epksi = CLHEP::MeV*G4Exp(lnepksi);
393 func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
394
395 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
396 } while(func2 < func1*G4UniformRand());
397
398 G4double gEnergy = epksi;
399
400 // ===== sample angle =====
401
402 //
403 // angles of the emitted gamma using general interface
404
405 G4ThreeVector gamDir =
406 GetAngularDistribution()->SampleDirection(dp, gEnergy, Z,
407 couple->GetMaterial());
408 // create G4DynamicParticle object for the Gamma
409 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma, gamDir, gEnergy);
410 vdp->push_back(gamma);
411 // compute post-interaction kinematics of primary e-/e+ based on
412 // energy-momentum conservation
413 const G4double totMomentum =
414 std::sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
415 G4ThreeVector dir =
416 (totMomentum*dp->GetMomentumDirection()-gEnergy*gamDir).unit();
417 const G4double finalE = kineticEnergy - gEnergy;
418 // if secondary gamma energy is higher than threshold(very high by default)
419 // then stop tracking the primary particle and create new secondary e-/e+
420 // instead of the primary one
421 if (gEnergy > SecondaryThreshold()) {
424 G4DynamicParticle* newdp = new G4DynamicParticle(particle, dir, finalE);
425 vdp->push_back(newdp);
426 } else { // continue tracking the primary e-/e+ otherwise
429 }
430}
431
432//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< G4Element * > G4ElementVector
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4fissionEvent * fe
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:130
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
virtual 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)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4ParticleDefinition * theGamma
G4ParticleChangeForLoss * fParticleChange
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
virtual 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 * > *)
Definition: G4VEmModel.hh:842
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:680
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void ProposeTrackStatus(G4TrackStatus status)
int G4lrint(double ad)
Definition: templates.hh:134