Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4mplIonisationModel Class Reference

#include <G4mplIonisationModel.hh>

+ Inheritance diagram for G4mplIonisationModel:

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 G4Material *, 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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
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 G4double SampleFluctuations (const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length, G4double &meanLoss)=0
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length)=0
 
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
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 57 of file G4mplIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4mplIonisationModel()

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

Definition at line 67 of file G4mplIonisationModel.cc.

69 magCharge(mCharge),
70 twoln10(log(100.0)),
71 betalow(0.01),
72 betalim(0.1),
73 beta2lim(betalim*betalim),
74 bg2lim(beta2lim*(1.0 + beta2lim))
75{
76 nmpl = G4int(abs(magCharge) * 2 * fine_structure_const + 0.5);
77 if(nmpl > 6) { nmpl = 6; }
78 else if(nmpl < 1) { nmpl = 1; }
79 pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
80 chargeSquare = magCharge * magCharge;
81 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
82 fParticleChange = 0;
83 monopole = 0;
84 mass = 0.0;
85}
int G4int
Definition: G4Types.hh:66
const G4double pi

◆ ~G4mplIonisationModel()

G4mplIonisationModel::~G4mplIonisationModel ( )
virtual

Definition at line 89 of file G4mplIonisationModel.cc.

90{}

Member Function Documentation

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 117 of file G4mplIonisationModel.cc.

121{
122 if(!monopole) { SetParticle(p); }
123 G4double tau = kineticEnergy / mass;
124 G4double gam = tau + 1.0;
125 G4double bg2 = tau * (tau + 2.0);
126 G4double beta2 = bg2 / (gam * gam);
127 G4double beta = sqrt(beta2);
128
129 // low-energy asymptotic formula
130 G4double dedx = dedxlim*beta*material->GetDensity();
131
132 // above asymptotic
133 if(beta > betalow) {
134
135 // high energy
136 if(beta >= betalim) {
137 dedx = ComputeDEDXAhlen(material, bg2);
138
139 } else {
140
141 G4double dedx1 = dedxlim*betalow*material->GetDensity();
142 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim);
143
144 // extrapolation between two formula
145 G4double kapa2 = beta - betalow;
146 G4double kapa1 = betalim - beta;
147 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
148 }
149 }
150 return dedx;
151}
double G4double
Definition: G4Types.hh:64
G4double GetDensity() const
Definition: G4Material.hh:179
void SetParticle(const G4ParticleDefinition *p)

◆ Dispersion()

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

Implements G4VEmFluctuationModel.

Definition at line 233 of file G4mplIonisationModel.cc.

237{
238 G4double siga = 0.0;
239 G4double tau = dp->GetKineticEnergy()/mass;
240 if(tau > 0.0) {
241 G4double electronDensity = material->GetElectronDensity();
242 G4double gam = tau + 1.0;
243 G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
244 siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
245 * electronDensity * chargeSquare;
246 }
247 return siga;
248}
G4double GetKineticEnergy() const
G4double GetElectronDensity() const
Definition: G4Material.hh:216

Referenced by SampleFluctuations().

◆ Initialise()

void G4mplIonisationModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 108 of file G4mplIonisationModel.cc.

110{
111 if(!monopole) { SetParticle(p); }
112 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
113}
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95

◆ SampleFluctuations()

G4double G4mplIonisationModel::SampleFluctuations ( const G4Material material,
const G4DynamicParticle dp,
G4double tmax,
G4double length,
G4double meanLoss 
)
virtual

Implements G4VEmFluctuationModel.

Definition at line 205 of file G4mplIonisationModel.cc.

211{
212 G4double siga = Dispersion(material,dp,tmax,length);
213 G4double loss = meanLoss;
214 siga = sqrt(siga);
215 G4double twomeanLoss = meanLoss + meanLoss;
216
217 if(twomeanLoss < siga) {
218 G4double x;
219 do {
220 loss = twomeanLoss*G4UniformRand();
221 x = (loss - meanLoss)/siga;
222 } while (1.0 - 0.5*x*x < G4UniformRand());
223 } else {
224 do {
225 loss = G4RandGauss::shoot(meanLoss,siga);
226 } while (0.0 > loss || loss > twomeanLoss);
227 }
228 return loss;
229}
#define G4UniformRand()
Definition: Randomize.hh:53
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length)

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 196 of file G4mplIonisationModel.cc.

201{}

◆ SetParticle()

void G4mplIonisationModel::SetParticle ( const G4ParticleDefinition p)

Definition at line 94 of file G4mplIonisationModel.cc.

95{
96 monopole = p;
97 mass = monopole->GetPDGMass();
98 G4double emin =
99 std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
100 G4double emax =
101 std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
102 SetLowEnergyLimit(emin);
103 SetHighEnergyLimit(emax);
104}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

Referenced by ComputeDEDXPerVolume(), and Initialise().


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