Geant4 10.7.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 &) 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 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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
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
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 
- Public Member Functions inherited from G4VEmFluctuationModel
 G4VEmFluctuationModel (const G4String &nam)
 
virtual ~G4VEmFluctuationModel ()
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, 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)
 
const G4StringGetName () const
 
G4VEmFluctuationModeloperator= (const G4VEmFluctuationModel &right)=delete
 
 G4VEmFluctuationModel (const G4VEmFluctuationModel &)=delete
 

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 G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 57 of file G4mplIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4mplIonisationModel()

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}
int G4int
Definition: G4Types.hh:85
const G4double pi

◆ ~G4mplIonisationModel()

G4mplIonisationModel::~G4mplIonisationModel ( )
virtual

Definition at line 93 of file G4mplIonisationModel.cc.

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

Member Function Documentation

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 144 of file G4mplIonisationModel.cc.

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

◆ Dispersion()

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

Implements G4VEmFluctuationModel.

Definition at line 264 of file G4mplIonisationModel.cc.

268{
269 G4double siga = 0.0;
270 G4double tau = dp->GetKineticEnergy()/mass;
271 if(tau > 0.0) {
272 G4double electronDensity = material->GetElectronDensity();
273 G4double gam = tau + 1.0;
274 G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
275 siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
276 * electronDensity * chargeSquare;
277 }
278 return siga;
279}
G4double GetKineticEnergy() const
G4double GetElectronDensity() const
Definition: G4Material.hh:215

Referenced by SampleFluctuations().

◆ Initialise()

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 G4Pow* g4calc = G4Pow::GetInstance();
128
129 // initialise vector assuming low conductivity
130 for(G4int i=0; i<numOfCouples; ++i) {
131
132 const G4Material* material =
133 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
134 G4double eDensity = material->GetElectronDensity();
135 G4double vF2 = 2*electron_Compton_length*g4calc->A13(3.*pi*pi*eDensity);
136 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
137 (G4Log(vF2/fine_structure_const) - 0.5)/vF2;
138 }
139 }
140}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
const G4Material * GetMaterial() const
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

◆ SampleFluctuations()

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

Implements G4VEmFluctuationModel.

Definition at line 234 of file G4mplIonisationModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 225 of file G4mplIonisationModel.cc.

230{}

◆ SetParticle()

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}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

Referenced by ComputeDEDXPerVolume(), and Initialise().


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