Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4mplIonisationWithDeltaModel Class Reference

#include <G4mplIonisationWithDeltaModel.hh>

Inheritance diagram for G4mplIonisationWithDeltaModel:
Collaboration diagram for G4mplIonisationWithDeltaModel:

Public Member Functions

 G4mplIonisationWithDeltaModel (G4double mCharge, const G4String &nam="mplIonisationWithDelta")
 
virtual ~G4mplIonisationWithDeltaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss) override
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length) override
 
void SetParticle (const G4ParticleDefinition *p)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
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 MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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)
 
- Public Member Functions inherited from G4VEmFluctuationModel
 G4VEmFluctuationModel (const G4String &nam)
 
virtual ~G4VEmFluctuationModel ()
 
virtual void InitialiseMe (const G4ParticleDefinition *)
 
virtual void SetParticleAndCharge (const G4ParticleDefinition *, G4double q2)
 
const G4StringGetName () const
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) override
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 58 of file G4mplIonisationWithDeltaModel.hh.

Constructor & Destructor Documentation

G4mplIonisationWithDeltaModel::G4mplIonisationWithDeltaModel ( G4double  mCharge,
const G4String nam = "mplIonisationWithDelta" 
)
explicit

Definition at line 73 of file G4mplIonisationWithDeltaModel.cc.

76  magCharge(mCharge),
77  twoln10(log(100.0)),
78  betalow(0.01),
79  betalim(0.1),
80  beta2lim(betalim*betalim),
81  bg2lim(beta2lim*(1.0 + beta2lim))
82 {
83  nmpl = G4lrint(std::fabs(magCharge) * 2 * fine_structure_const);
84  if(nmpl > 6) { nmpl = 6; }
85  else if(nmpl < 1) { nmpl = 1; }
86  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
87  chargeSquare = magCharge * magCharge;
88  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
89  fParticleChange = nullptr;
90  theElectron = G4Electron::Electron();
91  G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
92  << magCharge/eplus << G4endl;
93  monopole = nullptr;
94  mass = 0.0;
95 }
static constexpr double cm2
Definition: G4SIunits.hh:120
static constexpr double hbarc
static constexpr double g
Definition: G4SIunits.hh:183
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double electron_mass_c2
G4GLOB_DLL std::ostream G4cout
static constexpr double eplus
Definition: G4SIunits.hh:199
int G4lrint(double ad)
Definition: templates.hh:163
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
G4VEmFluctuationModel(const G4String &nam)
static constexpr double fine_structure_const

Here is the call graph for this function:

G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel ( )
virtual

Definition at line 99 of file G4mplIonisationWithDeltaModel.cc.

100 {
101  if(IsMaster()) { delete dedx0; }
102 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717

Here is the call graph for this function:

Member Function Documentation

G4double G4mplIonisationWithDeltaModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 248 of file G4mplIonisationWithDeltaModel.cc.

254 {
255  G4double cross =
256  Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
257  return cross;
258 }
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 228 of file G4mplIonisationWithDeltaModel.cc.

233 {
234  if(!monopole) { SetParticle(p); }
235  G4double cross = 0.0;
236  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
237  G4double maxEnergy = std::min(tmax,maxKinEnergy);
238  G4double cutEnergy = std::max(LowEnergyLimit(), cut);
239  if(cutEnergy < maxEnergy) {
240  cross = (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl;
241  }
242  return cross;
243 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
void SetParticle(const G4ParticleDefinition *p)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
T max(const T t1, const T t2)
brief Return the largest of the two arguments
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:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 150 of file G4mplIonisationWithDeltaModel.cc.

154 {
155  if(!monopole) { SetParticle(p); }
156  G4double tmax = MaxSecondaryEnergy(p,kineticEnergy);
157  G4double cutEnergy = std::min(tmax, maxEnergy);
158  cutEnergy = std::max(LowEnergyLimit(), cutEnergy);
159  G4double tau = kineticEnergy / mass;
160  G4double gam = tau + 1.0;
161  G4double bg2 = tau * (tau + 2.0);
162  G4double beta2 = bg2 / (gam * gam);
163  G4double beta = sqrt(beta2);
164 
165  // low-energy asymptotic formula
166  //G4double dedx = dedxlim*beta*material->GetDensity();
167  G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
168 
169  // above asymptotic
170  if(beta > betalow) {
171 
172  // high energy
173  if(beta >= betalim) {
174  dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
175 
176  } else {
177 
178  //G4double dedx1 = dedxlim*betalow*material->GetDensity();
179  G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
180  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
181 
182  // extrapolation between two formula
183  G4double kapa2 = beta - betalow;
184  G4double kapa1 = betalim - beta;
185  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
186  }
187  }
188  return dedx;
189 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
void SetParticle(const G4ParticleDefinition *p)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
T max(const T t1, const T t2)
brief Return the largest of the two arguments
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 G4mplIonisationWithDeltaModel::Dispersion ( const G4Material material,
const G4DynamicParticle dp,
G4double  tmax,
G4double  length 
)
overridevirtual

Implements G4VEmFluctuationModel.

Definition at line 352 of file G4mplIonisationWithDeltaModel.cc.

356 {
357  G4double siga = 0.0;
358  G4double tau = dp->GetKineticEnergy()/mass;
359  if(tau > 0.0) {
360  G4double electronDensity = material->GetElectronDensity();
361  G4double gam = tau + 1.0;
362  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
363  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
364  * electronDensity * chargeSquare;
365  }
366  return siga;
367 }
G4double GetKineticEnergy() const
static constexpr double twopi_mc2_rcl2
G4double GetElectronDensity() const
Definition: G4Material.hh:217
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4mplIonisationWithDeltaModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 121 of file G4mplIonisationWithDeltaModel.cc.

123 {
124  if(!monopole) { SetParticle(p); }
125  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
126  if(IsMaster()) {
127  if(!dedx0) { dedx0 = new std::vector<G4double>; }
128  G4ProductionCutsTable* theCoupleTable=
130  G4int numOfCouples = theCoupleTable->GetTableSize();
131  G4int n = dedx0->size();
132  if(n < numOfCouples) { dedx0->resize(numOfCouples); }
133 
134  // initialise vector
135  for(G4int i=0; i<numOfCouples; ++i) {
136 
137  const G4Material* material =
138  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
139  G4double eDensity = material->GetElectronDensity();
140  G4double vF = electron_Compton_length*pow(3.*pi*pi*eDensity,0.3333333333);
141  (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
142  (G4Log(2*vF/fine_structure_const) - 0.5)/vF;
143  }
144  }
145 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void SetParticle(const G4ParticleDefinition *p)
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
int G4int
Definition: G4Types.hh:78
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
static constexpr double fine_structure_const
static constexpr double electron_Compton_length
const G4Material * GetMaterial() const

Here is the call graph for this function:

G4double G4mplIonisationWithDeltaModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
)
overrideprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 372 of file G4mplIonisationWithDeltaModel.cc.

374 {
375  G4double tau = kinEnergy/mass;
376  return 2.0*electron_mass_c2*tau*(tau + 2.);
377 }
static constexpr double electron_mass_c2
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4mplIonisationWithDeltaModel::SampleFluctuations ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmax,
G4double  length,
G4double  meanLoss 
)
overridevirtual

Implements G4VEmFluctuationModel.

Definition at line 321 of file G4mplIonisationWithDeltaModel.cc.

327 {
328  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
329  G4double loss = meanLoss;
330  siga = sqrt(siga);
331  G4double twomeanLoss = meanLoss + meanLoss;
332 
333  if(twomeanLoss < siga) {
334  G4double x;
335  do {
336  loss = twomeanLoss*G4UniformRand();
337  x = (loss - meanLoss)/siga;
338  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
339  } while (1.0 - 0.5*x*x < G4UniformRand());
340  } else {
341  do {
342  loss = G4RandGauss::shoot(meanLoss,siga);
343  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
344  } while (0.0 > loss || loss > twomeanLoss);
345  }
346  return loss;
347 }
ThreeVector shoot(const G4int Ap, const G4int Af)
#define G4UniformRand()
Definition: Randomize.hh:97
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length) override
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const

Here is the call graph for this function:

void G4mplIonisationWithDeltaModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
overridevirtual

Implements G4VEmModel.

Definition at line 263 of file G4mplIonisationWithDeltaModel.cc.

268 {
269  G4double kineticEnergy = dp->GetKineticEnergy();
270  G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
271 
272  G4double maxKinEnergy = std::min(maxEnergy,tmax);
273  if(minKinEnergy >= maxKinEnergy) { return; }
274 
275  //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
276  // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
277  // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
278 
279  G4double totEnergy = kineticEnergy + mass;
280  G4double etot2 = totEnergy*totEnergy;
281  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
282 
283  // sampling without nuclear size effect
284  G4double q = G4UniformRand();
285  G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
286  /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
287 
288  // delta-electron is produced
289  G4double totMomentum = totEnergy*sqrt(beta2);
290  G4double deltaMomentum =
291  sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
292  G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
293  (deltaMomentum * totMomentum);
294  if(cost > 1.0) { cost = 1.0; }
295 
296  G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
297 
298  G4double phi = twopi * G4UniformRand() ;
299 
300  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
301  G4ThreeVector direction = dp->GetMomentumDirection();
302  deltaDirection.rotateUz(direction);
303 
304  // create G4DynamicParticle object for delta ray
305  G4DynamicParticle* delta =
306  new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
307 
308  vdp->push_back(delta);
309 
310  // Change kinematics of primary particle
311  kineticEnergy -= deltaKinEnergy;
312  G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
313  finalP = finalP.unit();
314 
315  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
316  fParticleChange->SetProposedMomentumDirection(finalP);
317 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
Hep3Vector unit() const
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:

void G4mplIonisationWithDeltaModel::SetParticle ( const G4ParticleDefinition p)

Definition at line 106 of file G4mplIonisationWithDeltaModel.cc.

107 {
108  monopole = p;
109  mass = monopole->GetPDGMass();
110  G4double emin =
111  std::min(LowEnergyLimit(),0.1*mass*(1./sqrt(1. - betalow*betalow) - 1.));
112  G4double emax =
113  std::max(HighEnergyLimit(),10*mass*(1./sqrt(1. - beta2lim) - 1.));
114  SetLowEnergyLimit(emin);
115  SetHighEnergyLimit(emax);
116 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
const char * p
Definition: xmltok.h:285
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
static const G4double emax
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731

Here is the call graph for this function:

Here is the caller graph for this function:


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