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

#include <G4mplIonisationModel.hh>

Inheritance diagram for G4mplIonisationModel:
G4VEmModel G4VEmFluctuationModel

Public Member Functions

 G4mplIonisationModel (G4double mCharge, const G4String &nam="mplIonisation")
 
virtual ~G4mplIonisationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- 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 G4mplIonisationModel.hh.

Constructor & Destructor Documentation

G4mplIonisationModel::G4mplIonisationModel ( G4double  mCharge,
const G4String nam = "mplIonisation" 
)

Definition at line 72 of file G4mplIonisationModel.cc.

References python.hepunit::cm2, python.hepunit::electron_mass_c2, python.hepunit::fine_structure_const, g(), python.hepunit::GeV, python.hepunit::hbarc, and python.hepunit::pi.

74  magCharge(mCharge),
75  twoln10(log(100.0)),
76  betalow(0.01),
77  betalim(0.1),
78  beta2lim(betalim*betalim),
79  bg2lim(beta2lim*(1.0 + beta2lim))
80 {
81  nmpl = G4int(abs(magCharge) * 2 * fine_structure_const + 0.5);
82  if(nmpl > 6) { nmpl = 6; }
83  else if(nmpl < 1) { nmpl = 1; }
84  pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
85  chargeSquare = magCharge * magCharge;
86  dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
87  fParticleChange = 0;
88  monopole = 0;
89  mass = 0.0;
90 }
int G4int
Definition: G4Types.hh:78
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
float electron_mass_c2
Definition: hepunit.py:274
G4VEmFluctuationModel(const G4String &nam)
G4mplIonisationModel::~G4mplIonisationModel ( )
virtual

Definition at line 94 of file G4mplIonisationModel.cc.

References G4VEmModel::IsMaster().

95 {
96  if(IsMaster()) { delete dedx0; }
97 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:676

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 143 of file G4mplIonisationModel.cc.

References G4VEmModel::CurrentCouple(), G4InuclParticleNames::gam, G4MaterialCutsCouple::GetIndex(), and SetParticle().

147 {
148  if(!monopole) { SetParticle(p); }
149  G4double tau = kineticEnergy / mass;
150  G4double gam = tau + 1.0;
151  G4double bg2 = tau * (tau + 2.0);
152  G4double beta2 = bg2 / (gam * gam);
153  G4double beta = sqrt(beta2);
154 
155  // low-energy asymptotic formula
156  //G4double dedx = dedxlim*beta*material->GetDensity();
157  G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
158 
159  // above asymptotic
160  if(beta > betalow) {
161 
162  // high energy
163  if(beta >= betalim) {
164  dedx = ComputeDEDXAhlen(material, bg2);
165 
166  } else {
167 
168  //G4double dedx1 = dedxlim*betalow*material->GetDensity();
169  G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
170  G4double dedx2 = ComputeDEDXAhlen(material, bg2lim);
171 
172  // extrapolation between two formula
173  G4double kapa2 = beta - betalow;
174  G4double kapa1 = betalim - beta;
175  dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
176  }
177  }
178  return dedx;
179 }
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
void SetParticle(const G4ParticleDefinition *p)
double G4double
Definition: G4Types.hh:76
G4double G4mplIonisationModel::Dispersion ( const G4Material material,
const G4DynamicParticle dp,
G4double  tmax,
G4double  length 
)
virtual

Implements G4VEmFluctuationModel.

Definition at line 261 of file G4mplIonisationModel.cc.

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

Referenced by SampleFluctuations().

265 {
266  G4double siga = 0.0;
267  G4double tau = dp->GetKineticEnergy()/mass;
268  if(tau > 0.0) {
269  G4double electronDensity = material->GetElectronDensity();
270  G4double gam = tau + 1.0;
271  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
272  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
273  * electronDensity * chargeSquare;
274  }
275  return siga;
276 }
G4double GetKineticEnergy() const
G4double GetElectronDensity() const
Definition: G4Material.hh:215
double G4double
Definition: G4Types.hh:76
void G4mplIonisationModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 115 of file G4mplIonisationModel.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().

117 {
118  if(!monopole) { SetParticle(p); }
119  if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
120  if(IsMaster()) {
121  if(!dedx0) { dedx0 = new std::vector<G4double>; }
122  G4ProductionCutsTable* theCoupleTable=
124  G4int numOfCouples = theCoupleTable->GetTableSize();
125  G4int n = dedx0->size();
126  if(n < numOfCouples) { dedx0->resize(numOfCouples); }
127 
128  // initialise vector
129  for(G4int i=0; i<numOfCouples; ++i) {
130 
131  const G4Material* material =
132  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
133  G4double eDensity = material->GetElectronDensity();
134  G4double vF = electron_Compton_length*pow(3*pi*pi*eDensity,0.3333333333);
135  (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
136  (G4Log(2*vF/fine_structure_const) - 0.5)/vF;
137  }
138  }
139 }
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
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
void SetParticle(const G4ParticleDefinition *p)
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 G4mplIonisationModel::SampleFluctuations ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmax,
G4double  length,
G4double  meanLoss 
)
virtual

Implements G4VEmFluctuationModel.

Definition at line 233 of file G4mplIonisationModel.cc.

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

239 {
240  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
241  G4double loss = meanLoss;
242  siga = sqrt(siga);
243  G4double twomeanLoss = meanLoss + meanLoss;
244 
245  if(twomeanLoss < siga) {
246  G4double x;
247  do {
248  loss = twomeanLoss*G4UniformRand();
249  x = (loss - meanLoss)/siga;
250  } while (1.0 - 0.5*x*x < G4UniformRand());
251  } else {
252  do {
253  loss = G4RandGauss::shoot(meanLoss,siga);
254  } while (0.0 > loss || loss > twomeanLoss);
255  }
256  return loss;
257 }
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 G4mplIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 224 of file G4mplIonisationModel.cc.

229 {}
void G4mplIonisationModel::SetParticle ( const G4ParticleDefinition p)

Definition at line 101 of file G4mplIonisationModel.cc.

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

Referenced by ComputeDEDXPerVolume(), and Initialise().

102 {
103  monopole = p;
104  mass = monopole->GetPDGMass();
105  G4double emin =
106  std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
107  G4double emax =
108  std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
109  SetLowEnergyLimit(emin);
110  SetHighEnergyLimit(emax);
111 }
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: