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

#include <G4mplIonisationWithDeltaModel.hh>

+ Inheritance diagram for G4mplIonisationWithDeltaModel:

Public Member Functions

 G4mplIonisationWithDeltaModel (G4double mCharge, const G4String &nam="mplIonisationWithDelta")
 
virtual ~G4mplIonisationWithDeltaModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss) override
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) 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 *, 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
 

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) override
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Additional Inherited Members

- 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 G4mplIonisationWithDeltaModel.hh.

Constructor & Destructor Documentation

◆ G4mplIonisationWithDeltaModel()

G4mplIonisationWithDeltaModel::G4mplIonisationWithDeltaModel ( G4double  mCharge,
const G4String nam = "mplIonisationWithDelta" 
)
explicit

Definition at line 73 of file G4mplIonisationWithDeltaModel.cc.

76 magCharge(mCharge),
77 twoln10(std::log(100.0)),
78 betalow(0.01),
79 betalim(0.1),
80 beta2lim(betalim*betalim),
81 bg2lim(beta2lim*(1.0 + beta2lim))
82{
83 nmpl = G4lrint(std::abs(magCharge) * 2 * fine_structure_const);
84 if(nmpl > 6) { nmpl = 6; }
85 else if(nmpl < 1) { nmpl = 1; }
86 pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
87 chargeSquare = magCharge * magCharge;
88 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
89 fParticleChange = nullptr;
90 theElectron = G4Electron::Electron();
91 G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
92 << magCharge/eplus << G4endl;
93 monopole = nullptr;
94 mass = 0.0;
95}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4double pi
int G4lrint(double ad)
Definition: templates.hh:134

◆ ~G4mplIonisationWithDeltaModel()

G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel ( )
virtual

Definition at line 99 of file G4mplIonisationWithDeltaModel.cc.

100{
101 if(IsMaster()) { delete dedx0; }
102}
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4mplIonisationWithDeltaModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 253 of file G4mplIonisationWithDeltaModel.cc.

259{
260 G4double cross =
261 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
262 return cross;
263}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Definition at line 235 of file G4mplIonisationWithDeltaModel.cc.

240{
241 if(!monopole) { SetParticle(p); }
242 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
243 G4double maxEnergy = std::min(tmax, maxKinEnergy);
244 G4double cutEnergy = std::max(LowEnergyLimit(), cut);
245 G4double cross = (cutEnergy < maxEnergy)
246 ? (0.5/cutEnergy - 0.5/maxEnergy)*pi_hbarc2_over_mc2 * nmpl * nmpl : 0.0;
247 return cross;
248}
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
void SetParticle(const G4ParticleDefinition *p)

Referenced by ComputeCrossSectionPerAtom().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 160 of file G4mplIonisationWithDeltaModel.cc.

164{
165 if(!monopole) { SetParticle(p); }
166 G4double tmax = MaxSecondaryEnergy(p,kineticEnergy);
167 G4double cutEnergy = std::min(tmax, maxEnergy);
168 cutEnergy = std::max(LowEnergyLimit(), cutEnergy);
169 G4double tau = kineticEnergy / mass;
170 G4double gam = tau + 1.0;
171 G4double bg2 = tau * (tau + 2.0);
172 G4double beta2 = bg2 / (gam * gam);
173 G4double beta = sqrt(beta2);
174
175 // low-energy asymptotic formula
176 G4double dedx = (*dedx0)[CurrentCouple()->GetIndex()]*beta;
177
178 // above asymptotic
179 if(beta > betalow) {
180
181 // high energy
182 if(beta >= betalim) {
183 dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
184
185 } else {
186 G4double dedx1 = (*dedx0)[CurrentCouple()->GetIndex()]*betalow;
187 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
188
189 // extrapolation between two formula
190 G4double kapa2 = beta - betalow;
191 G4double kapa1 = betalim - beta;
192 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
193 }
194 }
195 return dedx;
196}
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486

◆ Dispersion()

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

Implements G4VEmFluctuationModel.

Definition at line 358 of file G4mplIonisationWithDeltaModel.cc.

363{
364 G4double siga = 0.0;
365 G4double tau = dp->GetKineticEnergy()/mass;
366 if(tau > 0.0) {
367 const G4double beta = dp->GetBeta();
368 siga = (tmax/(beta*beta) - 0.5*tcut) * twopi_mc2_rcl2 * length
369 * material->GetElectronDensity() * chargeSquare;
370 }
371 return siga;
372}
G4double GetKineticEnergy() const
G4double GetBeta() const
G4double GetElectronDensity() const
Definition: G4Material.hh:212

Referenced by SampleFluctuations().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 121 of file G4mplIonisationWithDeltaModel.cc.

