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

#include <G4AtimaEnergyLossModel.hh>

Inheritance diagram for G4AtimaEnergyLossModel:

Public Member Functions

 G4AtimaEnergyLossModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="Atima")
 ~G4AtimaEnergyLossModel () override
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override
G4double GetChargeSquareRatio (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
void CorrectionsAlongStep (const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &) override
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
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 void StartTracking (G4Track *)
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
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 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 SetDeexcitationFlag (G4bool val)
void SetForceBuildTable (G4bool val)
void SetFluctuationFlag (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)
void SetLPMFlag (G4bool)
void SetMasterThread (G4bool)
G4VEmModeloperator= (const G4VEmModel &right)=delete
 G4VEmModel (const G4VEmModel &)=delete

Protected Member Functions

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) override
G4double GetChargeSquareRatio () const
void SetChargeSquareRatio (G4double val)
Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
G4ParticleChangeForGammaGetParticleChangeForGamma ()
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
std::size_t currentCoupleIndex = 0
std::size_t basedCoupleIndex = 0
G4bool lossFlucFlag = true

Detailed Description

Definition at line 57 of file G4AtimaEnergyLossModel.hh.

Constructor & Destructor Documentation

◆ G4AtimaEnergyLossModel()

G4AtimaEnergyLossModel::G4AtimaEnergyLossModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "Atima" )
explicit

Definition at line 74 of file G4AtimaEnergyLossModel.cc.

76 : G4VEmModel(nam),
77 particle(nullptr)
78{
79 g4calc = G4Pow::GetInstance();
80 fParticleChange = nullptr;
81 theElectron = G4Electron::Electron();
84 SetLowEnergyLimit(2.0*MeV);
85 MLN10 = 2.30258509299;
86 atomic_mass_unit = 931.4940954; // MeV/c^2
87 dedx_constant = 0.3070749187; //4*pi*Na*me*c^2*r_e^2 //MeV cm^2
88 electron_mass = 0.510998928; // MeV/c^2
89 fine_structure = 1.0/137.035999139;
90 domega2dx_constant = dedx_constant*electron_mass; //4*pi*Na*me*c^2*r_e^2 //MeV^2 cm^2
91 if(tableE[0] == 0.0) {
92 G4double logmin = 0.;
93 G4double logmax = 5.;
94 stepE = (logmax-logmin)/199.;
95 for(G4int i=0; i<200; ++i){
96 tableE[i] = G4Exp(MLN10*(logmin + ((G4double)i)*stepE));
97 }
98 }
99}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4NistManager * Instance()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

Referenced by ~G4AtimaEnergyLossModel().

◆ ~G4AtimaEnergyLossModel()

G4AtimaEnergyLossModel::~G4AtimaEnergyLossModel ( )
overridedefault

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 222 of file G4AtimaEnergyLossModel.cc.

228{
230 (p,kineticEnergy,cutEnergy,maxEnergy);
231 return cross;
232}
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4AtimaEnergyLossModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )

Definition at line 191 of file G4AtimaEnergyLossModel.cc.

195{
196 G4double cross = 0.0;
197 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
198 G4double maxEnergy = std::min(tmax,maxKinEnergy);
199 if(cutEnergy < maxEnergy) {
200
201 G4double totEnergy = kineticEnergy + mass;
202 G4double energy2 = totEnergy*totEnergy;
203 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
204
205 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
206 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
207
208 // +term for spin=1/2 particle
209 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
210
211 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
212 }
213
214 // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy
215 // << " tmax= " << tmax << " cross= " << cross << G4endl;
216
217 return cross;
218}
G4double G4Log(G4double x)
Definition G4Log.hh:169
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

G4double G4AtimaEnergyLossModel::ComputeDEDXPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double  )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 251 of file G4AtimaEnergyLossModel.cc.

255{
256 //Call to ATIMA Stopping Power function
257 G4double zt = material->GetIonisation()->GetZeffective();
258 zt = std::min(zt,93.);
259 G4double at = nist->GetAtomicMassAmu(G4lrint(zt));
260 G4double dedx = StoppingPower(p->GetPDGMass(), p->GetPDGCharge(),
261 kineticEnergy/(MeV), at, zt) *material->GetDensity()/(g/cm3);
262 dedx = std::max(dedx, 0.0);
263
264 // G4cout << "E(MeV)= " << kineticEnergy/MeV << " dedx= " << dedx
265 // << " " << material->GetName() << G4endl;
266 return dedx;
267}
G4double GetZeffective() const
G4double GetDensity() const
G4IonisParamMat * GetIonisation() const
int G4lrint(double ad)
Definition templates.hh:134

Referenced by CorrectionsAlongStep().

◆ CorrectionsAlongStep()

void G4AtimaEnergyLossModel::CorrectionsAlongStep ( const G4Material * mat,
const G4ParticleDefinition * p,
const G4double kinEnergy,
const G4double cutEnergy,
const G4double & length,
G4double & eloss )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 271 of file G4AtimaEnergyLossModel.cc.

277{
278 if (isIon) {
279 eloss = ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy)*length/(cm);
280 }
281}
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double) override

◆ CrossSectionPerVolume()

G4double G4AtimaEnergyLossModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 236 of file G4AtimaEnergyLossModel.cc.

242{
243 G4double eDensity = material->GetElectronDensity();
245 (p,kineticEnergy,cutEnergy,maxEnergy);
246 return cross;
247}
G4double GetElectronDensity() const

◆ GetChargeSquareRatio() [1/2]

G4double G4AtimaEnergyLossModel::GetChargeSquareRatio ( ) const
inlineprotected

Definition at line 209 of file G4AtimaEnergyLossModel.hh.

210{
211 return chargeSquare;
212}

◆ GetChargeSquareRatio() [2/2]

G4double G4AtimaEnergyLossModel::GetChargeSquareRatio ( const G4ParticleDefinition * p,
const G4Material * mat,
G4double kineticEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 130 of file G4AtimaEnergyLossModel.cc.

133{
134 // this method is called only for ions
135 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
136 corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
137 return corrFactor;
138}

◆ GetParticleCharge()

G4double G4AtimaEnergyLossModel::GetParticleCharge ( const G4ParticleDefinition * p,
const G4Material * mat,
G4double kineticEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 142 of file G4AtimaEnergyLossModel.cc.

145{
146 //G4cout<<"G4AtimaEnergyLossModel::GetParticleCharge e= "<<kineticEnergy <<
147 // " q= " << corr->GetParticleCharge(p,mat,kineticEnergy) <<G4endl;
148 // this method is called only for ions, so no check if it is an ion
149 return corr->GetParticleCharge(p,mat,kineticEnergy);
150}

◆ Initialise()

void G4AtimaEnergyLossModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 107 of file G4AtimaEnergyLossModel.cc.

109{
110 SetGenericIon(p);
111 SetParticle(p);
112
113 //G4cout << "G4AtimaEnergyLossModel::Initialise for " << p->GetParticleName()
114 // << " isIon= " << isIon
115 // << G4endl;
116
117 // always false before the run
118 SetDeexcitationFlag(false);
119
120 if(nullptr == fParticleChange) {
121 fParticleChange = GetParticleChangeForLoss();
123 SetAngularDistribution(new G4DeltaAngle());
124 }
125 }
126}
G4VEmAngularDistribution * GetAngularDistribution()
void SetDeexcitationFlag(G4bool val)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ MaxSecondaryEnergy()

G4double G4AtimaEnergyLossModel::MaxSecondaryEnergy ( const G4ParticleDefinition * pd,
G4double kinEnergy )
overrideprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 400 of file G4AtimaEnergyLossModel.cc.

402{
403 // here particle type is checked for any method
404 SetParticle(pd);
405 G4double tau = kinEnergy/mass;
406 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
407 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
408 return std::min(tmax,tlimit);
409}

Referenced by ComputeCrossSectionPerElectron(), and SampleSecondaries().

◆ MinEnergyCut()

G4double G4AtimaEnergyLossModel::MinEnergyCut ( const G4ParticleDefinition * ,
const G4MaterialCutsCouple * couple )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 182 of file G4AtimaEnergyLossModel.cc.

184{
186}
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 285 of file G4AtimaEnergyLossModel.cc.

290{
291 G4double kineticEnergy = dp->GetKineticEnergy();
292 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
293
294 G4double maxKinEnergy = std::min(maxEnergy,tmax);
295 if(minKinEnergy >= maxKinEnergy) { return; }
296
297 //G4cout << "G4AtimaEnergyLossModel::SampleSecondaries Emin= " << minKinEnergy
298 // << " Emax= " << maxKinEnergy << G4endl;
299
300 G4double totEnergy = kineticEnergy + mass;
301 G4double etot2 = totEnergy*totEnergy;
302 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
303
304 G4double deltaKinEnergy, f;
305 G4double f1 = 0.0;
306 G4double fmax = 1.0;
307 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
308
309 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
310 G4double rndm[2];
311
312 // sampling without nuclear size effect
313 do {
314 rndmEngineMod->flatArray(2, rndm);
315 deltaKinEnergy = minKinEnergy*maxKinEnergy
316 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
317
318 f = 1.0 - beta2*deltaKinEnergy/tmax;
319 if( 0.0 < spin ) {
320 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
321 f += f1;
322 }
323
324 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
325 } while( fmax*rndm[1] > f);
326
327 // projectile formfactor - suppresion of high energy
328 // delta-electron production at high energy
329
330 G4double x = formfact*deltaKinEnergy*(deltaKinEnergy + 2*electron_mass_c2);
331 if(x > 1.e-6) {
332
333 G4double x10 = 1.0 + x;
334 G4double grej = 1.0/(x10*x10);
335 if( 0.0 < spin ) {
336 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
337 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
338 }
339 if(grej > 1.1) {
340 G4cout << "### G4AtimaEnergyLossModel WARNING: grej= " << grej
341 << " " << dp->GetDefinition()->GetParticleName()
342 << " Ekin(MeV)= " << kineticEnergy
343 << " delEkin(MeV)= " << deltaKinEnergy
344 << G4endl;
345 }
346 if(rndmEngineMod->flat() > grej) { return; }
347 }
348
349 G4ThreeVector deltaDirection;
350
352
353 const G4Material* mat = couple->GetMaterial();
355
356 deltaDirection =
357 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
358
359 } else {
360
361 G4double deltaMomentum =
362 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
363 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
364 (deltaMomentum * dp->GetTotalMomentum());
365 if(cost > 1.0) { cost = 1.0; }
366 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
367
368 G4double phi = twopi*rndmEngineMod->flat();
369
370 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
371 deltaDirection.rotateUz(dp->GetMomentumDirection());
372 }
373 /*
374 G4cout << "### G4AtimaEnergyLossModel "
375 << dp->GetDefinition()->GetParticleName()
376 << " Ekin(MeV)= " << kineticEnergy
377 << " delEkin(MeV)= " << deltaKinEnergy
378 << " tmin(MeV)= " << minKinEnergy
379 << " tmax(MeV)= " << maxKinEnergy
380 << " dir= " << dp->GetMomentumDirection()
381 << " dirDelta= " << deltaDirection
382 << G4endl;
383 */
384 // create G4DynamicParticle object for delta ray
385 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
386
387 vdp->push_back(delta);
388
389 // Change kinematics of primary particle
390 kineticEnergy -= deltaKinEnergy;
391 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
392 finalP = finalP.unit();
393
394 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
395 fParticleChange->SetProposedMomentumDirection(finalP);
396}
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *) const

◆ SetChargeSquareRatio()

void G4AtimaEnergyLossModel::SetChargeSquareRatio ( G4double val)
inlineprotected

Definition at line 216 of file G4AtimaEnergyLossModel.hh.

217{
218 chargeSquare = val;
219}

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