Geant4 9.6.0
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=0, const G4String &nam="BetheHeitlerLPM")
 
virtual ~G4PairProductionRelModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double)
 
void SetCurrentElement (G4double)
 
void SetLPMconstant (G4double val)
 
G4double LPMconstant () const
 
void SetLPMflag (G4bool)
 
G4bool LPMflag () const
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

G4double Phi1 (G4double delta) const
 
G4double Phi2 (G4double delta) const
 
G4double DeltaMax () const
 
G4double DeltaMin (G4double) const
 
void CalcLPMFunctions (G4double k, G4double eplus)
 
G4double ScreenFunction1 (G4double ScreenVariable)
 
G4double ScreenFunction2 (G4double ScreenVariable)
 
G4double ComputeXSectionPerAtom (G4double totalEnergy, G4double Z)
 
G4double ComputeDXSectionPerAtom (G4double eplusEnergy, G4double totalEnergy, G4double Z)
 
G4double ComputeRelDXSectionPerAtom (G4double eplusEnergy, G4double totalEnergy, G4double Z)
 
G4PairProductionRelModeloperator= (const G4PairProductionRelModel &right)
 
 G4PairProductionRelModel (const G4PairProductionRelModel &)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4NistManagernist
 
G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleDefinitionthePositron
 
G4ParticleChangeForGammafParticleChange
 
G4double fLPMconstant
 
G4bool fLPMflag
 
G4double z13
 
G4double z23
 
G4double lnZ
 
G4double Fel
 
G4double Finel
 
G4double fCoulomb
 
G4double currentZ
 
G4double lpmEnergy
 
G4double xiLPM
 
G4double phiLPM
 
G4double gLPM
 
G4bool use_completescreening
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Static Protected Attributes

static const G4double xgi [8]
 
static const G4double wgi [8]
 
static const G4double Fel_light [5] = {0., 5.31 , 4.79 , 4.74 , 4.71}
 
static const G4double Finel_light [5] = {0., 6.144 , 5.621 , 5.805 , 5.924}
 
static const G4double facFel = log(184.15)
 
static const G4double facFinel = log(1194.)
 
static const G4double preS1 = 1./(184.15*184.15)
 
static const G4double logTwo = log(2.)
 

Detailed Description

Definition at line 61 of file G4PairProductionRelModel.hh.

Constructor & Destructor Documentation

◆ G4PairProductionRelModel() [1/2]

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

Definition at line 86 of file G4PairProductionRelModel.cc.

88 : G4VEmModel(nam),
89 fLPMconstant(fine_structure_const*electron_mass_c2*electron_mass_c2/(4.*pi*hbarc)*0.5),
90 fLPMflag(true),
91 lpmEnergy(0.),
93{
98
100
101 currentZ = z13 = z23 = lnZ = Fel = Finel = fCoulomb = phiLPM = gLPM = xiLPM = 0;
102
103}
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4NistManager * Instance()
G4ParticleDefinition * thePositron
G4ParticleDefinition * theElectron
G4ParticleChangeForGamma * fParticleChange
static G4Positron * Positron()
Definition: G4Positron.cc:94

◆ ~G4PairProductionRelModel()

G4PairProductionRelModel::~G4PairProductionRelModel ( )
virtual

Definition at line 107 of file G4PairProductionRelModel.cc.

108{}

◆ G4PairProductionRelModel() [2/2]

G4PairProductionRelModel::G4PairProductionRelModel ( const G4PairProductionRelModel )
protected

Member Function Documentation

◆ CalcLPMFunctions()

void G4PairProductionRelModel::CalcLPMFunctions ( G4double  k,
G4double  eplus 
)
protected

Definition at line 222 of file G4PairProductionRelModel.cc.

223{
224 // *** calculate lpm variable s & sprime ***
225 // Klein eqs. (78) & (79)
226 G4double sprime = sqrt(0.125*k*lpmEnergy/(eplusEnergy*(k-eplusEnergy)));
227
228 G4double s1 = preS1*z23;
229 G4double logS1 = 2./3.*lnZ-2.*facFel;
230 G4double logTS1 = logTwo+logS1;
231
232 xiLPM = 2.;
233
234 if (sprime>1)
235 xiLPM = 1.;
236 else if (sprime>sqrt(2.)*s1) {
237 G4double h = log(sprime)/logTS1;
238 xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
239 }
240
241 G4double s0 = sprime/sqrt(xiLPM);
242 // G4cout<<"k="<<k<<" y="<<eplusEnergy/k<<G4endl;
243 // G4cout<<"s0="<<s0<<G4endl;
244
245 // *** calculate supression functions phi and G ***
246 // Klein eqs. (77)
247 G4double s2=s0*s0;
248 G4double s3=s0*s2;
249 G4double s4=s2*s2;
250
251 if (s0<0.1) {
252 // high suppression limit
253 phiLPM = 6.*s0 - 18.84955592153876*s2 + 39.47841760435743*s3
254 - 57.69873135166053*s4;
255 gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
256 }
257 else if (s0<1.9516) {
258 // intermediate suppression
259 // using eq.77 approxim. valid s0<2.
260 phiLPM = 1.-exp(-6.*s0*(1.+(3.-pi)*s0)
261 +s3/(0.623+0.795*s0+0.658*s2));
262 if (s0<0.415827397755) {
263 // using eq.77 approxim. valid 0.07<s<2
264 G4double psiLPM = 1-exp(-4*s0-8*s2/(1+3.936*s0+4.97*s2-0.05*s3+7.50*s4));
265 gLPM = 3*psiLPM-2*phiLPM;
266 }
267 else {
268 // using alternative parametrisiation
269 G4double pre = -0.16072300849123999 + s0*3.7550300067531581 + s2*-1.7981383069010097
270 + s3*0.67282686077812381 + s4*-0.1207722909879257;
271 gLPM = tanh(pre);
272 }
273 }
274 else {
275 // low suppression limit valid s>2.
276 phiLPM = 1. - 0.0119048/s4;
277 gLPM = 1. - 0.0230655/s4;
278 }
279
280 // *** make sure suppression is smaller than 1 ***
281 // *** caused by Migdal approximation in xi ***
282 if (xiLPM*phiLPM>1. || s0>0.57) xiLPM=1./phiLPM;
283}
double G4double
Definition: G4Types.hh:64
T sqr(const T &x)
Definition: templates.hh:145

Referenced by ComputeRelDXSectionPerAtom(), and SampleSecondaries().

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 289 of file G4PairProductionRelModel.cc.

292{
293 // static const G4double gammaEnergyLimit = 1.5*MeV;
294 G4double crossSection = 0.0 ;
295 if ( Z < 0.9 ) return crossSection;
296 if ( gammaEnergy <= 2.0*electron_mass_c2 ) return crossSection;
297
299 // choose calculator according to parameters and switches
300 // in the moment only one calculator:
301 crossSection=ComputeXSectionPerAtom(gammaEnergy,Z);
302
303 G4double xi = Finel/(Fel - fCoulomb); // inelastic contribution
304 crossSection*=4.*Z*(Z+xi)*fine_structure_const*classic_electr_radius*classic_electr_radius;
305
306 return crossSection;
307}
G4double ComputeXSectionPerAtom(G4double totalEnergy, G4double Z)

◆ ComputeDXSectionPerAtom()

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

Definition at line 165 of file G4PairProductionRelModel.cc.

168{
169 // most simple case - complete screening:
170
171 // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
172 // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
173 // y = E+/k
174 G4double yp=eplusEnergy/totalEnergy;
175 G4double ym=1.-yp;
176
177 G4double cross = 0.;
179 cross = (yp*yp + ym*ym + 2./3.*ym*yp)*(Fel - fCoulomb) + 1./9.*yp*ym;
180 else {
181 G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
182 cross = (yp*yp + ym*ym)*(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
183 + 2./3.*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
184 }
185 return cross/totalEnergy;
186
187}
G4double Phi2(G4double delta) const
G4double DeltaMin(G4double) const
G4double Phi1(G4double delta) const

Referenced by ComputeXSectionPerAtom().

◆ ComputeRelDXSectionPerAtom()

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

Definition at line 191 of file G4PairProductionRelModel.cc.

194{
195 // most simple case - complete screening:
196
197 // dsig/dE+ = 4 * alpha * Z**2 * r0**2 / k
198 // * [ (y**2 + (1-y**2) + 2/3*y*(1-y) ) * ( log (183 * Z**-1/3) + 1/9 * y*(1-y) ]
199 // y = E+/k
200 G4double yp=eplusEnergy/totalEnergy;
201 G4double ym=1.-yp;
202
203 CalcLPMFunctions(totalEnergy,eplusEnergy); // gamma
204
205 G4double cross = 0.;
207 cross = xiLPM*(2./3.*phiLPM*(yp*yp + ym*ym) + gLPM)*(Fel - fCoulomb);
208 else {
209 G4double delta = 0.25*DeltaMin(totalEnergy)/(yp*ym);
210 cross = (1./3.*gLPM + 2./3.*phiLPM)*(yp*yp + ym*ym)
211 *(0.25*Phi1(delta) - lnZ/3. - fCoulomb)
212 + 2./3.*gLPM*ym*yp*(0.25*Phi2(delta) - lnZ/3. - fCoulomb);
213 cross *= xiLPM;
214 }
215 return cross/totalEnergy;
216
217}
void CalcLPMFunctions(G4double k, G4double eplus)

Referenced by ComputeXSectionPerAtom().

◆ ComputeXSectionPerAtom()

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

Definition at line 121 of file G4PairProductionRelModel.cc.

122{
123 G4double cross = 0.0;
124
125 // number of intervals and integration step
126 G4double vcut = electron_mass_c2/totalEnergy ;
127
128 // limits by the screening variable
129 G4double dmax = DeltaMax();
130 G4double dmin = min(DeltaMin(totalEnergy),dmax);
131 G4double vcut1 = 0.5 - 0.5*sqrt(1. - dmin/dmax) ;
132 vcut = max(vcut, vcut1);
133
134
135 G4double vmax = 0.5;
136 G4int n = 1; // needs optimisation
137
138 G4double delta = (vmax - vcut)*totalEnergy/G4double(n);
139
140 G4double e0 = vcut*totalEnergy;
141 G4double xs;
142
143 // simple integration
144 for(G4int l=0; l<n; l++,e0 += delta) {
145 for(G4int i=0; i<8; i++) {
146
147 G4double eg = (e0 + xgi[i]*delta);
148 if (fLPMflag && totalEnergy>100.*GeV)
149 xs = ComputeRelDXSectionPerAtom(eg,totalEnergy,Z);
150 else
151 xs = ComputeDXSectionPerAtom(eg,totalEnergy,Z);
152 cross += wgi[i]*xs;
153
154 }
155 }
156
157 cross *= delta*2.;
158
159 return cross;
160}
int G4int
Definition: G4Types.hh:66
static const G4double wgi[8]
G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)
static const G4double xgi[8]
G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z)

Referenced by ComputeCrossSectionPerAtom().

◆ DeltaMax()

G4double G4PairProductionRelModel::DeltaMax ( ) const
inlineprotected

Definition at line 274 of file G4PairProductionRelModel.hh.

275{
276 // k > 50 MeV
277 G4double FZ = 8.*(lnZ/3. + fCoulomb);
278 return std::exp( (42.24-FZ)/8.368 ) + 0.952;
279}

Referenced by ComputeXSectionPerAtom().

◆ DeltaMin()

G4double G4PairProductionRelModel::DeltaMin ( G4double  k) const
inlineprotected

Definition at line 281 of file G4PairProductionRelModel.hh.

282{
283 return 4.*136./z13*(CLHEP::electron_mass_c2/k);
284}

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 112 of file G4PairProductionRelModel.cc.

114{
117}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123

◆ LPMconstant()

G4double G4PairProductionRelModel::LPMconstant ( ) const
inline

Definition at line 166 of file G4PairProductionRelModel.hh.

167{
168 return fLPMconstant;
169}

◆ LPMflag()

G4bool G4PairProductionRelModel::LPMflag ( ) const
inline

Definition at line 182 of file G4PairProductionRelModel.hh.

183{
184 return fLPMflag;
185}

◆ operator=()

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

◆ Phi1()

G4double G4PairProductionRelModel::Phi1 ( G4double  delta) const
inlineprotected

Definition at line 213 of file G4PairProductionRelModel.hh.

214{
215 G4double screenVal;
216
217 if (delta > 1.)
218 screenVal = 21.12 - 4.184*std::log(delta+0.952);
219 else
220 screenVal = 20.868 - delta*(3.242 - 0.625*delta);
221
222 return screenVal;
223}

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

