Geant4 11.1.1
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")
 
 ~G4mplIonisationModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss) override
 
G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override
 
void SetParticle (const G4ParticleDefinition *p)
 
G4mplIonisationModeloperator= (const G4mplIonisationModel &right)=delete
 
 G4mplIonisationModel (const G4mplIonisationModel &)=delete
 
- 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
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)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 *)
 
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 *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss)=0
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const 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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 57 of file G4mplIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4mplIonisationModel() [1/2]

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

Definition at line 69 of file G4mplIonisationModel.cc.

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

◆ ~G4mplIonisationModel()

G4mplIonisationModel::~G4mplIonisationModel ( )
override

Definition at line 88 of file G4mplIonisationModel.cc.

89{
90 if(IsMaster()) { delete dedx0; }
91}
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4mplIonisationModel() [2/2]

G4mplIonisationModel::G4mplIonisationModel ( const G4mplIonisationModel )
delete

Member Function Documentation

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 139 of file G4mplIonisationModel.cc.

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

◆ Dispersion()

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

Implements G4VEmFluctuationModel.

Definition at line 260 of file G4mplIonisationModel.cc.

265{
266 G4double siga = 0.0;
267 G4double tau = dp->GetKineticEnergy()/mass;
268 if(tau > 0.0) {
269 const G4double beta = dp->GetBeta();
270 siga = (tmax/(beta*beta) - 0.5*tcut) * twopi_mc2_rcl2 * length
271 * material->GetElectronDensity() * chargeSquare;
272 }
273 return siga;
274}
G4double GetKineticEnergy() const
G4double GetBeta() const
G4double GetElectronDensity() const
Definition: G4Material.hh:212

Referenced by SampleFluctuations().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 109 of file G4mplIonisationModel.cc.

111{
112 if(nullptr == monopole) { SetParticle(p); }
113 if(nullptr == fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
114 if(IsMaster()) {
115 if(nullptr == dedx0) { dedx0 = new std::vector<G4double>; }
116 G4ProductionCutsTable* theCoupleTable=
118 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
119 G4int n = (G4int)dedx0->size();
120 if(n < numOfCouples) { dedx0->resize(numOfCouples); }
121
122 G4Pow* g4calc = G4Pow::GetInstance();
123
124 // initialise vector assuming low conductivity
125 for(G4int i=0; i<numOfCouples; ++i) {
126
127 const G4Material* material =
128 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
129 G4double eDensity = material->GetElectronDensity();
130 G4double vF2 = 2*electron_Compton_length*g4calc->A13(3.*pi*pi*eDensity);
131 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
132 (G4Log(vF2/fine_structure_const) - 0.5)/vF2;
133 }
134 }
135}
const G4Material * GetMaterial() const
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:116
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

◆ operator=()

G4mplIonisationModel & G4mplIonisationModel::operator= ( const G4mplIonisationModel right)
delete

◆ SampleFluctuations()

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

Implements G4VEmFluctuationModel.

Definition at line 229 of file G4mplIonisationModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 220 of file G4mplIonisationModel.cc.

225{}

◆ SetParticle()

void G4mplIonisationModel::SetParticle ( const G4ParticleDefinition p)

Definition at line 95 of file G4mplIonisationModel.cc.

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

Referenced by ComputeDEDXPerVolume(), and Initialise().


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