Geant4  10.02.p03
G4eBremsstrahlungRelModel Class Reference

#include <G4eBremsstrahlungRelModel.hh>

Inheritance diagram for G4eBremsstrahlungRelModel:
Collaboration diagram for G4eBremsstrahlungRelModel:

Public Member Functions

 G4eBremsstrahlungRelModel (const G4ParticleDefinition *p=0, const G4String &nam="eBremLPM")
 
virtual ~G4eBremsstrahlungRelModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut)
 
void SetLPMconstant (G4double val)
 
G4double LPMconstant () const
 
void SetLowestKinEnergy (G4double)
 
G4double LowestKinEnergy () const
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
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 *> *)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Member Functions

virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
void SetCurrentElement (const G4double)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4NistManagernist
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLoss * fParticleChange
 
G4double bremFactor
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double currentZ
 
G4double densityFactor
 
G4double densityCorr
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

void InitialiseConstants ()
 
void CalcLPMFunctions (G4double gammaEnergy)
 
G4double ComputeBremLoss (G4double cutEnergy)
 
G4double ComputeXSectionPerAtom (G4double cutEnergy)
 
G4double ComputeRelDXSectionPerAtom (G4double gammaEnergy)
 
void SetParticle (const G4ParticleDefinition *p)
 
G4double Phi1 (G4double, G4double)
 
G4double Phi1M2 (G4double, G4double)
 
G4double Psi1 (G4double, G4double)
 
G4double Psi1M2 (G4double, G4double)
 
G4eBremsstrahlungRelModeloperator= (const G4eBremsstrahlungRelModel &right)
 
 G4eBremsstrahlungRelModel (const G4eBremsstrahlungRelModel &)
 

Private Attributes

G4double lowestKinEnergy
 
G4double fMigdalConstant
 
G4double fLPMconstant
 
G4double energyThresholdLPM
 
G4double facFel
 
G4double facFinel
 
G4double preS1
 
G4double logTwo
 
G4double z13
 
G4double z23
 
G4double lnZ
 
G4double Fel
 
G4double Finel
 
G4double fCoulomb
 
G4double fMax
 
G4double lpmEnergy
 
G4PhysicsVectorfXiLPM
 
G4PhysicsVectorfPhiLPM
 
G4PhysicsVectorfGLPM
 
G4double xiLPM
 
G4double phiLPM
 
G4double gLPM
 
G4double klpm
 
G4double kp
 
G4bool use_completescreening
 

Static Private 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}
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 62 of file G4eBremsstrahlungRelModel.hh.

Constructor & Destructor Documentation

◆ G4eBremsstrahlungRelModel() [1/2]

G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eBremLPM" 
)

Definition at line 87 of file G4eBremsstrahlungRelModel.cc.

89  : G4VEmModel(nam),
90  particle(0),
92  isElectron(true),
95  fXiLPM(0), fPhiLPM(0), fGLPM(0),
97 {
98  fParticleChange = 0;
100 
101  lowestKinEnergy = 1.0*MeV;
103 
105 
106  SetLPMFlag(true);
107  //SetAngularDistribution(new G4ModifiedTsai());
109 
112  = xiLPM = phiLPM = gLPM = klpm = kp = 0.0;
113 
114  energyThresholdLPM = 1.e39;
115 
117  if(p) { SetParticle(p); }
118 }
static const double MeV
Definition: G4SIunits.hh:211
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
const G4ParticleDefinition * particle
int fine_structure_const
Definition: hepunit.py:287
float hbarc
Definition: hepunit.py:265
float electron_Compton_length
Definition: hepunit.py:289
float electron_mass_c2
Definition: hepunit.py:274
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
int classic_electr_radius
Definition: hepunit.py:288
static const double pi
Definition: G4SIunits.hh:74
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:774
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
G4ParticleChangeForLoss * fParticleChange
void SetParticle(const G4ParticleDefinition *p)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
Here is the call graph for this function:

◆ ~G4eBremsstrahlungRelModel()

G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel ( )
virtual

Definition at line 133 of file G4eBremsstrahlungRelModel.cc.

134 {
135 }

◆ G4eBremsstrahlungRelModel() [2/2]

G4eBremsstrahlungRelModel::G4eBremsstrahlungRelModel ( const G4eBremsstrahlungRelModel )
private

Member Function Documentation

◆ CalcLPMFunctions()

void G4eBremsstrahlungRelModel::CalcLPMFunctions ( G4double  gammaEnergy)
private

Definition at line 347 of file G4eBremsstrahlungRelModel.cc.

348 {
349  // *** calculate lpm variable s & sprime ***
350  // Klein eqs. (78) & (79)
351  G4double sprime = sqrt(0.125*k*lpmEnergy/(totalEnergy*(totalEnergy-k)));
352 
353  G4double s1 = preS1*z23;
354  G4double logS1 = 2./3.*lnZ-2.*facFel;
355  G4double logTS1 = logTwo+logS1;
356 
357  xiLPM = 2.;
358 
359  if (sprime>1)
360  xiLPM = 1.;
361  else if (sprime>sqrt(2.)*s1) {
362  G4double h = G4Log(sprime)/logTS1;
363  xiLPM = 1+h-0.08*(1-h)*(1-sqr(1-h))/logTS1;
364  }
365 
366  G4double s0 = sprime/sqrt(xiLPM);
367 
368  // *** merging with density effect*** should be only necessary in region
369  // "close to" kp, e.g. k<100*kp using Ter-Mikaelian eq. (20.9)
370  G4double k2 = k*k;
371  s0 *= (1 + (densityCorr/k2) );
372 
373  // recalculate Xi using modified s above
374  // Klein eq. (75)
375  xiLPM = 1.;
376  if (s0<=s1) xiLPM = 2.;
377  else if ( (s1<s0) && (s0<=1) ) { xiLPM = 1. + G4Log(s0)/logS1; }
378 
379 
380  // *** calculate supression functions phi and G ***
381  // Klein eqs. (77)
382  G4double s2=s0*s0;
383  G4double s3=s0*s2;
384  G4double s4=s2*s2;
385 
386  if (s0<0.1) {
387  // high suppression limit
388  phiLPM = 6.*s0 - 18.84955592153876*s2 + 39.47841760435743*s3
389  - 57.69873135166053*s4;
390  gLPM = 37.69911184307752*s2 - 236.8705056261446*s3 + 807.7822389*s4;
391  }
392  else if (s0<1.9516) {
393  // intermediate suppression
394  // using eq.77 approxim. valid s<2.
395  phiLPM = 1.-G4Exp(-6.*s0*(1.+(3.-pi)*s0)
396  +s3/(0.623+0.795*s0+0.658*s2));
397  if (s0<0.415827397755) {
398  // using eq.77 approxim. valid 0.07<s<2
399  G4double psiLPM = 1-G4Exp(-4*s0-8*s2/(1+3.936*s0+4.97*s2-0.05*s3+7.50*s4));
400  gLPM = 3*psiLPM-2*phiLPM;
401  }
402  else {
403  // using alternative parametrisiation
404  G4double pre = -0.16072300849123999 + s0*3.7550300067531581 + s2*-1.7981383069010097
405  + s3*0.67282686077812381 + s4*-0.1207722909879257;
406  gLPM = tanh(pre);
407  }
408  }
409  else {
410  // low suppression limit valid s>2.
411  phiLPM = 1. - 0.0119048/s4;
412  gLPM = 1. - 0.0230655/s4;
413  }
414 
415  // *** make sure suppression is smaller than 1 ***
416  // *** caused by Migdal approximation in xi ***
417  if (xiLPM*phiLPM>1. || s0>0.57) { xiLPM=1./phiLPM; }
418 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static const double pi
Definition: G4SIunits.hh:74
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeBremLoss()

G4double G4eBremsstrahlungRelModel::ComputeBremLoss ( G4double  cutEnergy)
private

Definition at line 244 of file G4eBremsstrahlungRelModel.cc.

245 {
246  G4double loss = 0.0;
247 
248  // number of intervals and integration step
249  G4double vcut = cut/totalEnergy;
250  G4int n = (G4int)(20*vcut) + 3;
251  G4double delta = vcut/G4double(n);
252 
253  G4double e0 = 0.0;
254  G4double xs;
255 
256  // integration
257  for(G4int l=0; l<n; l++) {
258 
259  for(G4int i=0; i<8; i++) {
260 
261  G4double eg = (e0 + xgi[i]*delta)*totalEnergy;
262 
265  } else {
266  xs = ComputeDXSectionPerAtom(eg);
267  }
268  loss += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
269  }
270  e0 += delta;
271  }
272 
273  loss *= delta*totalEnergy;
274 
275  return loss;
276 }
int G4int
Definition: G4Types.hh:78
Char_t n[5]
G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy)
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeCrossSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  tkin,
G4double  Z,
G4double  ,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 280 of file G4eBremsstrahlungRelModel.cc.

286 {
287  if(!particle) { SetParticle(p); }
288  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
289  G4double cut = std::min(cutEnergy, kineticEnergy);
290  G4double tmax = std::min(maxEnergy, kineticEnergy);
291 
292  if(cut >= tmax) { return 0.0; }
293 
295 
296  G4double cross = ComputeXSectionPerAtom(cut);
297 
298  // allow partial integration
299  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
300 
301  cross *= Z*Z*bremFactor;
302 
303  return cross;
304 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
const G4ParticleDefinition * particle
Float_t Z
G4double ComputeXSectionPerAtom(G4double cutEnergy)
void SetParticle(const G4ParticleDefinition *p)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ComputeDEDXPerVolume()

G4double G4eBremsstrahlungRelModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 210 of file G4eBremsstrahlungRelModel.cc.

215 {
216  if(!particle) { SetParticle(p); }
217  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
218  G4double cut = std::min(cutEnergy, kineticEnergy);
219  if(cut == 0.0) { return 0.0; }
220 
221  SetupForMaterial(particle, material,kineticEnergy);
222 
223  const G4ElementVector* theElementVector = material->GetElementVector();
224  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
225 
226  G4double dedx = 0.0;
227 
228  // loop for elements in the material
229  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
230 
231  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
232  SetCurrentElement((*theElementVector)[i]->GetZ());
233 
234  dedx += theAtomicNumDensityVector[i]*currentZ*currentZ*ComputeBremLoss(cut);
235  }
236  dedx *= bremFactor;
237 
238 
239  return dedx;
240 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
std::vector< G4Element * > G4ElementVector
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
const G4ParticleDefinition * particle
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
void SetParticle(const G4ParticleDefinition *p)
double G4double
Definition: G4Types.hh:76
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:459
G4double ComputeBremLoss(G4double cutEnergy)
Here is the call graph for this function:

◆ ComputeDXSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
protectedvirtual

Reimplemented in G4LivermoreBremsstrahlungModel, and G4SeltzerBergerModel.

Definition at line 449 of file G4eBremsstrahlungRelModel.cc.

454 {
455 
456  if(gammaEnergy < 0.0) { return 0.0; }
457 
458  G4double y = gammaEnergy/totalEnergy;
459 
460  G4double main=0.,secondTerm=0.;
461 
463  // ** form factors complete screening case **
464  main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currentZ );
465  secondTerm = (1.-y)/12.*(1.+1./currentZ);
466  }
467  else {
468  // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5**
469  G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy);
470  G4double gg=dd/z13;
471  G4double eps=dd/z23;
472  G4double phi1=Phi1(gg,currentZ), phi1m2=Phi1M2(gg,currentZ);
473  G4double psi1=Psi1(eps,currentZ), psi1m2=Psi1M2(eps,currentZ);
474 
475  main = (3./4.*y*y - y + 1.) * ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currentZ );
476  secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currentZ);
477  }
478  G4double cross = main+secondTerm;
479  return cross;
480 }
G4double Psi1(G4double, G4double)
int main(int argc, char **argv)
Definition: genwindef.cpp:30
static const G4double eps
G4double Psi1M2(G4double, G4double)
Double_t y
G4double Phi1(G4double, G4double)
float electron_mass_c2
Definition: hepunit.py:274
G4double Phi1M2(G4double, G4double)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeRelDXSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeRelDXSectionPerAtom ( G4double  gammaEnergy)
private

