Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4mplIonisationModel Class Reference

#include <G4mplIonisationModel.hh>

Inheritance diagram for G4mplIonisationModel:
Collaboration diagram for G4mplIonisationModel:

Public Member Functions

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

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
 
G4bool lossFlucFlag
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 58 of file G4mplIonisationModel.hh.

Constructor & Destructor Documentation

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

Definition at line 71 of file G4mplIonisationModel.cc.

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

Here is the call graph for this function:

G4mplIonisationModel::~G4mplIonisationModel ( )
virtual

Definition at line 93 of file G4mplIonisationModel.cc.

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

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 142 of file G4mplIonisationModel.cc.

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

Here is the call graph for this function:

G4double G4mplIonisationModel::Dispersion ( const G4Material material,
const G4DynamicParticle dp,
G4double  tmax,
G4double  length 
)
overridevirtual

Implements G4VEmFluctuationModel.

Definition at line 262 of file G4mplIonisationModel.cc.

266 {
267  G4double siga = 0.0;
268  G4double tau = dp->GetKineticEnergy()/mass;
269  if(tau > 0.0) {
270  G4double electronDensity = material->GetElectronDensity();
271  G4double gam = tau + 1.0;
272  G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
273  siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
274  * electronDensity * chargeSquare;
275  }
276  return siga;
277 }
G4double GetKineticEnergy() const
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 G4mplIonisationModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 114 of file G4mplIonisationModel.cc.

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

Here is the call graph for this function:

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

Implements G4VEmFluctuationModel.

Definition at line 232 of file G4mplIonisationModel.cc.

238 {
239  G4double siga = Dispersion(couple->GetMaterial(),dp,tmax,length);
240  G4double loss = meanLoss;
241  siga = sqrt(siga);
242  G4double twomeanLoss = meanLoss + meanLoss;
243 
244  if(twomeanLoss < siga) {
245  G4double x;
246  do {
247  loss = twomeanLoss*G4UniformRand();
248  x = (loss - meanLoss)/siga;
249  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
250  } while (1.0 - 0.5*x*x < G4UniformRand());
251  } else {
252  do {
253  loss = G4RandGauss::shoot(meanLoss,siga);
254  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
255  } while (0.0 > loss || loss > twomeanLoss);
256  }
257  return loss;
258 }
ThreeVector shoot(const G4int Ap, const G4int Af)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length) override
tuple x
Definition: test.py:50
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 223 of file G4mplIonisationModel.cc.

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

Definition at line 100 of file G4mplIonisationModel.cc.

101 {
102  monopole = p;
103  mass = monopole->GetPDGMass();
104  G4double emin =
105  std::min(LowEnergyLimit(),0.1*mass*(1./sqrt(1. - betalow*betalow) - 1.));
106  G4double emax =
107  std::max(HighEnergyLimit(),10.*mass*(1./sqrt(1. - beta2lim) - 1.));
108  SetLowEnergyLimit(emin);
109  SetHighEnergyLimit(emax);
110 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
const char * p
Definition: xmltok.h:285
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
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:739

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: