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

#include <G4BetheHeitlerModel.hh>

+ Inheritance diagram for G4BetheHeitlerModel:

Classes

struct  ElementData
 

Public Member Functions

 G4BetheHeitlerModel (const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitler")
 
virtual ~G4BetheHeitlerModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4double ScreenFunction1 (const G4double delta)
 
G4double ScreenFunction2 (const G4double delta)
 
void ScreenFunction12 (const G4double delta, G4double &f1, G4double &f2)
 
void InitialiseElementData ()
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4PowfG4Calc
 
G4ParticleDefinitionfTheGamma
 
G4ParticleDefinitionfTheElectron
 
G4ParticleDefinitionfThePositron
 
G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Static Protected Attributes

static const G4int gMaxZet = 120
 
static std::vector< ElementData * > gElementData
 

Detailed Description

Definition at line 60 of file G4BetheHeitlerModel.hh.

Constructor & Destructor Documentation

◆ G4BetheHeitlerModel()

G4BetheHeitlerModel::G4BetheHeitlerModel ( const G4ParticleDefinition p = 0,
const G4String nam = "BetheHeitler" 
)
explicit

Definition at line 60 of file G4BetheHeitlerModel.cc.

62: G4VEmModel(nam),
65 fParticleChange(nullptr)
66{
68}
G4ParticleDefinition * fTheElectron
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * fTheGamma
G4ParticleDefinition * fThePositron
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4BetheHeitlerModel()

G4BetheHeitlerModel::~G4BetheHeitlerModel ( )
virtual

Definition at line 70 of file G4BetheHeitlerModel.cc.

71{
72 if (IsMaster()) {
73 // clear ElementData container
74 for (size_t iz = 0; iz < gElementData.size(); ++iz) {
75 if (gElementData[iz]) delete gElementData[iz];
76 }
77 gElementData.clear();
78 }
79}
static std::vector< ElementData * > gElementData
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 106 of file G4BetheHeitlerModel.cc.

109{
110 G4double xSection = 0.0 ;
111 // short versions
112 static const G4double kMC2 = CLHEP::electron_mass_c2;
113 // zero cross section below the kinematical limit: Eg<2mc^2
114 if (Z < 0.9 || gammaEnergy <= 2.0*kMC2) { return xSection; }
115 //
116 static const G4double gammaEnergyLimit = 1.5*CLHEP::MeV;
117 // set coefficients a, b c
118 static const G4double a0 = 8.7842e+2*CLHEP::microbarn;
119 static const G4double a1 = -1.9625e+3*CLHEP::microbarn;
120 static const G4double a2 = 1.2949e+3*CLHEP::microbarn;
121 static const G4double a3 = -2.0028e+2*CLHEP::microbarn;
122 static const G4double a4 = 1.2575e+1*CLHEP::microbarn;
123 static const G4double a5 = -2.8333e-1*CLHEP::microbarn;
124
125 static const G4double b0 = -1.0342e+1*CLHEP::microbarn;
126 static const G4double b1 = 1.7692e+1*CLHEP::microbarn;
127 static const G4double b2 = -8.2381 *CLHEP::microbarn;
128 static const G4double b3 = 1.3063 *CLHEP::microbarn;
129 static const G4double b4 = -9.0815e-2*CLHEP::microbarn;
130 static const G4double b5 = 2.3586e-3*CLHEP::microbarn;
131
132 static const G4double c0 = -4.5263e+2*CLHEP::microbarn;
133 static const G4double c1 = 1.1161e+3*CLHEP::microbarn;
134 static const G4double c2 = -8.6749e+2*CLHEP::microbarn;
135 static const G4double c3 = 2.1773e+2*CLHEP::microbarn;
136 static const G4double c4 = -2.0467e+1*CLHEP::microbarn;
137 static const G4double c5 = 6.5372e-1*CLHEP::microbarn;
138 // check low energy limit of the approximation (1.5 MeV)
139 G4double gammaEnergyOrg = gammaEnergy;
140 if (gammaEnergy < gammaEnergyLimit) { gammaEnergy = gammaEnergyLimit; }
141 // compute gamma energy variables
142 const G4double x = G4Log(gammaEnergy/kMC2);
143 const G4double x2 = x *x;
144 const G4double x3 = x2*x;
145 const G4double x4 = x3*x;
146 const G4double x5 = x4*x;
147 //
148 const G4double F1 = a0 + a1*x + a2*x2 + a3*x3 + a4*x4 + a5*x5;
149 const G4double F2 = b0 + b1*x + b2*x2 + b3*x3 + b4*x4 + b5*x5;
150 const G4double F3 = c0 + c1*x + c2*x2 + c3*x3 + c4*x4 + c5*x5;
151 // compute the approximated cross section
152 xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
153 // check if we are below the limit of the approximation and apply correction
154 if (gammaEnergyOrg < gammaEnergyLimit) {
155 const G4double dum = (gammaEnergyOrg-2.*kMC2)/(gammaEnergyLimit-2.*kMC2);
156 xSection *= dum*dum;
157 }
158 // make sure that the cross section is never negative
159 xSection = std::max(xSection, 0.);
160 return xSection;
161}
const G4double a0
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83

◆ Initialise()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedGammaConversionModel.

Definition at line 81 of file G4BetheHeitlerModel.cc.

83{
84 if (IsMaster()) {
86 }
88 if (IsMaster()) {
90 }
91}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

Referenced by G4PolarizedGammaConversionModel::Initialise().

◆ InitialiseElementData()

void G4BetheHeitlerModel::InitialiseElementData ( )
protected

Definition at line 298 of file G4BetheHeitlerModel.cc.

299{
300 G4int size = gElementData.size();
301 if (size < gMaxZet+1) {
302 gElementData.resize(gMaxZet+1, nullptr);
303 }
304 // create for all elements that are in the detector
305 const G4ElementTable* elemTable = G4Element::GetElementTable();
306 size_t numElems = (*elemTable).size();
307 for (size_t ie = 0; ie < numElems; ++ie) {
308 const G4Element* elem = (*elemTable)[ie];
309 const G4int iz = std::min(gMaxZet, elem->GetZasInt());
310 if (!gElementData[iz]) { // create it if doesn't exist yet
311 G4double FZLow = 8.*elem->GetIonisation()->GetlogZ3();
312 G4double FZHigh = FZLow + 8.*elem->GetfCoulomb();
313 ElementData* elD = new ElementData();
314 elD->fDeltaMaxLow = G4Exp((42.038 - FZLow )/8.29) - 0.958;
315 elD->fDeltaMaxHigh = G4Exp((42.038 - FZHigh)/8.29) - 0.958;
316 gElementData[iz] = elD;
317 }
318 }
319}
std::vector< G4Element * > G4ElementTable
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
int G4int
Definition: G4Types.hh:85
static const G4int gMaxZet
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
G4int GetZasInt() const
Definition: G4Element.hh:131
G4double GetlogZ3() const

Referenced by Initialise().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 93 of file G4BetheHeitlerModel.cc.

95{
97}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedGammaConversionModel.

Definition at line 175 of file G4BetheHeitlerModel.cc.

179{
180 // set some constant values
181 const G4double gammaEnergy = aDynamicGamma->GetKineticEnergy();
182 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy;
183 //
184 // check kinematical limit: gamma energy(Eg) must be at least 2 e- rest mass
185 if (eps0 > 0.5) { return; }
186 //
187 // select target element of the material (probs. are based on partial x-secs)
188 const G4Element* anElement = SelectTargetAtom(couple, fTheGamma, gammaEnergy,
189 aDynamicGamma->GetLogKineticEnergy());
190
191 //
192 // get the random engine
193 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
194 //
195 // 'eps' is the total energy transferred to one of the e-/e+ pair in initial
196 // gamma energy units Eg. Since the corresponding DCS is symmetric on eps=0.5,
197 // the kinematical limits for eps0=mc^2/Eg <= eps <= 0.5
198 // 1. 'eps' is sampled uniformly on the [eps0, 0.5] inteval if Eg<Egsmall
199 // 2. otherwise, on the [eps_min, 0.5] interval according to the DCS (case 2.)
200 G4double eps;
201 // case 1.
202 static const G4double Egsmall = 2.*CLHEP::MeV;
203 if (gammaEnergy < Egsmall) {
204 eps = eps0 + (0.5-eps0)*rndmEngine->flat();
205 } else {
206 // case 2.
207 // get the Coulomb factor for the target element (Z) and gamma energy (Eg)
208 // F(Z) = 8*ln(Z)/3 if Eg <= 50 [MeV] => no Coulomb correction
209 // F(Z) = 8*ln(Z)/3 + 8*fc(Z) if Eg > 50 [MeV] => fc(Z) is the Coulomb cor.
210 //
211 // The screening variable 'delta(eps)' = 136*Z^{-1/3}*eps0/[eps(1-eps)]
212 // Due to the Coulomb correction, the DCS can go below zero even at
213 // kinematicaly allowed eps > eps0 values. In order to exclude this eps
214 // range with negative DCS, the minimum eps value will be set to eps_min =
215 // max[eps0, epsp] with epsp is the solution of SF(delta(epsp)) - F(Z)/2 = 0
216 // with SF being the screening function (SF1=SF2 at high value of delta).
217 // The solution is epsp = 0.5 - 0.5*sqrt[ 1 - 4*136*Z^{-1/3}eps0/deltap]
218 // with deltap = Exp[(42.038-F(Z))/8.29]-0.958. So the limits are:
219 // - when eps=eps_max = 0.5 => delta_min = 136*Z^{-1/3}*eps0/4
220 // - epsp = 0.5 - 0.5*sqrt[ 1 - delta_min/deltap]
221 // - and eps_min = max[eps0, epsp]
222 static const G4double midEnergy = 50.*CLHEP::MeV;
223 const G4int iZet = std::min(gMaxZet, anElement->GetZasInt());
224 const G4double deltaFactor = 136.*eps0/anElement->GetIonisation()->GetZ3();
225 G4double deltaMax = gElementData[iZet]->fDeltaMaxLow;
226 G4double FZ = 8.*anElement->GetIonisation()->GetlogZ3();
227 if (gammaEnergy > midEnergy) {
228 FZ += 8.*(anElement->GetfCoulomb());
229 deltaMax = gElementData[iZet]->fDeltaMaxHigh;
230 }
231 const G4double deltaMin = 4.*deltaFactor;
232 //
233 // compute the limits of eps
234 const G4double epsp = 0.5 - 0.5*std::sqrt(1. - deltaMin/deltaMax) ;
235 const G4double epsMin = std::max(eps0,epsp);
236 const G4double epsRange = 0.5 - epsMin;
237 //
238 // sample the energy rate (eps) of the created electron (or positron)
240 ScreenFunction12(deltaMin, F10, F20);
241 F10 -= FZ;
242 F20 -= FZ;
243 const G4double NormF1 = std::max(F10 * epsRange * epsRange, 0.);
244 const G4double NormF2 = std::max(1.5 * F20 , 0.);
245 const G4double NormCond = NormF1/(NormF1 + NormF2);
246 // we will need 3 uniform random number for each trial of sampling
247 G4double rndmv[3];
248 G4double greject = 0.;
249 do {
250 rndmEngine->flatArray(3, rndmv);
251 if (NormCond > rndmv[0]) {
252 eps = 0.5 - epsRange * fG4Calc->A13(rndmv[1]);
253 const G4double delta = deltaFactor/(eps*(1.-eps));
254 greject = (ScreenFunction1(delta)-FZ)/F10;
255 } else {
256 eps = epsMin + epsRange*rndmv[1];
257 const G4double delta = deltaFactor/(eps*(1.-eps));
258 greject = (ScreenFunction2(delta)-FZ)/F20;
259 }
260 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
261 } while (greject < rndmv[2]);
262 } // end of eps sampling
263 //
264 // select charges randomly
265 G4double eTotEnergy, pTotEnergy;
266 if (rndmEngine->flat() > 0.5) {
267 eTotEnergy = (1.-eps)*gammaEnergy;
268 pTotEnergy = eps*gammaEnergy;
269 } else {
270 pTotEnergy = (1.-eps)*gammaEnergy;
271 eTotEnergy = eps*gammaEnergy;
272 }
273 //
274 // sample pair kinematics
275 const G4double eKinEnergy = std::max(0.,eTotEnergy - CLHEP::electron_mass_c2);
276 const G4double pKinEnergy = std::max(0.,pTotEnergy - CLHEP::electron_mass_c2);
277 //
278 G4ThreeVector eDirection, pDirection;
279 //
281 eKinEnergy, pKinEnergy,
282 eDirection, pDirection);
283 // create G4DynamicParticle object for the particle1
284 G4DynamicParticle* aParticle1= new G4DynamicParticle(
285 fTheElectron,eDirection,eKinEnergy);
286 // create G4DynamicParticle object for the particle2
287 G4DynamicParticle* aParticle2= new G4DynamicParticle(
288 fThePositron,pDirection,pKinEnergy);
289 // Fill output vector
290 fvect->push_back(aParticle1);
291 fvect->push_back(aParticle2);
292 // kill incident photon
295}
#define F10
#define F20
@ fStopAndKill
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
void ScreenFunction12(const G4double delta, G4double &f1, G4double &f2)
G4double ScreenFunction1(const G4double delta)
G4double ScreenFunction2(const G4double delta)
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4double GetZ3() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double A13(G4double A) const
Definition: G4Pow.cc:120
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:611
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:587
void ProposeTrackStatus(G4TrackStatus status)

