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

#include <G4PairProductionRelModel.hh>

+ Inheritance diagram for G4PairProductionRelModel:

Public Member Functions

 G4PairProductionRelModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheHeitlerLPM")
 
 ~G4PairProductionRelModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
void SetLPMflag (G4bool val)
 
G4bool LPMflag () const
 
G4PairProductionRelModeloperator= (const G4PairProductionRelModel &right)=delete
 
 G4PairProductionRelModel (const G4PairProductionRelModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

void ComputePhi12 (const G4double delta, G4double &phi1, G4double &phi2)
 
G4double ScreenFunction1 (const G4double delta)
 
G4double ScreenFunction2 (const G4double delta)
 
void ScreenFunction12 (const G4double delta, G4double &f1, G4double &f2)
 
G4double ComputeParametrizedXSectionPerAtom (G4double gammaEnergy, G4double Z)
 
G4double ComputeXSectionPerAtom (G4double gammaEnergy, G4double Z)
 
G4double ComputeDXSectionPerAtom (G4double eplusEnergy, G4double gammaEnergy, G4double Z)
 
G4double ComputeRelDXSectionPerAtom (G4double eplusEnergy, G4double gammaEnergy, G4double Z)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4bool fIsUseLPMCorrection
 
G4bool fIsUseCompleteScreening
 
G4double fLPMEnergy
 
G4double fParametrizedXSectionThreshold
 
G4double fCoulombCorrectionThreshold
 
G4PowfG4Calc
 
G4ParticleDefinitionfTheGamma
 
G4ParticleDefinitionfTheElectron
 
G4ParticleDefinitionfThePositron
 
G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Static Protected Attributes

static const G4int gMaxZet = 120
 
static const G4double gLPMconstant
 
static const G4double gXGL [8]
 
static const G4double gWGL [8]
 
static const G4double gFelLowZet [8]
 
static const G4double gFinelLowZet [8]
 
static const G4double gXSecFactor
 
static const G4double gEgLPMActivation = 100.*CLHEP::GeV
 
static std::vector< ElementData * > gElementData
 
static LPMFuncs gLPMFuncs
 

Detailed Description

Definition at line 65 of file G4PairProductionRelModel.hh.

Constructor & Destructor Documentation

◆ G4PairProductionRelModel() [1/2]

G4PairProductionRelModel::G4PairProductionRelModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "BetheHeitlerLPM" 
)
explicit

Definition at line 114 of file G4PairProductionRelModel.cc.

119 fParticleChange(nullptr)
120{
121 // gamma energy below which the parametrized atomic x-section is used (80 GeV)
122 fParametrizedXSectionThreshold = 30.0*CLHEP::GeV;
123 // gamma energy below the Coulomb correction is turned off (50 MeV)
124 fCoulombCorrectionThreshold = 50.0*CLHEP::MeV;
125 // set angular generator used in the final state kinematics computation
127}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4ParticleDefinition * fThePositron
G4ParticleDefinition * fTheGamma
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * fTheElectron
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607

◆ ~G4PairProductionRelModel()

G4PairProductionRelModel::~G4PairProductionRelModel ( )
override

Definition at line 130 of file G4PairProductionRelModel.cc.

131{
132 if (IsMaster()) {
133 // clear ElementData container
134 for (std::size_t iz = 0; iz < gElementData.size(); ++iz) {
135 if (gElementData[iz]) delete gElementData[iz];
136 }
137 gElementData.clear();
138 // clear LPMFunctions (if any)
140 gLPMFuncs.fLPMFuncG.clear();
141 gLPMFuncs.fLPMFuncPhi.clear();
142 gLPMFuncs.fIsInitialized = false;
143 }
144 }
145}
static std::vector< ElementData * > gElementData
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4PairProductionRelModel() [2/2]