◆ Phi2()

G4double G4PairProductionRelModel::Phi2 ( G4double  delta) const
inlineprotected

Definition at line 227 of file G4PairProductionRelModel.hh.

228{
229 G4double screenVal;
230
231 if (delta > 1.)
232 screenVal = 21.12 - 4.184*std::log(delta+0.952);
233 else
234 screenVal = 20.209 - delta*(1.930 + 0.086*delta);
235
236 return screenVal;
237}

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 312 of file G4PairProductionRelModel.cc.

329{
330 const G4Material* aMaterial = couple->GetMaterial();
331
332 G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
333 G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
334
335 G4double epsil ;
336 G4double epsil0 = electron_mass_c2/GammaEnergy ;
337 if(epsil0 > 1.0) { return; }
338
339 // do it fast if GammaEnergy < 2. MeV
340 static const G4double Egsmall=2.*MeV;
341
342 // select randomly one element constituing the material
343 const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
344
345 if (GammaEnergy < Egsmall) {
346
347 epsil = epsil0 + (0.5-epsil0)*G4UniformRand();
348
349 } else {
350 // now comes the case with GammaEnergy >= 2. MeV
351
352 // Extract Coulomb factor for this Element
353 G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
354 if (GammaEnergy > 50.*MeV) FZ += 8.*(anElement->GetfCoulomb());
355
356 // limits of the screening variable
357 G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
358 G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
359 G4double screenmin = min(4.*screenfac,screenmax);
360
361 // limits of the energy sampling
362 G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
363 G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
364
365 //
366 // sample the energy rate of the created electron (or positron)
367 //
368 //G4double epsil, screenvar, greject ;
369 G4double screenvar, greject ;
370
371 G4double F10 = ScreenFunction1(screenmin) - FZ;
372 G4double F20 = ScreenFunction2(screenmin) - FZ;
373 G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
374 G4double NormF2 = max(1.5*F20,0.);
375
376 do {
377 if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) {
378 epsil = 0.5 - epsilrange*pow(G4UniformRand(), 0.333333);
379 screenvar = screenfac/(epsil*(1-epsil));
380 if (fLPMflag && GammaEnergy>100.*GeV) {
381 CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
382 greject = xiLPM*((gLPM+2.*phiLPM)*Phi1(screenvar) - gLPM*Phi2(screenvar) - phiLPM*FZ)/F10;
383 }
384 else {
385 greject = (ScreenFunction1(screenvar) - FZ)/F10;
386 }
387
388 } else {
389 epsil = epsilmin + epsilrange*G4UniformRand();
390 screenvar = screenfac/(epsil*(1-epsil));
391 if (fLPMflag && GammaEnergy>100.*GeV) {
392 CalcLPMFunctions(GammaEnergy,GammaEnergy*epsil);
393 greject = xiLPM*((0.5*gLPM+phiLPM)*Phi1(screenvar) + 0.5*gLPM*Phi2(screenvar) - 0.5*(gLPM+phiLPM)*FZ)/F20;
394 }
395 else {
396 greject = (ScreenFunction2(screenvar) - FZ)/F20;
397 }
398 }
399
400 } while( greject < G4UniformRand() );
401
402 } // end of epsil sampling
403
404 //
405 // fixe charges randomly
406 //
407
408 G4double ElectTotEnergy, PositTotEnergy;
409 if (G4UniformRand() > 0.5) {
410
411 ElectTotEnergy = (1.-epsil)*GammaEnergy;
412 PositTotEnergy = epsil*GammaEnergy;
413
414 } else {
415
416 PositTotEnergy = (1.-epsil)*GammaEnergy;
417 ElectTotEnergy = epsil*GammaEnergy;
418 }
419
420 //
421 // scattered electron (positron) angles. ( Z - axis along the parent photon)
422 //
423 // universal distribution suggested by L. Urban
424 // (Geant3 manual (1993) Phys211),
425 // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
426
427 G4double u;
428 const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
429
430 if (9./(9.+d) >G4UniformRand()) u= - log(G4UniformRand()*G4UniformRand())/a1;
431 else u= - log(G4UniformRand()*G4UniformRand())/a2;
432
433 G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
434 G4double TetPo = u*electron_mass_c2/PositTotEnergy;
435 G4double Phi = twopi * G4UniformRand();
436 G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
437 G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
438
439 //
440 // kinematic of the created pair
441 //
442 // the electron and positron are assumed to have a symetric
443 // angular distribution with respect to the Z axis along the parent photon.
444
445 G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
446
447 G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
448 ElectDirection.rotateUz(GammaDirection);
449
450 // create G4DynamicParticle object for the particle1
451 G4DynamicParticle* aParticle1= new G4DynamicParticle(
452 theElectron,ElectDirection,ElectKineEnergy);
453
454 // the e+ is always created (even with Ekine=0) for further annihilation.
455
456 G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
457
458 G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
459 PositDirection.rotateUz(GammaDirection);
460
461 // create G4DynamicParticle object for the particle2
462 G4DynamicParticle* aParticle2= new G4DynamicParticle(
463 thePositron,PositDirection,PositKineEnergy);
464
465 // Fill output vector
466 fvect->push_back(aParticle1);
467 fvect->push_back(aParticle2);
468
469 // kill incident photon
472}
#define F10
#define F20
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetfCoulomb() const
Definition: G4Element.hh:201
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:209
G4double GetlogZ3() const
G4double GetZ3() const
const G4Material * GetMaterial() const
G4double ScreenFunction2(G4double ScreenVariable)
G4double ScreenFunction1(G4double ScreenVariable)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void ProposeTrackStatus(G4TrackStatus status)

