Geant4 11.1.1
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 if(nullptr != p) { SetParticle(p); }
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98
101{
102 return minThreshold;
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106
109 G4double cut)
110{
111 return std::max(lowestKinEnergy, cut);
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115
117{
118 if(nullptr == particle) {
119 particle = p;
121 rmass = mass/CLHEP::electron_mass_c2 ;
122 cc = CLHEP::classic_electr_radius/rmass ;
123 coeff = 16.*CLHEP::fine_structure_const*cc*cc/3. ;
124 }
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
128
130 const G4DataVector& cuts)
131{
132 SetParticle(p);
133 if(nullptr == fParticleChange) {
135 }
136 if(IsMaster() && p == particle && lowestKinEnergy < HighEnergyLimit()) {
137 if(0.0 == fDN[1]) {
138 for(G4int i=1; i<93; ++i) {
139 G4double dn = 1.54*nist->GetA27(i);
140 fDN[i] = dn;
141 if(1 < i) {
142 fDN[i] /= std::pow(dn, 1./G4double(i));
143 }
144 }
145 }
147 }
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151
153 G4VEmModel* masterModel)
154{
157 }
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
163 const G4Material* material,
165 G4double kineticEnergy,
166 G4double cutEnergy)
167{
168 G4double dedx = 0.0;
169 if (kineticEnergy <= lowestKinEnergy) { return dedx; }
170
171 G4double cut = std::max(cutEnergy, minThreshold);
172 cut = std::min(cut, kineticEnergy);
173
174 const G4ElementVector* theElementVector = material->GetElementVector();
175 const G4double* theAtomicNumDensityVector =
176 material->GetAtomicNumDensityVector();
177
178 // loop for elements in the material
179 for (size_t i=0; i<material->GetNumberOfElements(); ++i) {
180 G4double loss =
181 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
182 dedx += loss*theAtomicNumDensityVector[i];
183 }
184 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
185 dedx = std::max(dedx, 0.);
186 return dedx;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
192 G4double tkin, G4double cut)
193{
194 G4double totalEnergy = mass + tkin;
195 static const G4double ak1 = 0.05;
196 static const G4int k2 = 5;
197 G4double loss = 0.;
198
199 G4double vcut = cut/totalEnergy;
200 G4int kkk = (G4int)(vcut/ak1) + k2;
201 if (kkk > 8) { kkk = 8; }
202 else if (kkk < 1) { kkk = 1; }
203 G4double hhh = vcut/(G4double)(kkk);
204
205 G4double aa = 0.;
206 for(G4int l=0; l<kkk; ++l) {
207 for(G4int i=0; i<6; ++i) {
208 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
209 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
210 }
211 aa += hhh;
212 }
213
214 loss *= hhh*totalEnergy;
215 return loss;
216}
217
218//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
219
221 G4double tkin,
222 G4double Z,
223 G4double cut)
224{
225 G4double totalEnergy = tkin + mass;
226 static const G4double ak1 = 2.3;
227 static const G4int k2 = 4;
228 G4double cross = 0.;
229
230 if(cut >= tkin) return cross;
231
232 G4double vcut = cut/totalEnergy;
233 G4double vmax = tkin/totalEnergy;
234
235 G4double aaa = G4Log(vcut);
236 G4double bbb = G4Log(vmax);
237 G4int kkk = (G4int)((bbb-aaa)/ak1) + k2 ;
238 if(kkk > 8) { kkk = 8; }
239 else if (kkk < 1) { kkk = 1; }
240 G4double hhh = (bbb-aaa)/(G4double)(kkk);
241 G4double aa = aaa;
242
243 for(G4int l=0; l<kkk; ++l) {
244 for(G4int i=0; i<6; ++i) {
245 G4double ep = G4Exp(aa + xgi[i]*hhh)*totalEnergy;
246 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
247 }
248 aa += hhh;
249 }
250
251 cross *= hhh;
252 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
253 return cross;
254}
255
256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
257
259 G4double tkin,
260 G4double Z,
261 G4double gammaEnergy)
262// differential cross section
263{
264 G4double dxsection = 0.;
265 if(gammaEnergy > tkin) { return dxsection; }
266
267 G4double E = tkin + mass ;
268 G4double v = gammaEnergy/E ;
269 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
270 G4double rab0 = delta*sqrte ;
271
272 G4int iz = G4lrint(Z);
273 if(iz < 1) { iz = 1; }
274 else if(iz > 92) { iz = 92; }
275
276 G4double z13 = 1.0/nist->GetZ13(iz);
277 G4double dnstar = fDN[iz];
278
279 G4double b,b1;
280 if(1 == iz) {
281 b = bh;
282 b1 = bh1;
283 } else {
284 b = btf;
285 b1 = btf1;
286 }
287
288 // nucleus contribution logarithm
289 G4double rab1 = b*z13;
290 G4double fn = G4Log(rab1/(dnstar*(CLHEP::electron_mass_c2+rab0*rab1))*
291 (mass + delta*(dnstar*sqrte-2.)));
292 fn = std::max(fn, 0.);
293 // electron contribution logarithm
294 G4double epmax1 = E/(1.+0.5*mass*rmass/E);
295 G4double fe = 0.;
296 if(gammaEnergy < epmax1) {
297 G4double rab2 = b1*z13*z13;
298 fe = G4Log(rab2*mass/((1.+delta*rmass/(CLHEP::electron_mass_c2*sqrte))*
299 (CLHEP::electron_mass_c2+rab0*rab2)));
300 fe = std::max(fe, 0.);
301 }
302
303 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
304 dxsection = std::max(dxsection, 0.0);
305 return dxsection;
306}
307
308//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
309
312 G4double kineticEnergy,
314 G4double cutEnergy,
315 G4double maxEnergy)
316{
317 G4double cross = 0.0;
318 if (kineticEnergy <= lowestKinEnergy) return cross;
319 G4double tmax = std::min(maxEnergy, kineticEnergy);
320 G4double cut = std::min(cutEnergy, kineticEnergy);
321 if (cut < minThreshold) cut = minThreshold;
322 if (cut >= tmax) return cross;
323
324 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
325 if(tmax < kineticEnergy) {
326 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
327 }
328 return cross;
329}
330
331//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
332
334 std::vector<G4DynamicParticle*>* vdp,
335 const G4MaterialCutsCouple* couple,
336 const G4DynamicParticle* dp,
337 G4double minEnergy,
338 G4double maxEnergy)
339{
340 G4double kineticEnergy = dp->GetKineticEnergy();
341 // check against insufficient energy
342 G4double tmax = std::min(kineticEnergy, maxEnergy);
343 G4double tmin = std::min(kineticEnergy, minEnergy);
344 tmin = std::max(tmin, minThreshold);
345 if(tmin >= tmax) return;
346
347 // ===== sampling of energy transfer ======
348
349 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
350
351 // select randomly one element constituing the material
352 const G4Element* anElement = SelectRandomAtom(couple,particle,kineticEnergy);
353 G4double Z = anElement->GetZ();
354 G4double func1 = tmin*
355 ComputeDMicroscopicCrossSection(kineticEnergy, Z, tmin);
356
357 G4double gEnergy;
358 G4double func2;
359
360 G4double xmin = G4Log(tmin/minThreshold);
361 G4double xmax = G4Log(tmax/tmin);
362
363 do {
364 gEnergy = minThreshold*G4Exp(xmin + G4UniformRand()*xmax);
365 func2 = gEnergy*ComputeDMicroscopicCrossSection(kineticEnergy, Z, gEnergy);
366
367 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
368 } while(func2 < func1*G4UniformRand());
369
370 // angles of the emitted gamma using general interface
371 G4ThreeVector gamDir =
373 couple->GetMaterial());
374 // create G4DynamicParticle object for the Gamma
375 G4DynamicParticle* gamma = new G4DynamicParticle(theGamma, gamDir, gEnergy);
376 vdp->push_back(gamma);
377
378 // compute post-interaction kinematics of primary e-/e+ based on
379 // energy-momentum conservation
380 const G4double totMomentum =
381 std::sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
382 G4ThreeVector dir =
383 (totMomentum*dp->GetMomentumDirection() - gEnergy*gamDir).unit();
384 const G4double finalE = kineticEnergy - gEnergy;
385
386 // if secondary gamma energy is higher than threshold(very high by default)
387 // then stop tracking the primary particle and create new secondary e-/e+
388 // instead of the primary one
389 if (gEnergy > SecondaryThreshold()) {
392 G4DynamicParticle* newdp = new G4DynamicParticle(particle, dir, finalE);
393 vdp->push_back(newdp);
394 } else {
395 // continue tracking the primary e-/e+ otherwise
398 }
399}
400
401//....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
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:211
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 * > *)
Definition: G4VEmModel.hh:831
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:561
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
void ProposeTrackStatus(G4TrackStatus status)
Definition: DoubConv.h:17
int G4lrint(double ad)
Definition: templates.hh:134