Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaConversionToMuons.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// ------------ G4GammaConversionToMuons physics process ------
28// by H.Burkhardt, S. Kelner and R. Kokoulin, April 2002
29//
30//
31// 07-08-02: missprint in OR condition in DoIt : f1<0 || f1>f1_max ..etc ...
32// 25-10-04: migrade to new interfaces of ParticleChange (vi)
33// ---------------------------------------------------------------------------
34
37#include "G4SystemOfUnits.hh"
38#include "G4UnitsTable.hh"
39#include "G4MuonPlus.hh"
40#include "G4MuonMinus.hh"
41#include "G4EmProcessSubType.hh"
42#include "G4EmParameters.hh"
43#include "G4LossTableManager.hh"
45#include "G4Gamma.hh"
46#include "G4Electron.hh"
47#include "G4Positron.hh"
48#include "G4NistManager.hh"
49#include "G4Log.hh"
50#include "G4Exp.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
54
55using namespace std;
56
57static const G4double sqrte=sqrt(exp(1.));
58static const G4double PowSat=-0.88;
59
61 G4ProcessType type)
62 : G4VDiscreteProcess (processName, type),
63 Mmuon(G4MuonPlus::MuonPlus()->GetPDGMass()),
64 Rc(CLHEP::elm_coupling/Mmuon),
65 LimitEnergy (5.*Mmuon),
66 LowestEnergyLimit (2.*Mmuon),
67 HighestEnergyLimit(1e12*CLHEP::GeV), // ok to 1e12GeV, then LPM suppression
68 theGamma(G4Gamma::Gamma()),
69 theMuonPlus(G4MuonPlus::MuonPlus()),
70 theMuonMinus(G4MuonMinus::MuonMinus())
71{
74 fManager->Register(this);
75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
78
80{
81 fManager->DeRegister(this);
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
85
87{
88 return (&part == theGamma);
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
94// Build cross section and mean free path tables
95{ //here no tables, just calling PrintInfoDefinition
97 if(Energy5DLimit > 0.0 && nullptr != f5Dmodel) {
98 f5Dmodel = new G4BetheHeitler5DModel();
99 f5Dmodel->SetLeptonPair(theMuonPlus, theMuonMinus);
100 const std::size_t numElems = G4ProductionCutsTable::GetProductionCutsTable()->GetTableSize();
101 const G4DataVector cuts(numElems);
102 f5Dmodel->Initialise(&p, cuts);
103 }
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
111
112// returns the photon mean free path in GEANT4 internal units
113// (MeanFreePath is a private member of the class)
114
115{
116 const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
117 G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
118 const G4Material* aMaterial = aTrack.GetMaterial();
119 return ComputeMeanFreePath(GammaEnergy, aMaterial);
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123
126 const G4Material* aMaterial)
127
128// computes and returns the photon mean free path in GEANT4 internal units
129{
130 if(GammaEnergy <= LowestEnergyLimit) { return DBL_MAX; }
131 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
132 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
133
134 G4double SIGMA = 0.0;
135 G4double fact = 1.0;
136 G4double e = GammaEnergy;
137 // low energy approximation as in Bethe-Heitler model
138 if(e < LimitEnergy) {
139 G4double y = (e - LowestEnergyLimit)/(LimitEnergy - LowestEnergyLimit);
140 fact = y*y;
141 e = LimitEnergy;
142 }
143
144 for ( std::size_t i=0 ; i < aMaterial->GetNumberOfElements(); ++i)
145 {
146 SIGMA += NbOfAtomsPerVolume[i] * fact *
147 ComputeCrossSectionPerAtom(e, (*theElementVector)[i]->GetZasInt());
148 }
149 return (SIGMA > 0.0) ? 1./SIGMA : DBL_MAX;
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153
155 const G4DynamicParticle* aDynamicGamma,
156 const G4Element* anElement)
157
158// gives the total cross section per atom in GEANT4 internal units
159{
160 return ComputeCrossSectionPerAtom(aDynamicGamma->GetKineticEnergy(),
161 anElement->GetZasInt());
162}
163
164//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
165
167 G4double Egam, G4int Z)
168
169// Calculates the microscopic cross section in GEANT4 internal units.
170// Total cross section parametrisation from H.Burkhardt
171// It gives a good description at any energy (from 0 to 10**21 eV)
172{
173 if(Egam <= LowestEnergyLimit) { return 0.0; }
174
176
177 G4double PowThres,Ecor,B,Dn,Zthird,Winfty,WMedAppr,
178 Wsatur,sigfac;
179
180 if(Z==1) // special case of Hydrogen
181 { B=202.4;
182 Dn=1.49;
183 }
184 else
185 { B=183.;
186 Dn=1.54*nist->GetA27(Z);
187 }
188 Zthird=1./nist->GetZ13(Z); // Z**(-1/3)
189 Winfty=B*Zthird*Mmuon/(Dn*electron_mass_c2);
190 WMedAppr=1./(4.*Dn*sqrte*Mmuon);
191 Wsatur=Winfty/WMedAppr;
192 sigfac=4.*fine_structure_const*Z*Z*Rc*Rc;
193 PowThres=1.479+0.00799*Dn;
194 Ecor=-18.+4347./(B*Zthird);
195
196 G4double CorFuc=1.+.04*G4Log(1.+Ecor/Egam);
197 //G4double Eg=pow(1.-4.*Mmuon/Egam,PowThres)*pow( pow(Wsatur,PowSat)+
198 // pow(Egam,PowSat),1./PowSat); // threshold and saturation
199 G4double Eg=G4Exp(G4Log(1.-4.*Mmuon/Egam)*PowThres)*
200 G4Exp(G4Log( G4Exp(G4Log(Wsatur)*PowSat)+G4Exp(G4Log(Egam)*PowSat))/PowSat);
201 G4double CrossSection=7./9.*sigfac*G4Log(1.+WMedAppr*CorFuc*Eg);
202 CrossSection *= CrossSecFactor; // increase the CrossSection by (by default 1)
203 return CrossSection;
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
207
209// Set the factor to artificially increase the cross section
210{
211 if(fac < 0.0) return;
212 CrossSecFactor=fac;
213 G4cout << "The cross section for GammaConversionToMuons is artificially "
214 << "increased by the CrossSecFactor=" << CrossSecFactor << G4endl;
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
218
220 const G4Track& aTrack,
221 const G4Step& aStep)
222//
223// generation of gamma->mu+mu-
224//
225{
227 const G4Material* aMaterial = aTrack.GetMaterial();
228
229 // current Gamma energy and direction, return if energy too low
230 const G4DynamicParticle *aDynamicGamma = aTrack.GetDynamicParticle();
231 G4double Egam = aDynamicGamma->GetKineticEnergy();
232 if (Egam <= LowestEnergyLimit) {
233 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
234 }
235 //
236 // Kill the incident photon
237 //
241
242 if (Egam <= Energy5DLimit) {
243 std::vector<G4DynamicParticle*> fvect;
244 f5Dmodel->SampleSecondaries(&fvect, aTrack.GetMaterialCutsCouple(),
245 aTrack.GetDynamicParticle(), 0.0, DBL_MAX);
247 for(auto dp : fvect) { aParticleChange.AddSecondary(dp); }
248 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
249 }
250
251 G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
252
253 // select randomly one element constituting the material
254 const G4Element* anElement = SelectRandomAtom(aDynamicGamma, aMaterial);
255 G4int Z = anElement->GetZasInt();
257
258 G4double B,Dn;
259 G4double A027 = nist->GetA27(Z);
260
261 if(Z==1) // special case of Hydrogen
262 { B=202.4;
263 Dn=1.49;
264 }
265 else
266 { B=183.;
267 Dn=1.54*A027;
268 }
269 G4double Zthird=1./nist->GetZ13(Z); // Z**(-1/3)
270 G4double Winfty=B*Zthird*Mmuon/(Dn*electron_mass_c2);
271
272 G4double C1Num=0.138*A027;
273 G4double C1Num2=C1Num*C1Num;
274 G4double C2Term2=electron_mass_c2/(183.*Zthird*Mmuon);
275
276 G4double GammaMuonInv=Mmuon/Egam;
277
278 // generate xPlus according to the differential cross section by rejection
279 G4double xmin=(Egam < LimitEnergy) ? GammaMuonInv : .5-sqrt(.25-GammaMuonInv);
280 G4double xmax=1.-xmin;
281
282 G4double Ds2=(Dn*sqrte-2.);
283 G4double sBZ=sqrte*B*Zthird/electron_mass_c2;
284 G4double LogWmaxInv=1./G4Log(Winfty*(1.+2.*Ds2*GammaMuonInv)
285 /(1.+2.*sBZ*Mmuon*GammaMuonInv));
286 G4double xPlus,xMinus,xPM,result,W;
287 G4int nn = 0;
288 const G4int nmax = 1000;
289 do {
290 xPlus=xmin+G4UniformRand()*(xmax-xmin);
291 xMinus=1.-xPlus;
292 xPM=xPlus*xMinus;
293 G4double del=Mmuon*Mmuon/(2.*Egam*xPM);
294 W=Winfty*(1.+Ds2*del/Mmuon)/(1.+sBZ*del);
295 G4double xxp=1.-4./3.*xPM; // the main xPlus dependence
296 result=(xxp > 0.) ? xxp*G4Log(W)*LogWmaxInv : 0.0;
297 if(result>1.) {
298 G4cout << "G4GammaConversionToMuons::PostStepDoIt WARNING:"
299 << " in dSigxPlusGen, result=" << result << " > 1" << G4endl;
300 }
301 ++nn;
302 if(nn >= nmax) { break; }
303 }
304 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
305 while (G4UniformRand() > result);
306
307 // now generate the angular variables via the auxilary variables t,psi,rho
308 G4double t;
309 G4double psi;
310 G4double rho;
311
312 G4double a3 = (GammaMuonInv/(2.*xPM));
313 G4double a33 = a3*a3;
314 G4double f1;
315 G4double b1 = 1./(4.*C1Num2);
316 G4double b3 = b1*b1*b1;
317 G4double a21 = a33 + b1;
318
319 G4double f1_max=-(1.-xPM)*(2.*b1+(a21+a33)*G4Log(a33/a21))/(2*b3);
320
321 G4double thetaPlus,thetaMinus,phiHalf; // final angular variables
322 nn = 0;
323 // t, psi, rho generation start (while angle < pi)
324 do {
325 //generate t by the rejection method
326 do {
327 ++nn;
328 t=G4UniformRand();
329 G4double a34=a33/(t*t);
330 G4double a22 = a34 + b1;
331 if(std::abs(b1)<0.0001*a34)
332 // special case of a34=a22 because of logarithm accuracy
333 {
334 f1=(1.-2.*xPM+4.*xPM*t*(1.-t))/(12.*a34*a34*a34*a34);
335 }
336 else
337 {
338 f1=-(1.-2.*xPM+4.*xPM*t*(1.-t))*(2.*b1+(a22+a34)*G4Log(a34/a22))/(2*b3);
339 }
340 if(f1<0.0 || f1> f1_max) // should never happend
341 {
342 G4cout << "G4GammaConversionToMuons::PostStepDoIt WARNING:"
343 << "outside allowed range f1=" << f1
344 << " is set to zero, a34 = "<< a34 << " a22 = "<<a22<<"."
345 << G4endl;
346 f1 = 0.0;
347 }
348 if(nn > nmax) { break; }
349 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
350 } while ( G4UniformRand()*f1_max > f1);
351 // generate psi by the rejection method
352 G4double f2_max=1.-2.*xPM*(1.-4.*t*(1.-t));
353 // long version
354 G4double f2;
355 do {
356 ++nn;
357 psi=twopi*G4UniformRand();
358 f2=1.-2.*xPM+4.*xPM*t*(1.-t)*(1.+cos(2.*psi));
359 if(f2<0 || f2> f2_max) // should never happend
360 {
361 G4cout << "G4GammaConversionToMuons::PostStepDoIt WARNING:"
362 << "outside allowed range f2=" << f2 << " is set to zero"
363 << G4endl;
364 f2 = 0.0;
365 }
366 if(nn >= nmax) { break; }
367 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
368 } while ( G4UniformRand()*f2_max > f2);
369
370 // generate rho by direct transformation
371 G4double C2Term1=GammaMuonInv/(2.*xPM*t);
372 G4double C22 = C2Term1*C2Term1+C2Term2*C2Term2;
373 G4double C2=4.*C22*C22/sqrt(xPM);
374 G4double rhomax=(1./t-1.)*1.9/A027;
375 G4double beta=G4Log( (C2+rhomax*rhomax*rhomax*rhomax)/C2 );
376 rho=G4Exp(G4Log(C2 *( G4Exp(beta*G4UniformRand())-1. ))*0.25);
377
378 //now get from t and psi the kinematical variables
379 G4double u=sqrt(1./t-1.);
380 G4double xiHalf=0.5*rho*cos(psi);
381 phiHalf=0.5*rho/u*sin(psi);
382
383 thetaPlus =GammaMuonInv*(u+xiHalf)/xPlus;
384 thetaMinus=GammaMuonInv*(u-xiHalf)/xMinus;
385
386 // protection against infinite loop
387 if(nn > nmax) {
388 if(std::abs(thetaPlus)>pi) { thetaPlus = 0.0; }
389 if(std::abs(thetaMinus)>pi) { thetaMinus = 0.0; }
390 }
391
392 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
393 } while ( std::abs(thetaPlus)>pi || std::abs(thetaMinus) >pi);
394
395 // now construct the vectors
396 // azimuthal symmetry, take phi0 at random between 0 and 2 pi
397 G4double phi0=twopi*G4UniformRand();
398 G4double EPlus=xPlus*Egam;
399 G4double EMinus=xMinus*Egam;
400
401 // mu+ mu- directions for gamma in z-direction
402 G4ThreeVector MuPlusDirection ( sin(thetaPlus) *cos(phi0+phiHalf),
403 sin(thetaPlus) *sin(phi0+phiHalf), cos(thetaPlus) );
404 G4ThreeVector MuMinusDirection (-sin(thetaMinus)*cos(phi0-phiHalf),
405 -sin(thetaMinus) *sin(phi0-phiHalf), cos(thetaMinus) );
406 // rotate to actual gamma direction
407 MuPlusDirection.rotateUz(GammaDirection);
408 MuMinusDirection.rotateUz(GammaDirection);
410 // create G4DynamicParticle object for the particle1
411 G4DynamicParticle* aParticle1 =
412 new G4DynamicParticle(theMuonPlus,MuPlusDirection,EPlus-Mmuon);
413 aParticleChange.AddSecondary(aParticle1);
414 // create G4DynamicParticle object for the particle2
415 G4DynamicParticle* aParticle2 =
416 new G4DynamicParticle(theMuonMinus,MuMinusDirection,EMinus-Mmuon);
417 aParticleChange.AddSecondary(aParticle2);
418 // Reset NbOfInteractionLengthLeft and return aParticleChange
419 return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep );
420}
421
422//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
423
424const G4Element* G4GammaConversionToMuons::SelectRandomAtom(
425 const G4DynamicParticle* aDynamicGamma,
426 const G4Material* aMaterial)
427{
428 // select randomly 1 element within the material, invoked by PostStepDoIt
429
430 const std::size_t NumberOfElements = aMaterial->GetNumberOfElements();
431 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
432 const G4Element* elm = (*theElementVector)[0];
433
434 if (NumberOfElements > 1) {
435 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
436
437 G4double PartialSumSigma = 0.;
438 G4double rval = G4UniformRand()/MeanFreePath;
439
440 for (std::size_t i=0; i<NumberOfElements; ++i)
441 {
442 elm = (*theElementVector)[i];
443 PartialSumSigma += NbOfAtomsPerVolume[i]
444 *GetCrossSectionPerAtom(aDynamicGamma, elm);
445 if (rval <= PartialSumSigma) { break; }
446 }
447 }
448 return elm;
449}
450
451//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
452
454{
455 G4String comments ="gamma->mu+mu- Bethe Heitler process, SubType= ";
456 G4cout << G4endl << GetProcessName() << ": " << comments
458 G4cout << " good cross section parametrization from "
459 << G4BestUnit(LowestEnergyLimit,"Energy")
460 << " to " << HighestEnergyLimit/GeV << " GeV for all Z." << G4endl;
461}
462
463//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double B(G4double temperature)
std::vector< const G4Element * > G4ElementVector
@ fGammaConversionToMuMu
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4ForceCondition
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const double C2
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *fvect, const G4MaterialCutsCouple *couple, const G4DynamicParticle *aDynamicGamma, G4double, G4double) override
void SetLeptonPair(const G4ParticleDefinition *p1, const G4ParticleDefinition *p2)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition: G4Element.hh:132
static G4EmParameters * Instance()
G4double MaxEnergyFor5DMuPair() const
G4double ComputeMeanFreePath(G4double GammaEnergy, const G4Material *aMaterial)
G4double ComputeCrossSectionPerAtom(G4double GammaEnergy, G4int Z)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
G4bool IsApplicable(const G4ParticleDefinition &) override
G4double GetCrossSectionPerAtom(const G4DynamicParticle *aDynamicGamma, const G4Element *anElement)
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4GammaConversionToMuons(const G4String &processName="GammaToMuPair", G4ProcessType type=fElectromagnetic)
static G4LossTableManager * Instance()
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
G4double GetA27(G4int Z) const
G4double GetZ13(G4double Z) const
static G4NistManager * Instance()
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
Definition: G4Step.hh:62
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define W
Definition: crc32.c:84
Definition: DoubConv.h:17
#define DBL_MAX
Definition: templates.hh:62