57static const G4double sqrte=sqrt(exp(1.));
64 Rc(elm_coupling/Mmuon),
65 LimitEnergy (5.*Mmuon),
66 LowestEnergyLimit (2.*Mmuon),
67 HighestEnergyLimit(1e12*GeV),
92 return (&part == theGamma);
101 if(Energy5DLimit > 0.0 && !f5Dmodel) {
124 MeanFreePath = (GammaEnergy <= LowestEnergyLimit)
138 if(GammaEnergy <= LowestEnergyLimit) {
return DBL_MAX; }
146 if(e < LimitEnergy) {
147 G4double y = (e - LowestEnergyLimit)/(LimitEnergy - LowestEnergyLimit);
154 SIGMA += NbOfAtomsPerVolume[i] * fact *
157 return (SIGMA > 0.0) ? 1./SIGMA :
DBL_MAX;
181 if(Egam <= LowestEnergyLimit) {
return 0.0; }
185 G4double PowThres,Ecor,
B,Dn,Zthird,Winfty,WMedAppr,
196 Zthird=1./nist->
GetZ13(Z);
197 Winfty=
B*Zthird*Mmuon/(Dn*electron_mass_c2);
198 WMedAppr=1./(4.*Dn*sqrte*Mmuon);
199 Wsatur=Winfty/WMedAppr;
200 sigfac=4.*fine_structure_const*Z*Z*Rc*Rc;
201 PowThres=1.479+0.00799*Dn;
202 Ecor=-18.+4347./(
B*Zthird);
209 G4double CrossSection=7./9.*sigfac*
G4Log(1.+WMedAppr*CorFuc*Eg);
210 CrossSection*=CrossSecFactor;
220 G4cout <<
"The cross section for GammaConversionToMuons is artificially "
221 <<
"increased by the CrossSecFactor=" << CrossSecFactor <<
G4endl;
239 if (Egam <= LowestEnergyLimit) {
249 if (Egam <= Energy5DLimit) {
250 std::vector<G4DynamicParticle*> fvect;
261 const G4Element* anElement = SelectRandomAtom(aDynamicGamma, aMaterial);
277 G4double Winfty=
B*Zthird*Mmuon/(Dn*electron_mass_c2);
281 G4double C2Term2=electron_mass_c2/(183.*Zthird*Mmuon);
286 G4double xmin=(Egam < LimitEnergy) ? GammaMuonInv : .5-sqrt(.25-GammaMuonInv);
290 G4double sBZ=sqrte*
B*Zthird/electron_mass_c2;
292 /(1.+2.*sBZ*Mmuon*GammaMuonInv));
295 const G4int nmax = 1000;
300 G4double del=Mmuon*Mmuon/(2.*Egam*xPM);
301 W=Winfty*(1.+Ds2*del/Mmuon)/(1.+sBZ*del);
303 result=(xxp > 0.) ? xxp*
G4Log(W)*LogWmaxInv : 0.0;
305 G4cout <<
"G4GammaConversionToMuons::PostStepDoIt WARNING:"
306 <<
" in dSigxPlusGen, result=" << result <<
" > 1" <<
G4endl;
309 if(nn >= nmax) {
break; }
319 G4double a3 = (GammaMuonInv/(2.*xPM));
326 G4double f1_max=-(1.-xPM)*(2.*b1+(a21+a33)*
G4Log(a33/a21))/(2*b3);
328 G4double thetaPlus,thetaMinus,phiHalf;
338 if(std::abs(b1)<0.0001*a34)
341 f1=(1.-2.*xPM+4.*xPM*t*(1.-t))/(12.*a34*a34*a34*a34);
345 f1=-(1.-2.*xPM+4.*xPM*t*(1.-t))*(2.*b1+(a22+a34)*
G4Log(a34/a22))/(2*b3);
347 if(f1<0.0 || f1> f1_max)
349 G4cout <<
"G4GammaConversionToMuons::PostStepDoIt WARNING:"
350 <<
"outside allowed range f1=" << f1
351 <<
" is set to zero, a34 = "<< a34 <<
" a22 = "<<a22<<
"."
355 if(nn > nmax) {
break; }
359 G4double f2_max=1.-2.*xPM*(1.-4.*t*(1.-t));
365 f2=1.-2.*xPM+4.*xPM*t*(1.-t)*(1.+cos(2.*psi));
366 if(f2<0 || f2> f2_max)
368 G4cout <<
"G4GammaConversionToMuons::PostStepDoIt WARNING:"
369 <<
"outside allowed range f2=" << f2 <<
" is set to zero"
373 if(nn >= nmax) {
break; }
378 G4double C2Term1=GammaMuonInv/(2.*xPM*t);
379 G4double C22 = C2Term1*C2Term1+C2Term2*C2Term2;
388 phiHalf=0.5*rho/u*sin(psi);
390 thetaPlus =GammaMuonInv*(u+xiHalf)/xPlus;
391 thetaMinus=GammaMuonInv*(u-xiHalf)/xMinus;
395 if(std::abs(thetaPlus)>pi) { thetaPlus = 0.0; }
396 if(std::abs(thetaMinus)>pi) { thetaMinus = 0.0; }
400 }
while ( std::abs(thetaPlus)>pi || std::abs(thetaMinus) >pi);
409 G4ThreeVector MuPlusDirection ( sin(thetaPlus) *cos(phi0+phiHalf),
410 sin(thetaPlus) *sin(phi0+phiHalf), cos(thetaPlus) );
411 G4ThreeVector MuMinusDirection (-sin(thetaMinus)*cos(phi0-phiHalf),
412 -sin(thetaMinus) *sin(phi0-phiHalf), cos(thetaMinus) );
414 MuPlusDirection.
rotateUz(GammaDirection);
415 MuMinusDirection.
rotateUz(GammaDirection);
431const G4Element* G4GammaConversionToMuons::SelectRandomAtom(
439 const G4Element* elm = (*theElementVector)[0];
441 if (NumberOfElements > 1) {
447 for (
G4int i=0; i<NumberOfElements; ++i)
449 elm = (*theElementVector)[i];
450 PartialSumSigma += NbOfAtomsPerVolume[i]
452 if (rval <= PartialSumSigma) {
break; }
462 G4String comments =
"gamma->mu+mu- Bethe Heitler process, SubType= ";
465 G4cout <<
" good cross section parametrization from "
467 <<
" to " << HighestEnergyLimit/GeV <<
" GeV for all Z." <<
G4endl;
double B(double temperature)
std::vector< G4Element * > G4ElementVector
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4double G4Log(G4double x)
G4GLOB_DLL std::ostream G4cout
Hep3Vector & rotateUz(const Hep3Vector &)
virtual 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
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
void PrintInfoDefinition()
void SetCrossSecFactor(G4double fac)
~G4GammaConversionToMuons() 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
size_t GetNumberOfElements() const
const G4double * GetVecNbOfAtomsPerVolume() const
G4double GetA27(G4int Z) const
G4double GetZ13(G4double Z) const
static G4NistManager * Instance()
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
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
void SetProcessSubType(G4int)
G4int GetProcessSubType() const
const G4String & GetProcessName() const