123{
124 if(!monopole) { SetParticle(p); }
125 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
126 if(IsMaster()) {
127 if(!dedx0) { dedx0 = new std::vector<G4double>; }
128 G4ProductionCutsTable* theCoupleTable=
130 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
131 G4int n = (G4int)dedx0->size();
132 if(n < numOfCouples) { dedx0->resize(numOfCouples); }
133 G4Pow* g4calc = G4Pow::GetInstance();
134
135 // initialise vector assuming low conductivity
136 for(G4int i=0; i<numOfCouples; ++i) {
137
138 const G4Material* material =
139 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
140 G4double eDensity = material->GetElectronDensity();
141 G4double vF2 = 2*electron_Compton_length*g4calc->A13(3.*pi*pi*eDensity);
142 (*dedx0)[i] = pi_hbarc2_over_mc2*eDensity*nmpl*nmpl*
143 (G4Log(vF2/fine_structure_const) - 0.5)/vF2;
144 }
145 }
146}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
int G4int
Definition: G4Types.hh:85
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

◆ MaxSecondaryEnergy()

G4double G4mplIonisationWithDeltaModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
)
overrideprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 377 of file G4mplIonisationWithDeltaModel.cc.

379{
380 G4double tau = kinEnergy/mass;
381 return 2.0*electron_mass_c2*tau*(tau + 2.);
382}

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), and SampleSecondaries().

◆ MinEnergyCut()

G4double G4mplIonisationWithDeltaModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple couple 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 151 of file G4mplIonisationWithDeltaModel.cc.

153{
155}
G4double GetMeanExcitationEnergy() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221

◆ SampleFluctuations()

G4double G4mplIonisationWithDeltaModel::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 326 of file G4mplIonisationWithDeltaModel.cc.

333{
334 G4double siga = Dispersion(couple->GetMaterial(),dp,tcut,tmax,length);
335 G4double loss = meanLoss;
336 siga = std::sqrt(siga);
337 G4double twomeanLoss = meanLoss + meanLoss;
338
339 if(twomeanLoss < siga) {
340 G4double x;
341 do {
342 loss = twomeanLoss*G4UniformRand();
343 x = (loss - meanLoss)/siga;
344 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
345 } while (1.0 - 0.5*x*x < G4UniformRand());
346 } else {
347 do {
348 loss = G4RandGauss::shoot(meanLoss,siga);
349 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
350 } while (0.0 > loss || loss > twomeanLoss);
351 }
352 return loss;
353}
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length) override

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 268 of file G4mplIonisationWithDeltaModel.cc.

273{
274 G4double kineticEnergy = dp->GetKineticEnergy();
275 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
276
277 G4double maxKinEnergy = std::min(maxEnergy,tmax);
278 if(minKinEnergy >= maxKinEnergy) { return; }
279
280 //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
281 // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
282 // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
283
284 G4double totEnergy = kineticEnergy + mass;
285 G4double etot2 = totEnergy*totEnergy;
286 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
287
288 // sampling without nuclear size effect
290 G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
291 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
292
293 // delta-electron is produced
294 G4double totMomentum = totEnergy*sqrt(beta2);
295 G4double deltaMomentum =
296 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
297 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
298 (deltaMomentum * totMomentum);
299 cost = std::min(cost, 1.0);
300
301 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
302
303 G4double phi = twopi * G4UniformRand() ;
304
305 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
306 G4ThreeVector direction = dp->GetMomentumDirection();
307 deltaDirection.rotateUz(direction);
308
309 // create G4DynamicParticle object for delta ray
310 G4DynamicParticle* delta =
311 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
312
313 vdp->push_back(delta);
314
315 // Change kinematics of primary particle
316 kineticEnergy -= deltaKinEnergy;
317 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
318 finalP = finalP.unit();
319
320 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
321 fParticleChange->SetProposedMomentumDirection(finalP);
322}
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)

◆ SetParticle()

void G4mplIonisationWithDeltaModel::SetParticle ( const G4ParticleDefinition p)

Definition at line 106 of file G4mplIonisationWithDeltaModel.cc.

107{
108 monopole = p;
109 mass = monopole->GetPDGMass();
110 G4double emin =
111 std::min(LowEnergyLimit(),0.1*mass*(1./sqrt(1. - betalow*betalow) - 1.));
112 G4double emax =
113 std::max(HighEnergyLimit(),10*mass*(1./sqrt(1. - beta2lim) - 1.));
114 SetLowEnergyLimit(emin);
115 SetHighEnergyLimit(emax);
116}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), Initialise(), and G4mplIonisation::InitialiseEnergyLossProcess().


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