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

#include <G4mplIonisationWithDeltaModel.hh>

Inheritance diagram for G4mplIonisationWithDeltaModel:
G4VEmModel G4VEmFluctuationModel

Public Member Functions

 G4mplIonisationWithDeltaModel (G4double mCharge, const G4String &nam="mplIonisationWithDelta")
 
virtual ~G4mplIonisationWithDeltaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
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)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss)
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)
 
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 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 *)
 
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
 
- Public Member Functions inherited from G4VEmFluctuationModel
 G4VEmFluctuationModel (const G4String &nam)
 
virtual ~G4VEmFluctuationModel ()
 
virtual void InitialiseMe (const G4ParticleDefinition *)
 
virtual void SetParticleAndCharge (const G4ParticleDefinition *, G4double q2)
 
G4String GetName () const
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
 
- 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
 

Detailed Description

Definition at line 58 of file G4mplIonisationWithDeltaModel.hh.

Constructor & Destructor Documentation

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

Definition at line 74 of file G4mplIonisationWithDeltaModel.cc.

References python.hepunit::cm2, G4Electron::Electron(), python.hepunit::electron_mass_c2, python.hepunit::eplus, python.hepunit::fine_structure_const, g(), G4cout, G4endl, G4lrint(), python.hepunit::GeV, python.hepunit::hbarc, and python.hepunit::pi.

