Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaConversionToMuons Class Reference

#include <G4GammaConversionToMuons.hh>

+ Inheritance diagram for G4GammaConversionToMuons:

Public Member Functions

 G4GammaConversionToMuons (const G4String &processName="GammaToMuPair", G4ProcessType type=fElectromagnetic)
 
 ~G4GammaConversionToMuons () override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void PrintInfoDefinition ()
 
void SetCrossSecFactor (G4double fac)
 
G4double GetCrossSecFactor () const
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double GetCrossSectionPerAtom (const G4DynamicParticle *aDynamicGamma, const G4Element *anElement)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4double ComputeCrossSectionPerAtom (G4double GammaEnergy, G4int Z)
 
G4double ComputeMeanFreePath (G4double GammaEnergy, const G4Material *aMaterial)
 
G4GammaConversionToMuonsoperator= (const G4GammaConversionToMuons &right)=delete
 
 G4GammaConversionToMuons (const G4GammaConversionToMuons &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 60 of file G4GammaConversionToMuons.hh.

Constructor & Destructor Documentation

◆ G4GammaConversionToMuons() [1/2]

G4GammaConversionToMuons::G4GammaConversionToMuons ( const G4String processName = "GammaToMuPair",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 60 of file G4GammaConversionToMuons.cc.

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}
@ fGammaConversionToMuMu
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4GammaConversionToMuons()

G4GammaConversionToMuons::~G4GammaConversionToMuons ( )
override

Definition at line 79 of file G4GammaConversionToMuons.cc.

80{
81 fManager->DeRegister(this);
82}
void DeRegister(G4VEnergyLossProcess *p)

◆ G4GammaConversionToMuons() [2/2]

G4GammaConversionToMuons::G4GammaConversionToMuons ( const G4GammaConversionToMuons )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4GammaConversionToMuons::BuildPhysicsTable ( const G4ParticleDefinition p)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 93 of file G4GammaConversionToMuons.cc.

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}
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetLeptonPair(const G4ParticleDefinition *p1, const G4ParticleDefinition *p2)
static G4EmParameters * Instance()
G4double MaxEnergyFor5DMuPair() const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()

Referenced by G4GammaGeneralProcess::BuildPhysicsTable().

◆ ComputeCrossSectionPerAtom()

G4double G4GammaConversionToMuons::ComputeCrossSectionPerAtom ( G4double  GammaEnergy,
G4int  Z 
)

Definition at line 166 of file G4GammaConversionToMuons.cc.

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}
G4double B(G4double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
G4double GetA27(G4int Z) const
G4double GetZ13(G4double Z) const
static G4NistManager * Instance()

Referenced by ComputeMeanFreePath(), and GetCrossSectionPerAtom().

◆ ComputeMeanFreePath()

G4double G4GammaConversionToMuons::ComputeMeanFreePath ( G4double  GammaEnergy,
const G4Material aMaterial 
)

Definition at line 125 of file G4GammaConversionToMuons.cc.

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}
std::vector< const G4Element * > G4ElementVector
G4double ComputeCrossSectionPerAtom(G4double GammaEnergy, G4int Z)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
#define DBL_MAX
Definition: templates.hh:62

Referenced by G4GammaGeneralProcess::BuildPhysicsTable(), and GetMeanFreePath().

◆ GetCrossSecFactor()

G4double G4GammaConversionToMuons::GetCrossSecFactor ( ) const
inline

Definition at line 85 of file G4GammaConversionToMuons.hh.

85{ return CrossSecFactor;}

◆ GetCrossSectionPerAtom()

G4double G4GammaConversionToMuons::GetCrossSectionPerAtom ( const G4DynamicParticle aDynamicGamma,
const G4Element anElement 
)

Definition at line 154 of file G4GammaConversionToMuons.cc.

159{
160 return ComputeCrossSectionPerAtom(aDynamicGamma->GetKineticEnergy(),
161 anElement->GetZasInt());
162}
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition: G4Element.hh:132

◆ GetMeanFreePath()

G4double G4GammaConversionToMuons::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 109 of file G4GammaConversionToMuons.cc.

115{
116 const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
117 G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
118 const G4Material* aMaterial = aTrack.GetMaterial();
119 return ComputeMeanFreePath(GammaEnergy, aMaterial);
120}
G4double ComputeMeanFreePath(G4double GammaEnergy, const G4Material *aMaterial)
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const

◆ IsApplicable()

G4bool G4GammaConversionToMuons::IsApplicable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 86 of file G4GammaConversionToMuons.cc.

87{
88 return (&part == theGamma);
89}

◆ operator=()

G4GammaConversionToMuons & G4GammaConversionToMuons::operator= ( const G4GammaConversionToMuons right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4GammaConversionToMuons::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 219 of file G4GammaConversionToMuons.cc.

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}
@ fStopAndKill
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const double C2
#define G4UniformRand()
Definition: Randomize.hh:52
void SampleSecondaries(std::vector< G4DynamicParticle * > *fvect, const G4MaterialCutsCouple *couple, const G4DynamicParticle *aDynamicGamma, G4double, G4double) override
const G4ThreeVector & GetMomentumDirection() const
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
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
#define W
Definition: crc32.c:84
const G4double pi

◆ PrintInfoDefinition()

void G4GammaConversionToMuons::PrintInfoDefinition ( )

Definition at line 453 of file G4GammaConversionToMuons.cc.

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}
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

Referenced by BuildPhysicsTable().

◆ SetCrossSecFactor()

void G4GammaConversionToMuons::SetCrossSecFactor ( G4double  fac)

Definition at line 208 of file G4GammaConversionToMuons.cc.

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}

Referenced by G4EmExtraPhysics::ConstructProcess().


The documentation for this class was generated from the following files: