Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double tkin, G4double Z, G4double, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double) override
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut) override
 
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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 * > *)
 
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)
 
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=nullptr)
 
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 (G4int)
 
- 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
 
const G4ParticleDefinitionparticle
 
G4ParticleDefinitiontheGamma
 
G4ParticleChangeForLossfParticleChange
 
G4double bremFactor
 
G4double particleMass
 
G4double kinEnergy
 
G4double totalEnergy
 
G4double densityFactor
 
G4double densityCorr
 
G4int currentZ
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

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::G4eBremsstrahlungRelModel ( const G4ParticleDefinition p = 0,
const G4String nam = "eBremLPM" 
)
explicit

Definition at line 90 of file G4eBremsstrahlungRelModel.cc.

92  : G4VEmModel(nam),
93  particle(0),
95  isElectron(true),
98  use_completescreening(false)
99 {
100  fParticleChange = nullptr;
102 
103  lowestKinEnergy = 1.0*MeV;
104  SetLowEnergyLimit(lowestKinEnergy);
105 
107 
108  SetLPMFlag(true);
109  //SetAngularDistribution(new G4ModifiedTsai());
111 
112  particleMass = kinEnergy = totalEnergy = z13 = z23 = lnZ = Fel
113  = Finel = fCoulomb = fMax = densityFactor = densityCorr = lpmEnergy
114  = xiLPM = phiLPM = gLPM = klpm = kp = 0.0;
115  currentZ = 0;
116  energyThresholdLPM = 1.e39;
117 
118  InitialiseConstants();
119  if(p) { SetParticle(p); }
120 }
static constexpr double hbarc
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
const G4ParticleDefinition * particle
static constexpr double electron_mass_c2
static constexpr double classic_electr_radius
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:773
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
static constexpr double MeV
Definition: G4SIunits.hh:214
G4ParticleChangeForLoss * fParticleChange
static constexpr double pi
Definition: G4SIunits.hh:75
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double fine_structure_const
static constexpr double electron_Compton_length

Here is the call graph for this function:

G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel ( )
virtual

Definition at line 135 of file G4eBremsstrahlungRelModel.cc.

136 {
137 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 282 of file G4eBremsstrahlungRelModel.cc.

288 {
289  if(!particle) { SetParticle(p); }
290  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
291  G4double cut = std::min(cutEnergy, kineticEnergy);
292  G4double tmax = std::min(maxEnergy, kineticEnergy);
293 
294  if(cut >= tmax) { return 0.0; }
295 
297 
298  G4double cross = ComputeXSectionPerAtom(cut);
299 
300  // allow partial integration
301  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
302 
303  cross *= Z*Z*bremFactor;
304 
305  return cross;
306 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
const G4ParticleDefinition * particle
int G4lrint(double ad)
Definition: templates.hh:163
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 212 of file G4eBremsstrahlungRelModel.cc.

217 {
218  if(!particle) { SetParticle(p); }
219  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
220  G4double cut = std::min(cutEnergy, kineticEnergy);
221  if(cut == 0.0) { return 0.0; }
222 
223  SetupForMaterial(particle, material,kineticEnergy);
224 
225  const G4ElementVector* theElementVector = material->GetElementVector();
226  const G4double* theAtomicNumDensityVector = material->GetAtomicNumDensityVector();
227 
228  G4double dedx = 0.0;
229 
230  // loop for elements in the material
231  for (size_t i=0; i<material->GetNumberOfElements(); i++) {
232 
233  G4VEmModel::SetCurrentElement((*theElementVector)[i]);
234  SetCurrentElement((*theElementVector)[i]->GetZasInt());
235 
236  dedx += theAtomicNumDensityVector[i]*(currentZ*currentZ)*ComputeBremLoss(cut);
237  }
238  dedx *= bremFactor;
239 
240 
241  return dedx;
242 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
const G4ParticleDefinition * particle
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:458

Here is the call graph for this function:

G4double G4eBremsstrahlungRelModel::ComputeDXSectionPerAtom ( G4double  gammaEnergy)
protectedvirtual

Reimplemented in G4LivermoreBremsstrahlungModel, and G4SeltzerBergerModel.

Definition at line 452 of file G4eBremsstrahlungRelModel.cc.

457 {
458 
459  if(gammaEnergy < 0.0) { return 0.0; }
460 
461  G4double y = gammaEnergy/totalEnergy;
462 
463  G4double main=0.,secondTerm=0.;
464 
465  G4double currZ = (G4double)currentZ;
466  if (use_completescreening|| currentZ<5) {
467  // ** form factors complete screening case **
468  main = (3./4.*y*y - y + 1.) * ( (Fel-fCoulomb) + Finel/currZ );
469  secondTerm = (1.-y)/12.*(1.+1./currZ);
470  }
471  else {
472  // ** intermediate screening using Thomas-Fermi FF from Tsai only valid for Z>=5**
473  G4double dd=100.*electron_mass_c2*y/(totalEnergy-gammaEnergy);
474  G4double gg=dd/z13;
475  G4double eps=dd/z23;
476  G4double phi1=Phi1(gg,currZ), phi1m2=Phi1M2(gg,currZ);
477  G4double psi1=Psi1(eps,currZ), psi1m2=Psi1M2(eps,currZ);
478 
479  main = (3./4.*y*y - y + 1.) *
480  ( (0.25*phi1-1./3.*lnZ-fCoulomb) + (0.25*psi1-2./3.*lnZ)/currZ );
481  secondTerm = (1.-y)/8.*(phi1m2+psi1m2/currZ);
482  }
483  G4double cross = main+secondTerm;
484  return cross;
485 }
static const G4double eps
static constexpr double electron_mass_c2
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel, G4SeltzerBergerModel, and G4ePolarizedBremsstrahlungModel.

Definition at line 176 of file G4eBremsstrahlungRelModel.cc.

178 {
179  if(p) { SetParticle(p); }
180 
181  currentZ = 0;
182 
183  if(IsMaster() && LowEnergyLimit() < HighEnergyLimit()) {
184  InitialiseElementSelectors(p, cuts);
185  }
186 
188 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
G4ParticleChangeForLoss * fParticleChange

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 192 of file G4eBremsstrahlungRelModel.cc.

194 {
195  if(LowEnergyLimit() < HighEnergyLimit()) {
197  }
198 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809

Here is the call graph for this function:

G4double G4eBremsstrahlungRelModel::LowestKinEnergy ( ) const
inline

Definition at line 266 of file G4eBremsstrahlungRelModel.hh.

267 {
268  return lowestKinEnergy;
269 }

Here is the caller graph for this function:

G4double G4eBremsstrahlungRelModel::LPMconstant ( ) const
inline

Definition at line 254 of file G4eBremsstrahlungRelModel.hh.

255 {
256  return fLPMconstant;
257 }
G4double G4eBremsstrahlungRelModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double  cut 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 203 of file G4eBremsstrahlungRelModel.cc.

206 {
207  return std::max(lowestKinEnergy, cut);
208 }
T max(const T t1, const T t2)
brief Return the largest of the two arguments

Here is the call graph for this function:

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel, G4SeltzerBergerModel, and G4ePolarizedBremsstrahlungModel.

Definition at line 489 of file G4eBremsstrahlungRelModel.cc.

495 {
496  G4double kineticEnergy = dp->GetKineticEnergy();
497  if(kineticEnergy < LowEnergyLimit()) { return; }
498  G4double cut = std::min(cutEnergy, kineticEnergy);
499  G4double emax = std::min(maxEnergy, kineticEnergy);
500  if(cut >= emax) { return; }
501 
502  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
503 
504  const G4Element* elm =
505  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
506  SetCurrentElement(elm->GetZasInt());
507 
508  kinEnergy = kineticEnergy;
509  totalEnergy = kineticEnergy + particleMass;
511 
512  //G4double fmax= fMax;
513  G4bool highe = true;
514  if(totalEnergy < energyThresholdLPM) { highe = false; }
515 
516  G4double xmin = G4Log(cut*cut + densityCorr);
517  G4double xmax = G4Log(emax*emax + densityCorr);
518  G4double gammaEnergy, f, x;
519 
520  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
521 
522  do {
523  x = G4Exp(xmin + rndmEngine->flat()*(xmax - xmin)) - densityCorr;
524  if(x < 0.0) { x = 0.0; }
525  gammaEnergy = sqrt(x);
526  if(highe) { f = ComputeRelDXSectionPerAtom(gammaEnergy); }
527  else { f = ComputeDXSectionPerAtom(gammaEnergy); }
528 
529  if ( f > fMax ) {
530  G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
531  << f << " > " << fMax
532  << " Egamma(MeV)= " << gammaEnergy
533  << " Ee(MeV)= " << kineticEnergy
534  << " " << GetName()
535  << G4endl;
536  }
537 
538  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
539  } while (f < fMax*rndmEngine->flat());
540 
541  //
542  // angles of the emitted gamma. ( Z - axis along the parent particle)
543  // use general interface
544  //
545 
546  G4ThreeVector gammaDirection =
547  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
548  currentZ,
549  couple->GetMaterial());
550 
551  // create G4DynamicParticle object for the Gamma
552  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
553  gammaEnergy);
554  vdp->push_back(gamma);
555 
556  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
557  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
558  - gammaEnergy*gammaDirection).unit();
559 
560  // energy of primary
561  G4double finalE = kineticEnergy - gammaEnergy;
562 
563  // stop tracking and create new secondary instead of primary
564  if(gammaEnergy > SecondaryThreshold()) {
567  G4DynamicParticle* el =
568  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
569  direction, finalE);
570  vdp->push_back(el);
571 
572  // continue tracking
573  } else {
576  }
577 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:668
G4double GetKineticEnergy() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
virtual double flat()=0
const G4ParticleDefinition * particle
static constexpr double electron_mass_c2
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double flat()
Definition: G4AblaRandom.cc:47
void SetProposedMomentumDirection(const G4ThreeVector &dir)
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
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double) override
virtual G4double ComputeDXSectionPerAtom(G4double gammaEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4VEmModel.hh:794
G4ParticleChangeForLoss * fParticleChange
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
const G4Material * GetMaterial() const

Here is the call graph for this function:

void G4eBremsstrahlungRelModel::SetCurrentElement ( G4int  Z)
inlineprotected

Definition at line 188 of file G4eBremsstrahlungRelModel.hh.

189 {
190  if(Z != currentZ) {
191  currentZ = Z;
192 
193  z13 = nist->GetZ13(Z);
194  z23 = z13*z13;
195  lnZ = nist->GetLOGZ(Z);
196 
197  if (Z <= 4) {
198  Fel = Fel_light[Z];
199  Finel = Finel_light[Z] ;
200  }
201  else {
202  G4double lnzt = lnZ/3.;
203  Fel = facFel - lnzt;
204  Finel = facFinel - 2*lnzt;
205  }
206 
207  fCoulomb = GetCurrentElement()->GetfCoulomb();
208  G4double xz = 1.0/(G4double)Z;
209  fMax = Fel-fCoulomb + Finel*xz + (1. + xz)/12.;
210  }
211 }
G4double GetfCoulomb() const
Definition: G4Element.hh:191
G4double GetZ13(G4double Z) const
G4double GetLOGZ(G4int Z) const
double G4double
Definition: G4Types.hh:76
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:466

Here is the call graph for this function:

Here is the caller graph for this function:

void G4eBremsstrahlungRelModel::SetLowestKinEnergy ( G4double  val)
inline

Definition at line 259 of file G4eBremsstrahlungRelModel.hh.

260 {
261  lowestKinEnergy = val;
262 }

Here is the caller graph for this function:

void G4eBremsstrahlungRelModel::SetLPMconstant ( G4double  val)
inline

Definition at line 246 of file G4eBremsstrahlungRelModel.hh.

247 {
248  fLPMconstant = val;
249 }
void G4eBremsstrahlungRelModel::SetupForMaterial ( const G4ParticleDefinition ,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 151 of file G4eBremsstrahlungRelModel.cc.

154 {
155  densityFactor = mat->GetElectronDensity()*fMigdalConstant;
156  lpmEnergy = mat->GetRadlen()*fLPMconstant;
157 
158  // Threshold for LPM effect (i.e. below which LPM hidden by density effect)
159  if (LPMFlag()) {
160  energyThresholdLPM=sqrt(densityFactor)*lpmEnergy;
161  } else {
162  energyThresholdLPM=1.e39; // i.e. do not use LPM effect
163  }
164  // calculate threshold for density effect
165  kinEnergy = kineticEnergy;
166  totalEnergy = kineticEnergy + particleMass;
168 
169  // define critical gamma energies (important for integration/dicing)
170  klpm=totalEnergy*totalEnergy/lpmEnergy;
171  kp=sqrt(densityCorr);
172 }
G4bool LPMFlag() const
Definition: G4VEmModel.hh:675
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double GetRadlen() const
Definition: G4Material.hh:220

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

G4double G4eBremsstrahlungRelModel::bremFactor
protected

Definition at line 145 of file G4eBremsstrahlungRelModel.hh.

G4int G4eBremsstrahlungRelModel::currentZ
protected

Definition at line 154 of file G4eBremsstrahlungRelModel.hh.

G4double G4eBremsstrahlungRelModel::densityCorr
protected

Definition at line 152 of file G4eBremsstrahlungRelModel.hh.

G4double G4eBremsstrahlungRelModel::densityFactor
protected

Definition at line 151 of file G4eBremsstrahlungRelModel.hh.

G4ParticleChangeForLoss* G4eBremsstrahlungRelModel::fParticleChange
protected

Definition at line 143 of file G4eBremsstrahlungRelModel.hh.

G4bool G4eBremsstrahlungRelModel::isElectron
protected

Definition at line 155 of file G4eBremsstrahlungRelModel.hh.

G4double G4eBremsstrahlungRelModel::kinEnergy
protected

Definition at line 149 of file G4eBremsstrahlungRelModel.hh.

G4NistManager* G4eBremsstrahlungRelModel::nist
protected

Definition at line 140 of file G4eBremsstrahlungRelModel.hh.

const G4ParticleDefinition* G4eBremsstrahlungRelModel::particle
protected

Definition at line 141 of file G4eBremsstrahlungRelModel.hh.

G4double G4eBremsstrahlungRelModel::particleMass
protected

Definition at line 148 of file G4eBremsstrahlungRelModel.hh.

G4ParticleDefinition* G4eBremsstrahlungRelModel::theGamma
protected

Definition at line 142 of file G4eBremsstrahlungRelModel.hh.

G4double G4eBremsstrahlungRelModel::totalEnergy
protected

Definition at line 150 of file G4eBremsstrahlungRelModel.hh.


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