Referenced by G4PolarizedGammaConversionModel::SampleSecondaries().

◆ ScreenFunction1()

G4double G4BetheHeitlerModel::ScreenFunction1 ( const G4double  delta)
inlineprotected

Definition at line 142 of file G4BetheHeitlerModel.hh.

143{
144 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
145 : 42.184 - delta*(7.444 - 1.623*delta);
146}

Referenced by SampleSecondaries().

◆ ScreenFunction12()

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

Definition at line 158 of file G4BetheHeitlerModel.hh.

160{
161 if (delta > 1.4) {
162 f1 = 42.038 - 8.29*G4Log(delta + 0.958);
163 f2 = f1;
164 } else {
165 f1 = 42.184 - delta*(7.444 - 1.623*delta);
166 f2 = 41.326 - delta*(5.848 - 0.902*delta);
167 }
168}

Referenced by SampleSecondaries().

◆ ScreenFunction2()

G4double G4BetheHeitlerModel::ScreenFunction2 ( const G4double  delta)
inlineprotected

Definition at line 150 of file G4BetheHeitlerModel.hh.

151{
152 return (delta > 1.4) ? 42.038 - 8.29*G4Log(delta + 0.958)
153 : 41.326 - delta*(5.848 - 0.902*delta);
154}

Referenced by SampleSecondaries().