Definition at line 423 of file G4eBremsstrahlungRelModel.cc.

427 {
428  if(gammaEnergy < 0.0) { return 0.0; }
429 
430  G4double y = gammaEnergy/totalEnergy;
431  G4double y2 = y*y*.25;
432  G4double yone2 = (1.-y+2.*y2);
433 
434  // ** form factors complete screening case **
435 
436  // ** calc LPM functions -- include ter-mikaelian merging with density effect **
437  // G4double xiLPM, gLPM, phiLPM; // to be made member variables !!!
438  CalcLPMFunctions(gammaEnergy);
439 
440  G4double mainLPM = xiLPM*(y2 * gLPM + yone2*phiLPM) * ( (Fel-fCoulomb) + Finel/currentZ );
441  G4double secondTerm = (1.-y)/12.*(1.+1./currentZ);
442 
443  G4double cross = mainLPM+secondTerm;
444  return cross;
445 }
Double_t y2[nxs]
Double_t y
double G4double
Definition: G4Types.hh:76
void CalcLPMFunctions(G4double gammaEnergy)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeXSectionPerAtom()

G4double G4eBremsstrahlungRelModel::ComputeXSectionPerAtom ( G4double  cutEnergy)
private

Definition at line 309 of file G4eBremsstrahlungRelModel.cc.

310 {
311  G4double cross = 0.0;
312 
313  // number of intervals and integration step
314  G4double vcut = G4Log(cut/totalEnergy);
316  G4int n = (G4int)(0.45*(vmax - vcut)) + 4;
317  // n=1; // integration test
318  G4double delta = (vmax - vcut)/G4double(n);
319 
320  G4double e0 = vcut;
321  G4double xs;
322 
323  // integration
324  for(G4int l=0; l<n; l++) {
325 
326  for(G4int i=0; i<8; i++) {
327 
328  G4double eg = G4Exp(e0 + xgi[i]*delta)*totalEnergy;
329 
332  } else {
333  xs = ComputeDXSectionPerAtom(eg);
334  }
335  cross += wgi[i]*xs/(1.0 + densityCorr/(eg*eg));
336  }
337  e0 += delta;
338  }
339 
340  cross *= delta;
341 
342  return cross;
343 }
int G4int
Definition: G4Types.hh:78
Char_t n[5]
G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel, G4SeltzerBergerModel, and G4ePolarizedBremsstrahlungModel.

Definition at line 174 of file G4eBremsstrahlungRelModel.cc.