G4PairProductionRelModel::G4PairProductionRelModel ( const G4PairProductionRelModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PairProductionRelModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0.,
G4double  cut = 0.,
G4double  emax = DBL_MAX 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 311 of file G4PairProductionRelModel.cc.

313{
314 G4double crossSection = 0.0 ;
315 // check kinematical limit
316 if ( gammaEnergy <= 2.0*electron_mass_c2 ) { return crossSection; }
317 // compute the atomic cross section either by using x-section parametrization
318 // or by numerically integrationg the DCS (with or without LPM)
319 if ( gammaEnergy < fParametrizedXSectionThreshold) {
320 // using the parametrized cross sections (max up to 80 GeV)
321 crossSection = ComputeParametrizedXSectionPerAtom(gammaEnergy, Z);
322 } else {
323 // by numerical integration of the DCS:
324 // Computes the cross section with or without LPM suppression depending on
325 // settings (by default with if the gamma energy is above a given threshold)
326 // and using or not using complete sreening approximation (by default not).
327 // Only the dependent part is computed in the numerical integration of the DCS
328 // i.e. the result must be multiplied here with 4 \alpha r_0^2 Z(Z+\eta(Z))
329 crossSection = ComputeXSectionPerAtom(gammaEnergy, Z);
330 // apply the constant factors:
331 // - eta(Z) is a correction to account interaction in the field of e-
332 // - gXSecFactor = 4 \alpha r_0^2
333 const G4int iz = std::min(gMaxZet, G4lrint(Z));
334 const G4double eta = gElementData[iz]->fEtaValue;
335 crossSection *= gXSecFactor*Z*(Z+eta);
336 }
337 // final protection
338 return std::max(crossSection, 0.);
339}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4double ComputeParametrizedXSectionPerAtom(G4double gammaEnergy, G4double Z)
G4double ComputeXSectionPerAtom(G4double gammaEnergy, G4double Z)
int G4lrint(double ad)
Definition: templates.hh:134

◆ ComputeDXSectionPerAtom()

G4double G4PairProductionRelModel::ComputeDXSectionPerAtom ( G4double  eplusEnergy,
G4double  gammaEnergy,
G4double  Z 
)
protected

Definition at line 228 of file G4PairProductionRelModel.cc.

231{
232 G4double xSection = 0.;
233 const G4int iz = std::min(gMaxZet, G4lrint(Z));
234 const G4double eps = pEnergy/gammaEnergy;
235 const G4double epsm = 1.-eps;
236 const G4double dum = eps*epsm;
238 // complete screening:
239 const G4double Lel = gElementData[iz]->fLradEl;
240 const G4double fc = gElementData[iz]->fCoulomb;
241 xSection = (eps*eps + epsm*epsm + 2.*dum/3.)*(Lel-fc) - dum/9.;
242 } else {
243 // normal case:
244 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
245 const G4double fc = gElementData[iz]->fCoulomb;
246 const G4double lnZ13 = gElementData[iz]->fLogZ13;
247 const G4double delta = gElementData[iz]->fDeltaFactor*eps0/dum;
248 G4double phi1, phi2;
249 ComputePhi12(delta, phi1, phi2);
250 xSection = (eps*eps + epsm*epsm)*(0.25*phi1-lnZ13-fc)
251 + 2.*dum*(0.25*phi2-lnZ13-fc)/3.;
252 }
253 // non-const. part of the DCS differential in total energy transfer not in eps
254 // ds/dEt=ds/deps deps/dEt with deps/dEt=1/Eg
255 return std::max(xSection, 0.0)/gammaEnergy;
256}
void ComputePhi12(const G4double delta, G4double &phi1, G4double &phi2)

Referenced by ComputeXSectionPerAtom().

◆ ComputeParametrizedXSectionPerAtom()

G4double G4PairProductionRelModel::ComputeParametrizedXSectionPerAtom ( G4double  gammaEnergy,
G4double  Z 
)
protected

Definition at line 654 of file G4PairProductionRelModel.cc.

