Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4eBremsstrahlungRelModel Class Reference

#include <G4eBremsstrahlungRelModel.hh>

Inheritance diagram for G4eBremsstrahlungRelModel:
G4VEmModel G4LivermoreBremsstrahlungModel G4SeltzerBergerModel G4ePolarizedBremsstrahlungModel

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 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 *)
 
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
 

Protected Member Functions

virtual G4double ComputeDXSectionPerAtom (G4double gammaEnergy)
 
void SetCurrentElement (const G4double)
 
- 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 currentZ
 
G4double densityFactor
 
G4double densityCorr
 
G4bool isElectron
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 62 of file G4eBremsstrahlungRelModel.hh.

Constructor & Destructor Documentation

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

Definition at line 87 of file G4eBremsstrahlungRelModel.cc.

References currentZ, densityCorr, densityFactor, fParticleChange, G4Gamma::Gamma(), G4NistManager::Instance(), kinEnergy, python.hepunit::MeV, nist, particleMass, G4VEmModel::SetAngularDistribution(), G4VEmModel::SetLowEnergyLimit(), G4VEmModel::SetLPMFlag(), theGamma, and totalEnergy.

89  : G4VEmModel(nam),
90  particle(0),
92  isElectron(true),
95  fXiLPM(0), fPhiLPM(0), fGLPM(0),
96  use_completescreening(false),isInitialised(false)
97 {
98  fParticleChange = 0;
100 
101  lowestKinEnergy = 1.0*MeV;
102  SetLowEnergyLimit(lowestKinEnergy);
103 
105 
106  SetLPMFlag(true);
107  //SetAngularDistribution(new G4ModifiedTsai());
109 
110  particleMass = kinEnergy = totalEnergy = currentZ = z13 = z23 = lnZ = Fel
111  = Finel = fCoulomb = fMax = densityFactor = densityCorr = lpmEnergy
112  = xiLPM = phiLPM = gLPM = klpm = kp = 0.0;
113 
114  energyThresholdLPM = 1.e39;
115 
116  InitialiseConstants();
117  if(p) { SetParticle(p); }
118 }
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
const G4ParticleDefinition * particle
float electron_mass_c2
Definition: hepunit.py:274
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
electron_Compton_length
Definition: hepunit.py:289
void SetLPMFlag(G4bool val)
Definition: G4VEmModel.hh:732
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
G4ParticleChangeForLoss * fParticleChange
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4eBremsstrahlungRelModel::~G4eBremsstrahlungRelModel ( )
virtual

Definition at line 133 of file G4eBremsstrahlungRelModel.cc.

134 {
135 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 279 of file G4eBremsstrahlungRelModel.cc.

References bremFactor, kinEnergy, G4VEmModel::LowEnergyLimit(), G4INCL::Math::min(), particle, and SetCurrentElement().

285 {
286  if(!particle) { SetParticle(p); }
287  if(kineticEnergy < LowEnergyLimit()) { return 0.0; }
288  G4double cut = std::min(cutEnergy, kineticEnergy);
289  G4double tmax = std::min(maxEnergy, kineticEnergy);
290 
291  if(cut >= tmax) { return 0.0; }
292 
294 
295  G4double cross = ComputeXSectionPerAtom(cut);
296 
297  // allow partial integration
298  if(tmax < kinEnergy) { cross -= ComputeXSectionPerAtom(tmax); }
299 
300  cross *= Z*Z*bremFactor;
301 
302  return cross;
303 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
const G4ParticleDefinition * particle
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
G4double G4eBremsstrahlungRelModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 209 of file G4eBremsstrahlungRelModel.cc.

References bremFactor, currentZ, G4Material::GetAtomicNumDensityVector(), G4Material::GetElementVector(), G4Material::GetNumberOfElements(), G4VEmModel::LowEnergyLimit(), G4INCL::Math::min(), particle, SetCurrentElement(), G4VEmModel::SetCurrentElement(), and SetupForMaterial().

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

Reimplemented in G4LivermoreBremsstrahlungModel, and G4SeltzerBergerModel.

Definition at line 448 of file G4eBremsstrahlungRelModel.cc.

References currentZ, python.hepunit::electron_mass_c2, main(), and totalEnergy.

Referenced by SampleSecondaries().

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

Implements G4VEmModel.

Reimplemented in G4LivermoreBremsstrahlungModel, G4SeltzerBergerModel, and G4ePolarizedBremsstrahlungModel.

Definition at line 175 of file G4eBremsstrahlungRelModel.cc.

References currentZ, fParticleChange, G4VEmModel::GetParticleChangeForLoss(), G4VEmModel::InitialiseElementSelectors(), and G4VEmModel::IsMaster().

Referenced by G4LivermoreBremsstrahlungModel::Initialise(), and G4SeltzerBergerModel::Initialise().

177 {
178  if(p) { SetParticle(p); }
179 
180  currentZ = 0.;
181 
182  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
183 
184  if(isInitialised) { return; }
186  isInitialised = true;
187 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4ParticleChangeForLoss * fParticleChange
void G4eBremsstrahlungRelModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 191 of file G4eBremsstrahlungRelModel.cc.

References G4VEmModel::GetElementSelectors(), and G4VEmModel::SetElementSelectors().

193 {
195 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
G4double G4eBremsstrahlungRelModel::LowestKinEnergy ( ) const
inline

Definition at line 267 of file G4eBremsstrahlungRelModel.hh.

Referenced by G4SeltzerBergerModel::G4SeltzerBergerModel().

268 {
269  return lowestKinEnergy;
270 }
G4double G4eBremsstrahlungRelModel::LPMconstant ( ) const
inline

Definition at line 255 of file G4eBremsstrahlungRelModel.hh.

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

Reimplemented from G4VEmModel.

Definition at line 200 of file G4eBremsstrahlungRelModel.cc.

References G4INCL::Math::max().

203 {
204  return std::max(lowestKinEnergy, cut);
205 }
T max(const T t1, const T t2)
brief Return the largest of the two arguments
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 483 of file G4eBremsstrahlungRelModel.cc.

References ComputeDXSectionPerAtom(), currentZ, densityCorr, densityFactor, python.hepunit::electron_mass_c2, fParticleChange, fStopAndKill, G4cout, G4endl, G4Exp(), G4Log(), G4lrint(), G4UniformRand, G4VEmModel::GetAngularDistribution(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4VEmModel::GetName(), kinEnergy, G4VEmModel::LowEnergyLimit(), G4INCL::Math::min(), particle, particleMass, G4VParticleChange::ProposeTrackStatus(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SecondaryThreshold(), G4VEmModel::SelectRandomAtom(), SetCurrentElement(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), SetupForMaterial(), theGamma, totalEnergy, and test::x.

489 {
490  G4double kineticEnergy = dp->GetKineticEnergy();
491  if(kineticEnergy < LowEnergyLimit()) { return; }
492  G4double cut = std::min(cutEnergy, kineticEnergy);
493  G4double emax = std::min(maxEnergy, kineticEnergy);
494  if(cut >= emax) { return; }
495 
496  SetupForMaterial(particle, couple->GetMaterial(), kineticEnergy);
497 
498  const G4Element* elm =
499  SelectRandomAtom(couple,particle,kineticEnergy,cut,emax);
500  SetCurrentElement(elm->GetZ());
501 
502  kinEnergy = kineticEnergy;
503  totalEnergy = kineticEnergy + particleMass;
505 
506  //G4double fmax= fMax;
507  G4bool highe = true;
508  if(totalEnergy < energyThresholdLPM) { highe = false; }
509 
510  G4double xmin = G4Log(cut*cut + densityCorr);
511  G4double xmax = G4Log(emax*emax + densityCorr);
512  G4double gammaEnergy, f, x;
513 
514  do {
515  x = G4Exp(xmin + G4UniformRand()*(xmax - xmin)) - densityCorr;
516  if(x < 0.0) { x = 0.0; }
517  gammaEnergy = sqrt(x);
518  if(highe) { f = ComputeRelDXSectionPerAtom(gammaEnergy); }
519  else { f = ComputeDXSectionPerAtom(gammaEnergy); }
520 
521  if ( f > fMax ) {
522  G4cout << "### G4eBremsstrahlungRelModel Warning: Majoranta exceeded! "
523  << f << " > " << fMax
524  << " Egamma(MeV)= " << gammaEnergy
525  << " Ee(MeV)= " << kineticEnergy
526  << " " << GetName()
527  << G4endl;
528  }
529 
530  } while (f < fMax*G4UniformRand());
531 
532  //
533  // angles of the emitted gamma. ( Z - axis along the parent particle)
534  // use general interface
535  //
536 
537  G4ThreeVector gammaDirection =
538  GetAngularDistribution()->SampleDirection(dp, totalEnergy-gammaEnergy,
539  G4lrint(currentZ),
540  couple->GetMaterial());
541 
542  // create G4DynamicParticle object for the Gamma
543  G4DynamicParticle* gamma = new G4DynamicParticle(theGamma,gammaDirection,
544  gammaEnergy);
545  vdp->push_back(gamma);
546 
547  G4double totMomentum = sqrt(kineticEnergy*(totalEnergy + electron_mass_c2));
548  G4ThreeVector direction = (totMomentum*dp->GetMomentumDirection()
549  - gammaEnergy*gammaDirection).unit();
550 
551  // energy of primary
552  G4double finalE = kineticEnergy - gammaEnergy;
553 
554  // stop tracking and create new secondary instead of primary
555  if(gammaEnergy > SecondaryThreshold()) {
558  G4DynamicParticle* el =
559  new G4DynamicParticle(const_cast<G4ParticleDefinition*>(particle),
560  direction, finalE);
561  vdp->push_back(el);
562 
563  // continue tracking
564  } else {
567  }
568 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4double SecondaryThreshold() const
Definition: G4VEmModel.hh:627
G4double GetKineticEnergy() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double)
const G4ParticleDefinition * particle
#define G4UniformRand()
Definition: Randomize.hh:87
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)
float electron_mass_c2
Definition: hepunit.py:274
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
int G4lrint(double ad)
Definition: templates.hh:163
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:753
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:510
const G4Material * GetMaterial() const
void G4eBremsstrahlungRelModel::SetCurrentElement ( const G4double  Z)
inlineprotected

Definition at line 190 of file G4eBremsstrahlungRelModel.hh.

References currentZ, G4VEmModel::GetCurrentElement(), G4Element::GetfCoulomb(), G4NistManager::GetLOGZ(), G4NistManager::GetZ13(), iz, and nist.

Referenced by ComputeCrossSectionPerAtom(), ComputeDEDXPerVolume(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), and SampleSecondaries().

191 {
192  if(Z != currentZ) {
193  currentZ = Z;
194 
195  G4int iz = G4int(Z);
196  z13 = nist->GetZ13(iz);
197  z23 = z13*z13;
198  lnZ = nist->GetLOGZ(iz);
199 
200  if (iz <= 4) {
201  Fel = Fel_light[iz];
202  Finel = Finel_light[iz] ;
203  }
204  else {
205  Fel = facFel - lnZ/3. ;
206  Finel = facFinel - 2.*lnZ/3. ;
207  }
208 
209  fCoulomb = GetCurrentElement()->GetfCoulomb();
210  fMax = Fel-fCoulomb + Finel/currentZ + (1.+1./currentZ)/12.;
211  }
212 }
G4double GetZ13(G4double Z)
G4double GetfCoulomb() const
Definition: G4Element.hh:190
int G4int
Definition: G4Types.hh:78
G4double iz
Definition: TRTMaterials.hh:39
G4double GetLOGZ(G4int Z)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440
void G4eBremsstrahlungRelModel::SetLowestKinEnergy ( G4double  val)
inline

Definition at line 260 of file G4eBremsstrahlungRelModel.hh.

Referenced by G4SeltzerBergerModel::G4SeltzerBergerModel().

261 {
262  lowestKinEnergy = val;
263 }
void G4eBremsstrahlungRelModel::SetLPMconstant ( G4double  val)
inline

Definition at line 247 of file G4eBremsstrahlungRelModel.hh.

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

Reimplemented from G4VEmModel.

Definition at line 149 of file G4eBremsstrahlungRelModel.cc.

References densityCorr, densityFactor, G4Material::GetElectronDensity(), G4Material::GetRadlen(), kinEnergy, G4VEmModel::LPMFlag(), particleMass, and totalEnergy.

Referenced by ComputeDEDXPerVolume(), G4SeltzerBergerModel::SampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), and SampleSecondaries().

152 {
153  densityFactor = mat->GetElectronDensity()*fMigdalConstant;
154  lpmEnergy = mat->GetRadlen()*fLPMconstant;
155 
156  // Threshold for LPM effect (i.e. below which LPM hidden by density effect)
157  if (LPMFlag()) {
158  energyThresholdLPM=sqrt(densityFactor)*lpmEnergy;
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 }
G4bool LPMFlag() const
Definition: G4VEmModel.hh:634
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double GetRadlen() const
Definition: G4Material.hh:218

Field Documentation

G4double G4eBremsstrahlungRelModel::bremFactor
protected
G4double G4eBremsstrahlungRelModel::currentZ
protected
G4double G4eBremsstrahlungRelModel::densityCorr
protected
G4double G4eBremsstrahlungRelModel::densityFactor
protected
G4ParticleChangeForLoss* G4eBremsstrahlungRelModel::fParticleChange
protected
G4bool G4eBremsstrahlungRelModel::isElectron
protected
G4double G4eBremsstrahlungRelModel::kinEnergy
protected
G4NistManager* G4eBremsstrahlungRelModel::nist
protected

Definition at line 140 of file G4eBremsstrahlungRelModel.hh.

Referenced by G4eBremsstrahlungRelModel(), and SetCurrentElement().

const G4ParticleDefinition* G4eBremsstrahlungRelModel::particle
protected
G4double G4eBremsstrahlungRelModel::particleMass
protected
G4ParticleDefinition* G4eBremsstrahlungRelModel::theGamma
protected
G4double G4eBremsstrahlungRelModel::totalEnergy
protected

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