176 {
177  if(p) { SetParticle(p); }
178 
179  currentZ = 0.;
180 
181  if(IsMaster() && LowEnergyLimit() < HighEnergyLimit()) {
182  InitialiseElementSelectors(p, cuts);
183  }
184 
186 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4ParticleChangeForLoss * fParticleChange
void SetParticle(const G4ParticleDefinition *p)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseConstants()

void G4eBremsstrahlungRelModel::InitialiseConstants ( )
private

Definition at line 122 of file G4eBremsstrahlungRelModel.cc.

123 {
124  facFel = G4Log(184.15);
125  facFinel = G4Log(1194.);
126 
127  preS1 = 1./(184.15*184.15);
128  logTwo = G4Log(2.);
129 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseLocal()

void G4eBremsstrahlungRelModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 190 of file G4eBremsstrahlungRelModel.cc.

192 {
193  if(LowEnergyLimit() < HighEnergyLimit()) {
195  }
196 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
Here is the call graph for this function:

◆ LowestKinEnergy()

G4double G4eBremsstrahlungRelModel::LowestKinEnergy ( ) const
inline

Definition at line 267 of file G4eBremsstrahlungRelModel.hh.

268 {
269  return lowestKinEnergy;
270 }
Here is the caller graph for this function:

◆ LPMconstant()

G4double G4eBremsstrahlungRelModel::LPMconstant ( ) const
inline

Definition at line 255 of file G4eBremsstrahlungRelModel.hh.

256 {
257  return fLPMconstant;
258 }

◆ MinPrimaryEnergy()

G4double G4eBremsstrahlungRelModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double  cut 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 201 of file G4eBremsstrahlungRelModel.cc.

204 {
205  return std::max(lowestKinEnergy, cut);
206 }

◆ operator=()

G4eBremsstrahlungRelModel& G4eBremsstrahlungRelModel::operator= ( const G4eBremsstrahlungRelModel right)
private

◆ Phi1()

G4double G4eBremsstrahlungRelModel::Phi1 ( G4double  gg,
G4double   
)
inlineprivate

Definition at line 217 of file G4eBremsstrahlungRelModel.hh.

218 {
219  // Thomas-Fermi FF from Tsai, eq.(3.38) for Z>=5
220  return 20.863 - 2.*G4Log(1. + sqr(0.55846*gg) )
221  - 4.*( 1. - 0.6*G4Exp(-0.9*gg) - 0.4*G4Exp(-1.5*gg) );
222 }
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T sqr(const T &x)
Definition: templates.hh:145
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Phi1M2()

G4double G4eBremsstrahlungRelModel::Phi1M2 ( G4double  gg,
G4double   
)
inlineprivate

Definition at line 224 of file G4eBremsstrahlungRelModel.hh.

225 {
226  // Thomas-Fermi FF from Tsai, eq. (3.39) for Z>=5
227  // return Phi1(gg,Z) -
228  return 2./(3.*(1. + 6.5*gg +6.*gg*gg) );
229 }
Here is the caller graph for this function:

◆ Psi1()

G4double G4eBremsstrahlungRelModel::Psi1 ( G4double  eps,
G4double   
)
inlineprivate

Definition at line 231 of file G4eBremsstrahlungRelModel.hh.

232 {
233  // Thomas-Fermi FF from Tsai, eq.(3.40) for Z>=5
234  return 28.340 - 2.*G4Log(1. + sqr(3.621*eps) )
235  - 4.*( 1. - 0.7*G4Exp(-8*eps) - 0.3*G4Exp(-29.*eps) );
236 }
static const G4double eps
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T sqr(const T &x)
Definition: templates.hh:145
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Psi1M2()

G4double G4eBremsstrahlungRelModel::Psi1M2 ( G4double  eps,
G4double   
)
inlineprivate

Definition at line 238 of file G4eBremsstrahlungRelModel.hh.

239 {
240  // Thomas-Fermi FF from Tsai, eq. (3.41) for Z>=5
241  return 2./(3.*(1. + 40.*eps +400.*eps*eps) );
242 }
static const G4double eps
Here is the caller graph for this function:

◆ SampleSecondaries()

void G4eBremsstrahlungRelModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel, G4SeltzerBergerModel, and G4ePolarizedBremsstrahlungModel.

Definition at line 484 of file G4eBremsstrahlungRelModel.cc.

490 {
491  G4double kineticEnergy = dp->GetKineticEnergy();
492  if(kineticEnergy < LowEnergyLimit()) { return; }
493  G4double cut = std::min(cutEnergy, kineticEnergy);
494  G4double emax = std::min(maxEnergy, kineticEnergy);
495  if(cut >= emax) { return; }
496 
497  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
498 
499  const G4Element* elm =
500  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
501  SetCurrentElement(elm->GetZ());
502 
503  kinEnergy = kineticEnergy;
504  totalEnergy = kineticEnergy + particleMass;
506 
507  //G4double fmax= fMax;
508  G4bool highe = true;
509  if(totalEnergy < energyThresholdLPM) { highe = false; }
510 
511  G4double xmin = G4Log(cut*cut + densityCorr);
512  G4double xmax = G4Log(emax*emax + densityCorr);
513  G4double gammaEnergy, f, x;
514 
515  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
516 
517  do {
518  x = G4Exp(xmin + rndmEngine->flat()*(xmax - xmin)) - densityCorr;
519  if(x < 0.0) { x = 0.0; }
520  gammaEnergy = sqrt(x);
521  if(highe) { f = ComputeRelDXSectionPerAtom(gammaEnergy); }
522  else { f = ComputeDXSectionPerAtom(gammaEnergy); }
523 
524  if ( f > fMax ) {
525  G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
526  << f << " > " << fMax
527  << " Egamma(MeV)= " << gammaEnergy
528  << " Ee(MeV)= " << kineticEnergy
529  << " " << GetName()
530  << G4endl;
531  }
532 
533  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
534  } while (f < fMax*rndmEngine->flat());
535 
536  //
537  // angles of the emitted gamma. ( Z - axis along the parent particle)
538  // use general interface
539  //
540 
541  G4ThreeVector gammaDirection =
542  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
543  G4lrint(currentZ),
544  couple->GetMaterial());
545 
546  // create G4DynamicParticle object for the Gamma
547  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
548  gammaEnergy);
549  vdp->push_back(gamma);
550 
551  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
552  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
553  - gammaEnergy*gammaDirection).unit();
554 
555  // energy of primary
556  G4double finalE = kineticEnergy - gammaEnergy;
557 
558  // stop tracking and create new secondary instead of primary
559  if(gammaEnergy > SecondaryThreshold()) {
560  fParticleChange->ProposeTrackStatus(fStopAndKill);
561  fParticleChange->SetProposedKineticEnergy(0.0);
562  G4DynamicParticle* el =
563  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
564  direction, finalE);
565  vdp->push_back(el);
566 
567  // continue tracking
568  } else {
569  fParticleChange->SetProposedMomentumDirection(direction);
570  fParticleChange->SetProposedKineticEnergy(finalE);
571  }
572 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:669
const G4String & GetName() const
Definition: G4VEmModel.hh:795
const G4Material * GetMaterial() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
virtual double flat()=0
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
const G4ParticleDefinition * particle
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
bool G4bool
Definition: G4Types.hh:79
G4double ComputeRelDXSectionPerAtom(G4double gammaEnergy)
float electron_mass_c2
Definition: hepunit.py:274
double flat()
Definition: G4AblaRandom.cc:47
static const G4double emax
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
const G4ThreeVector & GetMomentumDirection() const
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForLoss * fParticleChange
double G4double
Definition: G4Types.hh:76
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
Here is the call graph for this function:

◆ SetCurrentElement()

void G4eBremsstrahlungRelModel::SetCurrentElement ( const G4double  Z)
inlineprotected

Definition at line 189 of file G4eBremsstrahlungRelModel.hh.

190 {
191  if(Z != currentZ) {
192  currentZ = Z;
193 
194  G4int iz = G4lrint(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  G4double lnzt = lnZ/3.;
205  Fel = facFel - lnzt;
206  Finel = facFinel - 2*lnzt;
207  }
208 
210  fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
211  }
212 }
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:467
G4double GetZ13(G4double Z)
int G4int
Definition: G4Types.hh:78
static const G4double Finel_light[5]
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
G4double GetfCoulomb() const
Definition: G4Element.hh:190
int G4lrint(double ad)
Definition: templates.hh:163
static const G4double Fel_light[5]
G4double GetLOGZ(G4int Z)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetLowestKinEnergy()

void G4eBremsstrahlungRelModel::SetLowestKinEnergy ( G4double  val)
inline

Definition at line 260 of file G4eBremsstrahlungRelModel.hh.

261 {
262  lowestKinEnergy = val;
263 }
Here is the caller graph for this function:

◆ SetLPMconstant()

void G4eBremsstrahlungRelModel::SetLPMconstant ( G4double  val)
inline

Definition at line 247 of file G4eBremsstrahlungRelModel.hh.

248 {
249  fLPMconstant = val;
250 }

◆ SetParticle()

void G4eBremsstrahlungRelModel::SetParticle ( const G4ParticleDefinition p)
private

Definition at line 139 of file G4eBremsstrahlungRelModel.cc.

140 {
141  particle = p;
142  particleMass = p->GetPDGMass();
143  if(p == G4Electron::Electron()) { isElectron = true; }
144  else { isElectron = false;}
145 }
const G4ParticleDefinition * particle
static G4Electron * Electron()
Definition: G4Electron.cc:94
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetupForMaterial()

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

Reimplemented from G4VEmModel.

Definition at line 149 of file G4eBremsstrahlungRelModel.cc.

152 {
155 
156  // Threshold for LPM effect (i.e. below which LPM hidden by density effect)
157  if (LPMFlag()) {
159  } else {
160  energyThresholdLPM=1.e39; // i.e. do not use LPM effect
161  }
162  // calculate threshold for density effect
163  kinEnergy = kineticEnergy;
164  totalEnergy = kineticEnergy + particleMass;
166 
167  // define critical gamma energies (important for integration/dicing)
168  klpm=totalEnergy*totalEnergy/lpmEnergy;
169  kp=sqrt(densityCorr);
170 }
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double GetRadlen() const
Definition: G4Material.hh:220
G4bool LPMFlag() const
Definition: G4VEmModel.hh:676
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ bremFactor

G4double G4eBremsstrahlungRelModel::bremFactor
protected

Definition at line 145 of file G4eBremsstrahlungRelModel.hh.

◆ currentZ

G4double G4eBremsstrahlungRelModel::currentZ
protected

Definition at line 151 of file G4eBremsstrahlungRelModel.hh.

◆ densityCorr

G4double G4eBremsstrahlungRelModel::densityCorr
protected

Definition at line 153 of file G4eBremsstrahlungRelModel.hh.

◆ densityFactor

G4double G4eBremsstrahlungRelModel::densityFactor
protected

Definition at line 152 of file G4eBremsstrahlungRelModel.hh.

◆ energyThresholdLPM

G4double G4eBremsstrahlungRelModel::energyThresholdLPM
private

Definition at line 167 of file G4eBremsstrahlungRelModel.hh.

◆ facFel

G4double G4eBremsstrahlungRelModel::facFel
private

Definition at line 168 of file G4eBremsstrahlungRelModel.hh.

◆ facFinel

G4double G4eBremsstrahlungRelModel::facFinel
private

Definition at line 168 of file G4eBremsstrahlungRelModel.hh.

◆ fCoulomb

G4double G4eBremsstrahlungRelModel::fCoulomb
private

Definition at line 173 of file G4eBremsstrahlungRelModel.hh.

◆ Fel

G4double G4eBremsstrahlungRelModel::Fel
private

Definition at line 173 of file G4eBremsstrahlungRelModel.hh.

◆ Fel_light

const G4double G4eBremsstrahlungRelModel::Fel_light = {0., 5.31 , 4.79 , 4.74 , 4.71}
staticprivate

Definition at line 160 of file G4eBremsstrahlungRelModel.hh.

◆ fGLPM

G4PhysicsVector * G4eBremsstrahlungRelModel::fGLPM
private

Definition at line 177 of file G4eBremsstrahlungRelModel.hh.

◆ Finel

G4double G4eBremsstrahlungRelModel::Finel
private

Definition at line 173 of file G4eBremsstrahlungRelModel.hh.

◆ Finel_light

const G4double G4eBremsstrahlungRelModel::Finel_light = {0., 6.144 , 5.621 , 5.805 , 5.924}
staticprivate

Definition at line 161 of file G4eBremsstrahlungRelModel.hh.

◆ fLPMconstant

G4double G4eBremsstrahlungRelModel::fLPMconstant
private

Definition at line 166 of file G4eBremsstrahlungRelModel.hh.

◆ fMax

G4double G4eBremsstrahlungRelModel::fMax
private

Definition at line 173 of file G4eBremsstrahlungRelModel.hh.

◆ fMigdalConstant

G4double G4eBremsstrahlungRelModel::fMigdalConstant
private

Definition at line 165 of file G4eBremsstrahlungRelModel.hh.

◆ fParticleChange

G4ParticleChangeForLoss* G4eBremsstrahlungRelModel::fParticleChange
protected

Definition at line 143 of file G4eBremsstrahlungRelModel.hh.

◆ fPhiLPM