◆ ScreenFunction1()

G4double G4PairProductionRelModel::ScreenFunction1 ( G4double  ScreenVariable)
inlineprotected

Definition at line 239 of file G4PairProductionRelModel.hh.

243{
244 G4double screenVal;
245
246 if (ScreenVariable > 1.)
247 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
248 else
249 screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable);
250
251 return screenVal;
252}

Referenced by SampleSecondaries().

◆ ScreenFunction2()

G4double G4PairProductionRelModel::ScreenFunction2 ( G4double  ScreenVariable)
inlineprotected

Definition at line 256 of file G4PairProductionRelModel.hh.

260{
261 G4double screenVal;
262
263 if (ScreenVariable > 1.)
264 screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952);
265 else
266 screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable);
267
268 return screenVal;
269}

Referenced by SampleSecondaries().

◆ SetCurrentElement()

void G4PairProductionRelModel::SetCurrentElement ( G4double  Z)
inline

Definition at line 189 of file G4PairProductionRelModel.hh.

190{
191 if(Z != currentZ) {
192 currentZ = Z;
193
194 G4int iz = G4int(Z);
195 z13 = nist->GetZ13(iz);
196 z23 = z13*z13;
197 lnZ = nist->GetLOGZ(iz);
198
199 if (iz <= 4) {
200 Fel = Fel_light[iz];
201 Finel = Finel_light[iz] ;
202 }
203 else {
204 Fel = facFel - lnZ/3. ;
205 Finel = facFinel - 2.*lnZ/3. ;
206 }
208 }
209}
G4double GetZ13(G4double Z)
G4double GetLOGZ(G4int Z)
static const G4double Finel_light[5]
static const G4double Fel_light[5]
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:391

Referenced by ComputeCrossSectionPerAtom().

◆ SetLPMconstant()

void G4PairProductionRelModel::SetLPMconstant ( G4double  val)
inline

Definition at line 158 of file G4PairProductionRelModel.hh.

159{
160 fLPMconstant = val;
161}

◆ SetLPMflag()

void G4PairProductionRelModel::SetLPMflag ( G4bool  val)
inline

Definition at line 174 of file G4PairProductionRelModel.hh.

175{
176 fLPMflag = val;
177}

◆ SetupForMaterial()

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

Reimplemented from G4VEmModel.

Definition at line 477 of file G4PairProductionRelModel.cc.

479{
481 // G4cout<<" lpmEnergy="<<lpmEnergy<<G4endl;
482}
G4double GetRadlen() const
Definition: G4Material.hh:219

Member Data Documentation

◆ currentZ

G4double G4PairProductionRelModel::currentZ
protected

Definition at line 135 of file G4PairProductionRelModel.hh.

Referenced by G4PairProductionRelModel(), and SetCurrentElement().

◆ facFel

const G4double G4PairProductionRelModel::facFel = log(184.15)
staticprotected

Definition at line 147 of file G4PairProductionRelModel.hh.

Referenced by CalcLPMFunctions(), and SetCurrentElement().

◆ facFinel

const G4double G4PairProductionRelModel::facFinel = log(1194.)
staticprotected

Definition at line 148 of file G4PairProductionRelModel.hh.

Referenced by SetCurrentElement().

◆ fCoulomb

◆ Fel

◆ Fel_light

const G4double G4PairProductionRelModel::Fel_light = {0., 5.31 , 4.79 , 4.74 , 4.71}
staticprotected

Definition at line 145 of file G4PairProductionRelModel.hh.

Referenced by SetCurrentElement().

◆ Finel

G4double G4PairProductionRelModel::Finel
protected

◆ Finel_light

const G4double G4PairProductionRelModel::Finel_light = {0., 6.144 , 5.621 , 5.805 , 5.924}
staticprotected

Definition at line 146 of file G4PairProductionRelModel.hh.

Referenced by SetCurrentElement().

◆ fLPMconstant

G4double G4PairProductionRelModel::fLPMconstant
protected

Definition at line 129 of file G4PairProductionRelModel.hh.

Referenced by LPMconstant(), SetLPMconstant(), and SetupForMaterial().

◆ fLPMflag

G4bool G4PairProductionRelModel::fLPMflag
protected

◆ fParticleChange

G4ParticleChangeForGamma* G4PairProductionRelModel::fParticleChange
protected

◆ gLPM

G4double G4PairProductionRelModel::gLPM
protected

◆ lnZ

◆ logTwo

const G4double G4PairProductionRelModel::logTwo = log(2.)
staticprotected

Definition at line 150 of file G4PairProductionRelModel.hh.

Referenced by CalcLPMFunctions().

◆ lpmEnergy

G4double G4PairProductionRelModel::lpmEnergy
protected

Definition at line 138 of file G4PairProductionRelModel.hh.

Referenced by CalcLPMFunctions(), and SetupForMaterial().

◆ nist

G4NistManager* G4PairProductionRelModel::nist
protected

Definition at line 122 of file G4PairProductionRelModel.hh.

Referenced by G4PairProductionRelModel(), and SetCurrentElement().

◆ phiLPM

G4double G4PairProductionRelModel::phiLPM
protected

◆ preS1

const G4double G4PairProductionRelModel::preS1 = 1./(184.15*184.15)
staticprotected

Definition at line 150 of file G4PairProductionRelModel.hh.

Referenced by CalcLPMFunctions().

◆ theElectron

G4ParticleDefinition* G4PairProductionRelModel::theElectron
protected

Definition at line 125 of file G4PairProductionRelModel.hh.

Referenced by G4PairProductionRelModel(), and SampleSecondaries().

◆ theGamma

G4ParticleDefinition* G4PairProductionRelModel::theGamma
protected

Definition at line 124 of file G4PairProductionRelModel.hh.

Referenced by G4PairProductionRelModel(), and SampleSecondaries().

◆ thePositron

G4ParticleDefinition* G4PairProductionRelModel::thePositron
protected

Definition at line 126 of file G4PairProductionRelModel.hh.

Referenced by G4PairProductionRelModel(), and SampleSecondaries().

◆ use_completescreening

G4bool G4PairProductionRelModel::use_completescreening
protected

◆ wgi

const G4double G4PairProductionRelModel::wgi
staticprotected
Initial value:
={ 0.0506, 0.1112, 0.1569, 0.1813,
0.1813, 0.1569, 0.1112, 0.0506 }

Definition at line 144 of file G4PairProductionRelModel.hh.

Referenced by ComputeXSectionPerAtom().

◆ xgi

const G4double G4PairProductionRelModel::xgi
staticprotected
Initial value:
={ 0.0199, 0.1017, 0.2372, 0.4083,
0.5917, 0.7628, 0.8983, 0.9801 }

Definition at line 144 of file G4PairProductionRelModel.hh.

Referenced by ComputeXSectionPerAtom().

◆ xiLPM

G4double G4PairProductionRelModel::xiLPM
protected

◆ z13

G4double G4PairProductionRelModel::z13
protected

◆ z23

G4double G4PairProductionRelModel::z23
protected

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