656{
657 G4double xSection = 0.0 ;
658 // short versions
659 static const G4double kMC2 = CLHEP::electron_mass_c2;
660 // zero cross section below the kinematical limit: Eg<2mc^2
661 if (Z < 0.9 || gammaE <= 2.0*kMC2) { return xSection; }
662 //
663 static const G4double gammaEnergyLimit = 1.5*CLHEP::MeV;
664 // set coefficients a, b c
665 static const G4double a0 = 8.7842e+2*CLHEP::microbarn;
666 static const G4double a1 = -1.9625e+3*CLHEP::microbarn;
667 static const G4double a2 = 1.2949e+3*CLHEP::microbarn;
668 static const G4double a3 = -2.0028e+2*CLHEP::microbarn;
669 static const G4double a4 = 1.2575e+1*CLHEP::microbarn;
670 static const G4double a5 = -2.8333e-1*CLHEP::microbarn;
671
672 static const G4double b0 = -1.0342e+1*CLHEP::microbarn;
673 static const G4double b1 = 1.7692e+1*CLHEP::microbarn;
674 static const G4double b2 = -8.2381 *CLHEP::microbarn;
675 static const G4double b3 = 1.3063 *CLHEP::microbarn;
676 static const G4double b4 = -9.0815e-2*CLHEP::microbarn;
677 static const G4double b5 = 2.3586e-3*CLHEP::microbarn;
678
679 static const G4double c0 = -4.5263e+2*CLHEP::microbarn;
680 static const G4double c1 = 1.1161e+3*CLHEP::microbarn;
681 static const G4double c2 = -8.6749e+2*CLHEP::microbarn;
682 static const G4double c3 = 2.1773e+2*CLHEP::microbarn;
683 static const G4double c4 = -2.0467e+1*CLHEP::microbarn;
684 static const G4double c5 = 6.5372e-1*CLHEP::microbarn;
685 // check low energy limit of the approximation (1.5 MeV)
686 G4double gammaEnergyOrg = gammaE;
687 if (gammaE < gammaEnergyLimit) { gammaE = gammaEnergyLimit; }
688 // compute gamma energy variables
689 const G4double x = G4Log(gammaE/kMC2);
690 const G4double x2 = x *x;
691 const G4double x3 = x2*x;
692 const G4double x4 = x3*x;
693 const G4double x5 = x4*x;
694 //
695 const G4double F1 = a0 + a1*x + a2*x2 + a3*x3 + a4*x4 + a5*x5;
696 const G4double F2 = b0 + b1*x + b2*x2 + b3*x3 + b4*x4 + b5*x5;
697 const G4double F3 = c0 + c1*x + c2*x2 + c3*x3 + c4*x4 + c5*x5;
698 // compute the approximated cross section
699 xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
700 // check if we are below the limit of the approximation and apply correction
701 if (gammaEnergyOrg < gammaEnergyLimit) {
702 const G4double dum = (gammaEnergyOrg-2.*kMC2)/(gammaEnergyLimit-2.*kMC2);
703 xSection *= dum*dum;
704 }
705 return xSection;
706}
const G4double a0
G4double G4Log(G4double x)
Definition: G4Log.hh:227

Referenced by ComputeCrossSectionPerAtom().

◆ ComputePhi12()

void G4PairProductionRelModel::ComputePhi12 ( const G4double  delta,
G4double phi1,
G4double phi2 
)
inlineprotected

Definition at line 200 of file G4PairProductionRelModel.hh.

203{
204 if (delta > 1.4) {
205 phi1 = 21.0190 - 4.145*G4Log(delta + 0.958);
206 phi2 = phi1;
207 } else {
208 phi1 = 20.806 - delta*(3.190 - 0.5710*delta);
209 phi2 = 20.234 - delta*(2.126 - 0.0903*delta);
210 }
211}

Referenced by ComputeDXSectionPerAtom(), ComputeRelDXSectionPerAtom(), and SampleSecondaries().

◆ ComputeRelDXSectionPerAtom()

G4double G4PairProductionRelModel::ComputeRelDXSectionPerAtom ( G4double  eplusEnergy,
G4double  gammaEnergy,
G4double  Z 
)
protected

Definition at line 277 of file G4PairProductionRelModel.cc.

280{
281 G4double xSection = 0.;
282 const G4int iz = std::min(gMaxZet, G4lrint(Z));
283 const G4double eps = pEnergy/gammaEnergy;
284 const G4double epsm = 1.-eps;
285 const G4double dum = eps*epsm;
286 // evaluate LPM suppression functions
287 G4double fXiS, fGS, fPhiS;
288 ComputeLPMfunctions(fXiS, fGS, fPhiS, eps, gammaEnergy, iz);
290 // complete screening:
291 const G4double Lel = gElementData[iz]->fLradEl;
292 const G4double fc = gElementData[iz]->fCoulomb;
293 xSection = (Lel-fc)*((eps*eps+epsm*epsm)*2.*fPhiS + fGS)/3. - dum*fGS/9.;
294 } else {
295 // normal case:
296 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
297 const G4double fc = gElementData[iz]->fCoulomb;
298 const G4double lnZ13 = gElementData[iz]->fLogZ13;
299 const G4double delta = gElementData[iz]->fDeltaFactor*eps0/dum;
300 G4double phi1, phi2;
301 ComputePhi12(delta, phi1, phi2);
302 xSection = (eps*eps + epsm*epsm)*(2.*fPhiS+fGS)*(0.25*phi1-lnZ13-fc)/3.
303 + 2.*dum*fGS*(0.25*phi2-lnZ13-fc)/3.;
304 }
305 // non-const. part of the DCS differential in total energy transfer not in eps
306 // ds/dEt=ds/deps deps/dEt with deps/dEt=1/Eg
307 return std::max(fXiS*xSection, 0.0)/gammaEnergy;
308}

Referenced by ComputeXSectionPerAtom().

◆ ComputeXSectionPerAtom()

G4double G4PairProductionRelModel::ComputeXSectionPerAtom ( G4double  gammaEnergy,
G4double  Z 
)
protected

Definition at line 173 of file G4PairProductionRelModel.cc.

175{
176 G4double xSection = 0.0;
177 // check if LPM suppression needs to be used
178 const G4bool isLPM = (fIsUseLPMCorrection && gammaEnergy>gEgLPMActivation);
179 // determine the kinematical limits (taken into account the correction due to
180 // the way in which the Coulomb correction is applied i.e. avoid negative DCS)
181 const G4int iz = std::min(gMaxZet, G4lrint(Z));
182 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
183 // Coulomb correction is always included in the DCS even below 50 MeV (note:
184 // that this DCS is only used to get the integrated x-section)
185 const G4double dmax = gElementData[iz]->fDeltaMaxHigh;
186 const G4double dmin = 4.*eps0*gElementData[iz]->fDeltaFactor;
187 const G4double eps1 = 0.5 - 0.5*std::sqrt(1.-dmin/dmax);
188 const G4double epsMin = std::max(eps0, eps1);
189 const G4double epsMax = 0.5; // DCS is symmetric around eps=0.5
190 // let Et be the total energy transferred to the e- or to the e+
191 // the [Et-min, Et-max] interval will be divided into i=1,2,..,n subintervals
192 // with width of dInterv = (Et-max - Et-min)/n and numerical integration will
193 // be done in each sub-inteval using the xi = (Et - Et_i-min)/dInterv variable
194 // that is in [0,1]. The 8-point GL q. is used for the integration on [0,1].
195 const G4int numSub = 2;
196 const G4double dInterv= (epsMax - epsMin)*gammaEnergy/G4double(numSub);
197 G4double minEti = epsMin*gammaEnergy; // Et-min i.e. Et_0-min
198 for (G4int i = 0; i < numSub; ++i) {
199 for (G4int ngl = 0; ngl < 8; ++ngl) {
200 const G4double Et = (minEti + gXGL[ngl]*dInterv);
201 const G4double xs = isLPM ? ComputeRelDXSectionPerAtom(Et, gammaEnergy, Z)
202 : ComputeDXSectionPerAtom(Et, gammaEnergy, Z);
203 xSection += gWGL[ngl]*xs;
204 }
205 // update minimum Et of the sub-inteval
206 minEti += dInterv;
207 }
208 // apply corrections of variable transformation and half interval integration
209 xSection = std::max(2.*xSection*dInterv, 0.);
210 return xSection;
211}
bool G4bool
Definition: G4Types.hh:86
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
static const G4double gEgLPMActivation
static const G4double gXGL[8]
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double gammaEnergy, G4double Z)
static const G4double gWGL[8]

Referenced by ComputeCrossSectionPerAtom().

◆ Initialise()

void G4PairProductionRelModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 147 of file G4PairProductionRelModel.cc.

149{
150 if (IsMaster()) {
151 // init element data and LPM funcs
152 if (IsMaster()) {
153 InitialiseElementData();
155 InitLPMFunctions();
156 }
157 }
158 }
162 }
163}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139

Referenced by G4LivermoreGammaConversionModel::Initialise(), and G4BetheHeitler5DModel::Initialise().

◆ InitialiseLocal()

void G4PairProductionRelModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 165 of file G4PairProductionRelModel.cc.

167{
170 }
171}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823

◆ LPMflag()

G4bool G4PairProductionRelModel::LPMflag ( ) const
inline

Definition at line 97 of file G4PairProductionRelModel.hh.

97{ return fIsUseLPMCorrection; }

◆ operator=()

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

◆ SampleSecondaries()

void G4PairProductionRelModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Definition at line 348 of file G4PairProductionRelModel.cc.

365{
366 const G4Material* mat = couple->GetMaterial();
367 const G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
368 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy ;
369 //
370 // check kinematical limit: gamma energy(Eg) must be at least 2 e- rest mass
371 // (but the model should be used at higher energies above 100 MeV)
372 if (eps0 > 0.5) { return; }
373 //
374 // select target atom of the material
375 const G4Element* anElement = SelectTargetAtom(couple, fTheGamma, gammaEnergy,
376 aDynamicGamma->GetLogKineticEnergy());
377 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
378 //
379 // 'eps' is the total energy transferred to one of the e-/e+ pair in initial
380 // gamma energy units Eg. Since the corresponding DCS is symmetric on eps=0.5,
381 // the kinematical limits for eps0=mc^2/Eg <= eps <= 0.5
382 // 1. 'eps' is sampled uniformly on the [eps0, 0.5] inteval if Eg<Egsmall
383 // 2. otherwise, on the [eps_min, 0.5] interval according to the DCS (case 2.)
384 G4double eps;
385 // case 1.
386 static const G4double Egsmall = 2.*CLHEP::MeV;
387 if (gammaEnergy < Egsmall) {
388 eps = eps0 + (0.5-eps0)*rndmEngine->flat();
389 } else {
390 // case 2.
391 // get the Coulomb factor for the target element (Z) and gamma energy (Eg)
392 // F(Z) = 8*ln(Z)/3 if Eg <= 50 [MeV] => no Coulomb correction
393 // F(Z) = 8*ln(Z)/3 + 8*fc(Z) if Eg > 50 [MeV] => fc(Z) is the Coulomb cor.
394 //
395 // The screening variable 'delta(eps)' = 136*Z^{-1/3}*eps0/[eps(1-eps)]
396 // Due to the Coulomb correction, the DCS can go below zero even at
397 // kinematicaly allowed eps > eps0 values. In order to exclude this eps
398 // range with negative DCS, the minimum eps value will be set to eps_min =
399 // max[eps0, epsp] with epsp is the solution of SF(delta(epsp)) - F(Z)/2 = 0
400 // with SF being the screening function (SF1=SF2 at high value of delta).
401 // The solution is epsp = 0.5 - 0.5*sqrt[ 1 - 4*136*Z^{-1/3}eps0/deltap]
402 // with deltap = Exp[(42.038-F(Z))/8.29]-0.958. So the limits are:
403 // - when eps=eps_max = 0.5 => delta_min = 136*Z^{-1/3}*eps0/4
404 // - epsp = 0.5 - 0.5*sqrt[ 1 - delta_min/deltap]
405 // - and eps_min = max[eps0, epsp]
406 const G4int iZet = std::min(gMaxZet, anElement->GetZasInt());
407 const G4double deltaFactor = gElementData[iZet]->fDeltaFactor*eps0;
408 const G4double deltaMin = 4.*deltaFactor;
409 G4double deltaMax = gElementData[iZet]->fDeltaMaxLow;
410 G4double FZ = 8.*gElementData[iZet]->fLogZ13;
411 if ( gammaEnergy > fCoulombCorrectionThreshold ) { // Eg > 50 MeV ?
412 FZ += 8.*gElementData[iZet]->fCoulomb;
413 deltaMax = gElementData[iZet]->fDeltaMaxHigh;
414 }
415 // compute the limits of eps
416 const G4double epsp = 0.5 - 0.5*std::sqrt(1. - deltaMin/deltaMax) ;
417 const G4double epsMin = std::max(eps0,epsp);
418 const G4double epsRange = 0.5 - epsMin;
419 //
420 // sample the energy rate (eps) of the created electron (or positron)
422 ScreenFunction12(deltaMin, F10, F20);
423 F10 -= FZ;
424 F20 -= FZ;
425 const G4double NormF1 = std::max(F10 * epsRange * epsRange, 0.);
426 const G4double NormF2 = std::max(1.5 * F20 , 0.);
427 const G4double NormCond = NormF1/(NormF1 + NormF2);
428 // check if LPM correction is active
429 const G4bool isLPM = (fIsUseLPMCorrection && gammaEnergy>gEgLPMActivation);
431 // we will need 3 uniform random number for each trial of sampling
432 G4double rndmv[3];
433 G4double greject = 0.;
434 do {
435 rndmEngine->flatArray(3, rndmv);
436 if (NormCond > rndmv[0]) {
437 eps = 0.5 - epsRange * fG4Calc->A13(rndmv[1]);
438 const G4double delta = deltaFactor/(eps*(1.-eps));
439 if (isLPM) {
440 G4double lpmXiS, lpmGS, lpmPhiS, phi1, phi2;
441 ComputePhi12(delta, phi1, phi2);
442 ComputeLPMfunctions(lpmXiS, lpmGS, lpmPhiS, eps, gammaEnergy, iZet);
443 greject = lpmXiS*((2.*lpmPhiS+lpmGS)*phi1-lpmGS*phi2-lpmPhiS*FZ)/F10;
444 } else {
445 greject = (ScreenFunction1(delta)-FZ)/F10;
446 }
447 } else {
448 eps = epsMin + epsRange*rndmv[1];
449 const G4double delta = deltaFactor/(eps*(1.-eps));
450 if (isLPM) {
451 G4double lpmXiS, lpmGS, lpmPhiS, phi1, phi2;
452 ComputePhi12(delta, phi1, phi2);
453 ComputeLPMfunctions(lpmXiS, lpmGS, lpmPhiS, eps, gammaEnergy, iZet);
454 greject = lpmXiS*( (lpmPhiS+0.5*lpmGS)*phi1 + 0.5*lpmGS*phi2
455 -0.5*(lpmGS+lpmPhiS)*FZ )/F20;
456 } else {
457 greject = (ScreenFunction2(delta)-FZ)/F20;
458 }
459 }
460 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
461 } while (greject < rndmv[2]);
462 // end of eps sampling
463 }
464 //
465 // select charges randomly
466 G4double eTotEnergy, pTotEnergy;
467 if (rndmEngine->flat() > 0.5) {
468 eTotEnergy = (1.-eps)*gammaEnergy;
469 pTotEnergy = eps*gammaEnergy;
470 } else {
471 pTotEnergy = (1.-eps)*gammaEnergy;
472 eTotEnergy = eps*gammaEnergy;
473 }
474 //
475 // sample pair kinematics
476 //
477 const G4double eKinEnergy = std::max(0.,eTotEnergy - CLHEP::electron_mass_c2);
478 const G4double pKinEnergy = std::max(0.,pTotEnergy - CLHEP::electron_mass_c2);
479 //
480 G4ThreeVector eDirection, pDirection;
481 //
483 eKinEnergy, pKinEnergy, eDirection, pDirection);
484 // create G4DynamicParticle object for the particle1
485 auto aParticle1 = new G4DynamicParticle(fTheElectron,eDirection,eKinEnergy);
486
487 // create G4DynamicParticle object for the particle2
488 auto aParticle2 = new G4DynamicParticle(fThePositron,pDirection,pKinEnergy);
489 // Fill output vector
490 fvect->push_back(aParticle1);
491 fvect->push_back(aParticle2);
492 // kill incident photon
495}
#define F10
#define F20
@ fStopAndKill
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition: G4Element.hh:132
const G4Material * GetMaterial() const
G4double GetRadlen() const
Definition: G4Material.hh:215
void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2)
G4double ScreenFunction1(const G4double delta)
static const G4double gLPMconstant
G4double ScreenFunction2(const G4double delta)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double A13(G4double A) const
Definition: G4Pow.cc:116
virtual void SamplePairDirections(const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
void ProposeTrackStatus(G4TrackStatus status)

◆ ScreenFunction1()

G4double G4PairProductionRelModel::ScreenFunction1 ( const G4double  delta)
inlineprotected

Definition at line 214 of file G4PairProductionRelModel.hh.

215{
216 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
217 : 42.184 - delta*(7.444 - 1.623*delta);
218}

Referenced by SampleSecondaries().

◆ ScreenFunction12()

void G4PairProductionRelModel::ScreenFunction12 ( const G4double  delta,
G4double f1,
G4double f2 
)
inlineprotected

Definition at line 228 of file G4PairProductionRelModel.hh.

230{
231 if (delta > 1.4) {
232 f1 = 42.038 - 8.29*G4Log(delta + 0.958);
233 f2 = f1;
234 } else {
235 f1 = 42.184 - delta*(7.444 - 1.623*delta);
236 f2 = 41.326 - delta*(5.848 - 0.902*delta);
237 }
238}

Referenced by SampleSecondaries().

◆ ScreenFunction2()

G4double G4PairProductionRelModel::ScreenFunction2 ( const G4double  delta)
inlineprotected

Definition at line 221 of file G4PairProductionRelModel.hh.

222{
223 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
224 : 41.326 - delta*(5.848 - 0.902*delta);
225}

Referenced by SampleSecondaries().

◆ SetLPMflag()

void G4PairProductionRelModel::SetLPMflag ( G4bool  val)
inline

Definition at line 96 of file G4PairProductionRelModel.hh.

96{ fIsUseLPMCorrection = val; }

◆ SetupForMaterial()

void G4PairProductionRelModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material mat,
G4double   
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 341 of file G4PairProductionRelModel.cc.

343{
345}

Member Data Documentation

◆ fCoulombCorrectionThreshold

G4double G4PairProductionRelModel::fCoulombCorrectionThreshold
protected

Definition at line 175 of file G4PairProductionRelModel.hh.

Referenced by G4PairProductionRelModel(), and SampleSecondaries().

◆ fG4Calc

G4Pow* G4PairProductionRelModel::fG4Calc
protected

Definition at line 177 of file G4PairProductionRelModel.hh.

Referenced by SampleSecondaries().

◆ fIsUseCompleteScreening

G4bool G4PairProductionRelModel::fIsUseCompleteScreening
protected

◆ fIsUseLPMCorrection

G4bool G4PairProductionRelModel::fIsUseLPMCorrection
protected

◆ fLPMEnergy

G4double G4PairProductionRelModel::fLPMEnergy
protected

Definition at line 172 of file G4PairProductionRelModel.hh.

Referenced by SampleSecondaries(), and SetupForMaterial().

◆ fParametrizedXSectionThreshold

G4double G4PairProductionRelModel::fParametrizedXSectionThreshold
protected

◆ fParticleChange

G4ParticleChangeForGamma* G4PairProductionRelModel::fParticleChange
protected

◆ fTheElectron

G4ParticleDefinition* G4PairProductionRelModel::fTheElectron
protected

◆ fTheGamma

G4ParticleDefinition* G4PairProductionRelModel::fTheGamma
protected

◆ fThePositron

G4ParticleDefinition* G4PairProductionRelModel::fThePositron
protected

Definition at line 180 of file G4PairProductionRelModel.hh.

Referenced by SampleSecondaries().

◆ gEgLPMActivation

const G4double G4PairProductionRelModel::gEgLPMActivation = 100.*CLHEP::GeV
staticprotected

Definition at line 164 of file G4PairProductionRelModel.hh.

Referenced by ComputeXSectionPerAtom(), and SampleSecondaries().

◆ gElementData

std::vector< G4PairProductionRelModel::ElementData * > G4PairProductionRelModel::gElementData
staticprotected

◆ gFelLowZet

const G4double G4PairProductionRelModel::gFelLowZet
staticprotected
Initial value:
= {
0.0, 5.3104, 4.7935, 4.7402, 4.7112, 4.6694, 4.6134, 4.5520
}

Definition at line 160 of file G4PairProductionRelModel.hh.

◆ gFinelLowZet

const G4double G4PairProductionRelModel::gFinelLowZet
staticprotected
Initial value:
= {
0.0, 5.9173, 5.6125, 5.5377, 5.4728, 5.4174, 5.3688, 5.3236
}

Definition at line 161 of file G4PairProductionRelModel.hh.

◆ gLPMconstant

const G4double G4PairProductionRelModel::gLPMconstant
staticprotected
Initial value:
=
CLHEP::fine_structure_const*CLHEP::electron_mass_c2*CLHEP::electron_mass_c2
/(4.*CLHEP::pi*CLHEP::hbarc)

Definition at line 156 of file G4PairProductionRelModel.hh.

Referenced by SampleSecondaries(), and SetupForMaterial().

◆ gLPMFuncs

G4PairProductionRelModel::LPMFuncs G4PairProductionRelModel::gLPMFuncs
staticprotected

Definition at line 167 of file G4PairProductionRelModel.hh.

Referenced by ~G4PairProductionRelModel().

◆ gMaxZet

const G4int G4PairProductionRelModel::gMaxZet = 120
staticprotected

◆ gWGL

const G4double G4PairProductionRelModel::gWGL
staticprotected
Initial value:
= {
5.06142681e-02, 1.11190517e-01, 1.56853323e-01, 1.81341892e-01,
1.81341892e-01, 1.56853323e-01, 1.11190517e-01, 5.06142681e-02
}

Definition at line 159 of file G4PairProductionRelModel.hh.

Referenced by ComputeXSectionPerAtom().

◆ gXGL

const G4double G4PairProductionRelModel::gXGL
staticprotected
Initial value:
= {
1.98550718e-02, 1.01666761e-01, 2.37233795e-01, 4.08282679e-01,
5.91717321e-01, 7.62766205e-01, 8.98333239e-01, 9.80144928e-01
}

Definition at line 158 of file G4PairProductionRelModel.hh.

Referenced by ComputeXSectionPerAtom().

◆ gXSecFactor

const G4double G4PairProductionRelModel::gXSecFactor
staticprotected
Initial value:
=
4.*CLHEP::fine_structure_const*CLHEP::classic_electr_radius
*CLHEP::classic_electr_radius

Definition at line 163 of file G4PairProductionRelModel.hh.

Referenced by ComputeCrossSectionPerAtom().


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