Member Data Documentation

◆ fG4Calc

G4Pow* G4BetheHeitlerModel::fG4Calc
protected

Definition at line 114 of file G4BetheHeitlerModel.hh.

Referenced by SampleSecondaries().

◆ fParticleChange

G4ParticleChangeForGamma* G4BetheHeitlerModel::fParticleChange
protected

Definition at line 118 of file G4BetheHeitlerModel.hh.

Referenced by Initialise(), and SampleSecondaries().

◆ fTheElectron

G4ParticleDefinition* G4BetheHeitlerModel::fTheElectron
protected

Definition at line 116 of file G4BetheHeitlerModel.hh.

Referenced by SampleSecondaries().

◆ fTheGamma

G4ParticleDefinition* G4BetheHeitlerModel::fTheGamma
protected

Definition at line 115 of file G4BetheHeitlerModel.hh.

Referenced by SampleSecondaries().

◆ fThePositron

G4ParticleDefinition* G4BetheHeitlerModel::fThePositron
protected

Definition at line 117 of file G4BetheHeitlerModel.hh.

Referenced by SampleSecondaries().

◆ gElementData

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

◆ gMaxZet

const G4int G4BetheHeitlerModel::gMaxZet = 120
staticprotected

Definition at line 112 of file G4BetheHeitlerModel.hh.

Referenced by InitialiseElementData(), and SampleSecondaries().


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