77  magCharge(mCharge),
78  twoln10(log(100.0)),
79  betalow(0.01),
80  betalim(0.1),
81  beta2lim(betalim*betalim),
82  bg2lim(beta2lim*(1.0 + beta2lim))
83 {
84  nmpl = G4lrint(std::fabs(magCharge) * 2 * fine_structure_const);
85  if(nmpl > 6) { nmpl = 6; }
86  else if(nmpl < 1) { nmpl = 1; }
87  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
88  chargeSquare = magCharge * magCharge;
89  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
90  fParticleChange = 0;
91  theElectron = G4Electron::Electron();
92  G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
93  << magCharge/eplus << G4endl;
94  monopole = 0;
95  mass = 0.0;
96 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
float electron_mass_c2
Definition: hepunit.py:274
int G4lrint(double ad)
Definition: templates.hh:163
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4VEmFluctuationModel(const G4String &nam)
G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel ( )
virtual

Definition at line 100 of file G4mplIonisationWithDeltaModel.cc.

References G4VEmModel::IsMaster().

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

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 249 of file G4mplIonisationWithDeltaModel.cc.

References ComputeCrossSectionPerElectron().

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

Definition at line 229 of file G4mplIonisationWithDeltaModel.cc.

References G4VEmModel::LowEnergyLimit(), G4INCL::Math::max(), MaxSecondaryEnergy(), G4INCL::Math::min(), and SetParticle().

Referenced by ComputeCrossSectionPerAtom().

234 {
235  if(!monopole) { SetParticle(p); }
236  G4double cross = 0.0;
237  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
238  G4double maxEnergy = std::min(tmax,maxKinEnergy);
239  G4double cutEnergy = std::max(LowEnergyLimit(), cut);
240  if(cutEnergy < maxEnergy) {
241  cross = (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl;
242  }
243  return cross;
244 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
void SetParticle(const G4ParticleDefinition *p)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
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
G4double G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 151 of file G4mplIonisationWithDeltaModel.cc.

References G4VEmModel::CurrentCouple(), G4InuclParticleNames::gam, G4MaterialCutsCouple::GetIndex(), G4VEmModel::LowEnergyLimit(), G4INCL::Math::max(), MaxSecondaryEnergy(), G4INCL::Math::min(), and SetParticle().

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

Implements G4VEmFluctuationModel.

Definition at line 351 of file G4mplIonisationWithDeltaModel.cc.

References G4InuclParticleNames::gam, G4Material::GetElectronDensity(), G4DynamicParticle::GetKineticEnergy(), and python.hepunit::twopi_mc2_rcl2.

Referenced by SampleFluctuations().

355 {
356  G4double siga = 0.0;
357  G4double tau = dp->GetKineticEnergy()/mass;
358  if(tau > 0.0) {
359  G4double electronDensity = material->GetElectronDensity();
360  G4double gam = tau + 1.0;
361  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
362  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
363  * electronDensity * chargeSquare;
364  }
365  return siga;
366 }
G4double GetKineticEnergy() const
G4double GetElectronDensity() const
Definition: G4Material.hh:215
double G4double
Definition: G4Types.hh:76
void G4mplIonisationWithDeltaModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 122 of file G4mplIonisationWithDeltaModel.cc.

References python.hepunit::electron_Compton_length, python.hepunit::fine_structure_const, G4Log(), G4Material::GetElectronDensity(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4VEmModel::GetParticleChangeForLoss(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4VEmModel::IsMaster(), eplot::material, n, python.hepunit::pi, and SetParticle().

124 {
125  if(!monopole) { SetParticle(p); }
126  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
127  if(IsMaster()) {
128  if(!dedx0) { dedx0 = new std::vector<G4double>; }
129  G4ProductionCutsTable* theCoupleTable=
131  G4int numOfCouples = theCoupleTable->GetTableSize();
132  G4int n = dedx0->size();
133  if(n < numOfCouples) { dedx0->resize(numOfCouples); }
134 
135  // initialise vector
136  for(G4int i=0; i<numOfCouples; ++i) {
137 
138  const G4Material* material =
139  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
140  G4double eDensity = material->GetElectronDensity();
141  G4double vF = electron_Compton_length*pow(3*pi*pi*eDensity,0.3333333333);
142  (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
143  (G4Log(2*vF/fine_structure_const) - 0.5)/vF;
144  }
145  }
146 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
void SetParticle(const G4ParticleDefinition *p)
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
G4double GetElectronDensity() const
Definition: G4Material.hh:215
const G4int n
G4double G4Log(G4double x)
Definition: G4Log.hh:227
electron_Compton_length
Definition: hepunit.py:289
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const
G4double G4mplIonisationWithDeltaModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
)
protectedvirtual

Reimplemented from G4VEmModel.

Definition at line 371 of file G4mplIonisationWithDeltaModel.cc.

References python.hepunit::electron_mass_c2.

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), and SampleSecondaries().

373 {
374  G4double tau = kinEnergy/mass;
375  return 2.0*electron_mass_c2*tau*(tau + 2.);
376 }
float electron_mass_c2
Definition: hepunit.py:274
double G4double
Definition: G4Types.hh:76
G4double G4mplIonisationWithDeltaModel::SampleFluctuations ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmax,
G4double  length,
G4double  meanLoss 
)
virtual

Implements G4VEmFluctuationModel.

Definition at line 322 of file G4mplIonisationWithDeltaModel.cc.

References Dispersion(), G4UniformRand, G4MaterialCutsCouple::GetMaterial(), G4INCL::DeJongSpin::shoot(), and test::x.

328 {
329  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
330  G4double loss = meanLoss;
331  siga = sqrt(siga);
332  G4double twomeanLoss = meanLoss + meanLoss;
333 
334  if(twomeanLoss < siga) {
335  G4double x;
336  do {
337  loss = twomeanLoss*G4UniformRand();
338  x = (loss - meanLoss)/siga;
339  } while (1.0 - 0.5*x*x < G4UniformRand());
340  } else {
341  do {
342  loss = G4RandGauss::shoot(meanLoss,siga);
343  } while (0.0 > loss || loss > twomeanLoss);
344  }
345  return loss;
346 }
ThreeVector shoot(const G4int Ap, const G4int Af)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)
#define G4UniformRand()
Definition: Randomize.hh:87
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const
void G4mplIonisationWithDeltaModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 264 of file G4mplIonisationWithDeltaModel.cc.

References python.hepunit::electron_mass_c2, G4UniformRand, G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), MaxSecondaryEnergy(), G4INCL::Math::min(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), python.hepunit::twopi, and CLHEP::Hep3Vector::unit().

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

Definition at line 107 of file G4mplIonisationWithDeltaModel.cc.

References G4ParticleDefinition::GetPDGMass(), G4VEmModel::HighEnergyLimit(), G4VEmModel::LowEnergyLimit(), G4INCL::Math::max(), G4INCL::Math::min(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), Initialise(), and G4mplIonisation::InitialiseEnergyLossProcess().

108 {
109  monopole = p;
110  mass = monopole->GetPDGMass();
111  G4double emin =
112  std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
113  G4double emax =
114  std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
115  SetLowEnergyLimit(emin);
116  SetHighEnergyLimit(emax);
117 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
const char * p
Definition: xmltok.h:285
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
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:690

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