G4PhysicsVector * G4eBremsstrahlungRelModel::fPhiLPM
private

Definition at line 177 of file G4eBremsstrahlungRelModel.hh.

◆ fXiLPM

G4PhysicsVector* G4eBremsstrahlungRelModel::fXiLPM
private

Definition at line 177 of file G4eBremsstrahlungRelModel.hh.

◆ gLPM

G4double G4eBremsstrahlungRelModel::gLPM
private

Definition at line 178 of file G4eBremsstrahlungRelModel.hh.

◆ isElectron

G4bool G4eBremsstrahlungRelModel::isElectron
protected

Definition at line 155 of file G4eBremsstrahlungRelModel.hh.

◆ kinEnergy

G4double G4eBremsstrahlungRelModel::kinEnergy
protected

Definition at line 149 of file G4eBremsstrahlungRelModel.hh.

◆ klpm

G4double G4eBremsstrahlungRelModel::klpm
private

Definition at line 181 of file G4eBremsstrahlungRelModel.hh.

◆ kp

G4double G4eBremsstrahlungRelModel::kp
private

Definition at line 181 of file G4eBremsstrahlungRelModel.hh.

◆ lnZ

G4double G4eBremsstrahlungRelModel::lnZ
private

Definition at line 172 of file G4eBremsstrahlungRelModel.hh.

◆ logTwo

G4double G4eBremsstrahlungRelModel::logTwo
private

Definition at line 169 of file G4eBremsstrahlungRelModel.hh.

◆ lowestKinEnergy

G4double G4eBremsstrahlungRelModel::lowestKinEnergy
private

Definition at line 164 of file G4eBremsstrahlungRelModel.hh.

◆ lpmEnergy

G4double G4eBremsstrahlungRelModel::lpmEnergy
private

Definition at line 176 of file G4eBremsstrahlungRelModel.hh.

◆ nist

G4NistManager* G4eBremsstrahlungRelModel::nist
protected

Definition at line 140 of file G4eBremsstrahlungRelModel.hh.

◆ particle

const G4ParticleDefinition* G4eBremsstrahlungRelModel::particle
protected

Definition at line 141 of file G4eBremsstrahlungRelModel.hh.

◆ particleMass

G4double G4eBremsstrahlungRelModel::particleMass
protected

Definition at line 148 of file G4eBremsstrahlungRelModel.hh.

◆ phiLPM

G4double G4eBremsstrahlungRelModel::phiLPM
private

Definition at line 178 of file G4eBremsstrahlungRelModel.hh.

◆ preS1

G4double G4eBremsstrahlungRelModel::preS1
private

Definition at line 169 of file G4eBremsstrahlungRelModel.hh.

◆ theGamma

G4ParticleDefinition* G4eBremsstrahlungRelModel::theGamma
protected

Definition at line 142 of file G4eBremsstrahlungRelModel.hh.

◆ totalEnergy

G4double G4eBremsstrahlungRelModel::totalEnergy
protected

Definition at line 150 of file G4eBremsstrahlungRelModel.hh.

◆ use_completescreening

G4bool G4eBremsstrahlungRelModel::use_completescreening
private

Definition at line 184 of file G4eBremsstrahlungRelModel.hh.

◆ wgi

const G4double G4eBremsstrahlungRelModel::wgi
staticprivate
Initial value:
={ 0.0506, 0.1112, 0.1569, 0.1813,
0.1813, 0.1569, 0.1112, 0.0506 }

Definition at line 159 of file G4eBremsstrahlungRelModel.hh.

◆ xgi

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

Definition at line 159 of file G4eBremsstrahlungRelModel.hh.

◆ xiLPM

G4double G4eBremsstrahlungRelModel::xiLPM
private

Definition at line 178 of file G4eBremsstrahlungRelModel.hh.

◆ z13

G4double G4eBremsstrahlungRelModel::z13
private

Definition at line 172 of file G4eBremsstrahlungRelModel.hh.

◆ z23

G4double G4eBremsstrahlungRelModel::z23
private

Definition at line 172 of file G4eBremsstrahlungRelModel.hh.


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