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

#include <G4GoudsmitSaundersonMscModel.hh>

Inheritance diagram for G4GoudsmitSaundersonMscModel:

Public Member Functions

 G4GoudsmitSaundersonMscModel (const G4String &nam="GoudsmitSaunderson")
 ~G4GoudsmitSaundersonMscModel () override
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal (const G4ParticleDefinition *p, G4VEmModel *masterModel) override
G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety) override
G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &currentMinimalStep) override
G4double ComputeGeomPathLength (G4double truePathLength) override
G4double ComputeTrueStepLength (G4double geomStepLength) override
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
void StartTracking (G4Track *) override
void SampleMSC ()
G4double GetTransportMeanFreePath (const G4ParticleDefinition *, G4double)
void SetOptionPWACorrection (G4bool opt)
G4bool GetOptionPWACorrection () const
void SetOptionMottCorrection (G4bool opt)
G4bool GetOptionMottCorrection () const
void SetOptionOptimisation (G4bool opt)
G4bool GetOptionOptimisation () const
G4GoudsmitSaundersonTableGetGSTable ()
G4GSPWACorrectionsGetPWACorrection ()
G4GoudsmitSaundersonMscModeloperator= (const G4GoudsmitSaundersonMscModel &right)=delete
 G4GoudsmitSaundersonMscModel (const G4GoudsmitSaundersonMscModel &)=delete
Public Member Functions inherited from G4VMscModel
 G4VMscModel (const G4String &nam)
 ~G4VMscModel () override
void InitialiseParameters (const G4ParticleDefinition *)
void DumpParameters (std::ostream &out) const
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
void SetStepLimitType (G4MscStepLimitType)
void SetLateralDisplasmentFlag (G4bool val)
void SetRangeFactor (G4double)
void SetGeomFactor (G4double)
void SetSkin (G4double)
void SetLambdaLimit (G4double)
void SetSafetyFactor (G4double)
void SetSampleZ (G4bool)
G4VEnergyLossProcessGetIonisation () const
void SetIonisation (G4VEnergyLossProcess *, const G4ParticleDefinition *part)
G4double ComputeSafety (const G4ThreeVector &position, G4double limit=DBL_MAX)
G4double ComputeGeomLimit (const G4Track &, G4double &presafety, G4double limit)
G4double GetDEDX (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetDEDX (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
G4double GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKineticEnergy)
G4double GetEnergy (const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy)
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy, G4double logKinEnergy)
G4VMscModeloperator= (const G4VMscModel &right)=delete
 G4VMscModel (const G4VMscModel &)=delete
Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
virtual ~G4VEmModel ()
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 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 CorrectionsAlongStep (const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, 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 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

Additional Inherited Members

Protected Member Functions inherited from G4VMscModel
G4ParticleChangeForMSCGetParticleChangeForMSC (const G4ParticleDefinition *p=nullptr)
G4double ConvertTrueToGeom (G4double &tLength, G4double &gLength)
void SetUseSplineForMSC (G4bool val)
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 G4VMscModel
G4double facrange = 0.04
G4double facgeom = 2.5
G4double facsafety = 0.6
G4double skin = 1.0
G4double dtrl = 0.05
G4double lambdalimit
G4double geomMin
G4double geomMax
G4ThreeVector fDisplacement
G4MscStepLimitType steppingAlgorithm
G4bool samplez = false
G4bool latDisplasment = true
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 144 of file G4GoudsmitSaundersonMscModel.hh.

Constructor & Destructor Documentation

◆ G4GoudsmitSaundersonMscModel() [1/2]

G4GoudsmitSaundersonMscModel::G4GoudsmitSaundersonMscModel ( const G4String & nam = "GoudsmitSaunderson")

Definition at line 148 of file G4GoudsmitSaundersonMscModel.cc.

149 : G4VMscModel(nam) {
150 currentMaterialIndex = -1;
151 presafety = 0.*mm;
152 //
153 particle = nullptr;
154 currentKinEnergy = 0.0;
155 currentRange = 0.0;
156 currentCouple = nullptr;
157 fParticleChange = nullptr;
158 //
159 // Moliere screeing parameter will be used and (by default) corrections are
160 // appalied to the integrated quantities (screeing parameter, elastic mfp, first
161 // and second moments) derived from the corresponding PWA quantities
162 // this PWA correction is ignored if Mott-correction is set to true because
163 // Mott-correction contains all these corrections as well
164 fIsUsePWACorrection = true;
165 fIsUseMottCorrection = false;
166 // Use the "range-rejection" like optimisation?
167 fIsUseOptimisation = true;
168 //
169 fLambda0 = 0.0; // elastic mean free path
170 fLambda1 = 0.0; // first transport mean free path
171 fScrA = 0.0; // screening parameter
172 fG1 = 0.0; // first transport coef.
173 //
174 fMCtoScrA = 1.0;
175 fMCtoQ1 = 1.0;
176 fMCtoG2PerG1 = 1.0;
177 //
178 fTheTrueStepLenght = 0.;
179 fTheZPathLenght = 0.;
180
181 fTheDisplacementVector.set(0.,0.,0.);
182 fTheNewDirection.set(0.,0.,1.);
183
184 fIsEndedUpOnBoundary = false;
185 fIsMultipleScattering = false;
186 fIsSingleScattering = false;
187 fIsNoScatteringInMSC = false;
188 fIsSimplified = false;
189
190 fGSTable = nullptr;
191 fPWACorrection = nullptr;
192}
G4VMscModel(const G4String &nam)

Referenced by G4GoudsmitSaundersonMscModel(), InitialiseLocal(), and operator=().

◆ ~G4GoudsmitSaundersonMscModel()

G4GoudsmitSaundersonMscModel::~G4GoudsmitSaundersonMscModel ( )
override

Definition at line 195 of file G4GoudsmitSaundersonMscModel.cc.

195 {
196 if (IsMaster()) {
197 if (fGSTable) {
198 delete fGSTable;
199 fGSTable = nullptr;
200 }
201 if (fPWACorrection) {
202 delete fPWACorrection;
203 fPWACorrection = nullptr;
204 }
205 }
206}
G4bool IsMaster() const

◆ G4GoudsmitSaundersonMscModel() [2/2]

G4GoudsmitSaundersonMscModel::G4GoudsmitSaundersonMscModel ( const G4GoudsmitSaundersonMscModel & )
delete

Member Function Documentation

◆ ComputeGeomPathLength()

G4double G4GoudsmitSaundersonMscModel::ComputeGeomPathLength ( G4double truePathLength)
overridevirtual

Implements G4VMscModel.

Definition at line 517 of file G4GoudsmitSaundersonMscModel.cc.

517 {
518 // convert true -> geom (before transportation still at the pre-step point)
519 // (called from the step limitation ComputeTruePathLengthLimit)
520 // "range-rejection" like optimisation (just give whatever geometrical length)
521 if (fIsSimplified) {
522 fTheZPathLenght = std::min(fTheTrueStepLenght, fLambda1);
523 }
524 // Note: as we computed the post-step point, transport vector and its projection
525 // along the pre-step point direction at the step limit phase in
526 // `ComputeTruePathLengthLimit` above, we just return that projection here
527 return fTheZPathLenght;
528}

◆ ComputeTruePathLengthLimit()

G4double G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit ( const G4Track & track,
G4double & currentMinimalStep )
overridevirtual

Implements G4VMscModel.

Definition at line 397 of file G4GoudsmitSaundersonMscModel.cc.

398 {
399 const G4DynamicParticle* dp = track.GetDynamicParticle();
400 G4StepPoint* sp = track.GetStep()->GetPreStepPoint();
401 G4StepStatus stepStatus = sp->GetStepStatus();
402 currentCouple = track.GetMaterialCutsCouple();
403 SetCurrentCouple(currentCouple);
404 currentMaterialIndex = (G4int)currentCouple->GetMaterial()->GetIndex();
405 currentKinEnergy = dp->GetKineticEnergy();
406 currentRange = GetRange(particle,currentKinEnergy,currentCouple,
407 dp->GetLogKineticEnergy());
408 // elastic and first transport mfp, screening parameter and G1 are also set
409 // (Mott-correction will be used if it was requested by the user)
410 fLambda1 = GetTransportMeanFreePath(particle,currentKinEnergy);
411 // Set initial values:
412 // : lengths are initialised to currentMinimalStep which is the true, minimum
413 // step length from all other physics
414 fTheTrueStepLenght = currentMinimalStep;
415 fTheZPathLenght = currentMinimalStep; // will need to be converted
416 fTheDisplacementVector.set(0.0, 0.0, 0.0);
417 fTheNewDirection.set(0.0, 0.0, 1.0);
418
419 // Multiple scattering needs to be sampled ?
420 fIsMultipleScattering = false;
421 // Single scattering needs to be sampled ?
422 fIsSingleScattering = false;
423 // Was zero deflection in multiple scattering sampling ?
424 fIsNoScatteringInMSC = false;
425
426 // === The error-free stepping and boundary crossing
427 presafety = ComputeSafety(sp->GetPosition(), fTheTrueStepLenght);
428 // "range-rejection" like optimisation (i.e. deep inside the volume)
429 fIsSimplified = false;
430 if (fIsUseOptimisation && currentRange < presafety ) {
431 // it's assumed that the e-/e+ never leaves the volume
432 // --> simplified tracking, no MSC step limit, no displacement only deflection
433 fIsSimplified = true;
434 SampleMSC();
435 return ConvertTrueToGeom(fTheTrueStepLenght, currentMinimalStep);
436 } else {
437 // the step limit that corresponds to the accurate stepping and boundary
438 // crossing algorithm
439 // Set skin depth = skin x elastic_mean_free_path
440 const G4double skindepth = skin*fLambda0;
441 // Check if we can try Single Scattering because we are within skindepth
442 // distance from/to a boundary OR the current minimum true-step-length is
443 // shorter than skindepth. NOTICE: the latest has only efficieny reasons
444 // because the MSC angular sampling is fine for any short steps but much
445 // faster to try single scattering in case of short steps.
446 if ((stepStatus == fGeomBoundary) || (presafety < skindepth) || (fTheTrueStepLenght < skindepth)) {
447 //Try single scattering:
448 // - sample distance to next single scattering interaction (sslimit)
449 // - compare to current minimum length
450 // == if sslimit is the shorter:
451 // - set the step length to sslimit
452 // - indicate that single scattering needs to be done
453 // == else : nothing to do
454 //- in both cases, the step length was very short so geometrical and
455 // true path length are the same
456 const G4double sslimit = -1.*fLambda0*G4Log(G4UniformRand());
457 // compare to current minimum step length
458 if (sslimit < fTheTrueStepLenght) {
459 fTheTrueStepLenght = sslimit;
460 fIsSingleScattering = true;
461 }
462 // short step -> true step length equal to geometrical path length
463 fTheZPathLenght = fTheTrueStepLenght;
464 // Set that everything was done in step-limit phase (so no MSC call later)
465 // We will check if we need to perform the single-scattering angular
466 // sampling i.e. if single elastic scattering was the winer!
467 } else {
468 // After checking we know that we cannot try single scattering so we will
469 // need to make an MSC step
470 // Indicate that we need to make and MSC step.
471 fIsMultipleScattering = true;
472 // limit from range factor (keep this to allow stricter control)
473 fTheTrueStepLenght = std::min(fTheTrueStepLenght, facrange*currentRange);
474 // never let the particle go further than the safety if we are out of the skin
475 // if we are here we are out of the skin, presafety > 0.
476 if (fTheTrueStepLenght > presafety) {
477 fTheTrueStepLenght = std::min(fTheTrueStepLenght, presafety);
478 }
479 // make sure that we are still within the aplicability of condensed histry model
480 // i.e. true step length is not longer than 1/2 first transport mean free path.
481 // We schould take into account energy loss along the 1/2*lambda_transport1
482 // step length as well but we don't. So let it just 0.5*lambda_transport1
483 fTheTrueStepLenght = std::min(fTheTrueStepLenght, fLambda1*0.5);
484 }
485 }
486 // === end of step limit
487 // performe single scattering/multiple scattering
488 if (fIsMultipleScattering) {
489 // sample multiple scattering (including zero, one, few and multiple)
490 // the final position will also be calculated using the accurate LLCA
491 // leading to the transport vector, its projection to the current direction
492 // (`fTheZPathLenght`) and the perpendicular plane (`fTheDisplacementVector`)
493 // The MSC scattering able is also used to set the `fTheNewDirection`
494 SampleMSC();
495 } else if (fIsSingleScattering) {
496 // Sample single scattering angular deflection and claculate the new
497 // direction that is set in `fTheNewDirection`
498 // Note: `fTheZPathLenght` set to `true-step-length` while `fTheDisplacementVector`
499 // stays `(0, 0, 0)` vector in this case.
500 const G4double lekin = G4Log(currentKinEnergy);
501 const G4double pt2 = currentKinEnergy*(currentKinEnergy+2.0*CLHEP::electron_mass_c2);
502 const G4double beta2 = pt2/(pt2+CLHEP::electron_mass_c2*CLHEP::electron_mass_c2);
503 G4double cost = fGSTable->SingleScattering(1.0, fScrA, lekin, beta2, currentMaterialIndex);
504 cost = std::max(-1.0, std::min(cost, 1.0));
505 // compute sint
506 const G4double dum = 1.0 - cost;
507 const G4double sint = std::sqrt(dum*(2.0 - dum));
508 const G4double phi = CLHEP::twopi*G4UniformRand();
509 const G4double sinPhi = std::sin(phi);
510 const G4double cosPhi = std::cos(phi);
511 fTheNewDirection.set(sint*cosPhi, sint*sinPhi, cost);
512 } // otherwise: single scattering case because within skin but not elastic
513 // scattering proposed the shortest step (but most likely geometry)
514 return ConvertTrueToGeom(fTheTrueStepLenght, currentMinimalStep);
515}
G4double G4Log(G4double x)
Definition G4Log.hh:169
G4StepStatus
@ fGeomBoundary
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4double GetTransportMeanFreePath(const G4ParticleDefinition *, G4double)
G4StepPoint * GetPreStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Step * GetStep() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4double facrange
G4double skin
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double ConvertTrueToGeom(G4double &tLength, G4double &gLength)
G4double ComputeSafety(const G4ThreeVector &position, G4double limit=DBL_MAX)

◆ ComputeTrueStepLength()

G4double G4GoudsmitSaundersonMscModel::ComputeTrueStepLength ( G4double geomStepLength)
overridevirtual

Implements G4VMscModel.

Definition at line 530 of file G4GoudsmitSaundersonMscModel.cc.

530 {
531 // convert geom -> true (after transportation already at the post-step point)
532 // Note: we computed the post-step point, transport vector and its projection
533 // along the pre-step point direction at the step limit phase in
534 // `ComputeTruePathLengthLimit` above. That computation was bsed on the
535 // actual true step length that was stored. Morover, we ensured that
536 // that post-step point will be either within the safety or used single
537 // scattering within the skin. Therefore, the post-step point is either
538 // exactly as we expected (`geomStepLength == fTheZPathLenght` as
539 // transportation could move the track with `fTheZPathLenght`) or the
540 // boundary was reached (`geomStepLength < fTheZPathLenght`) but then
541 // it was a single scattering step. The true step length is the one
542 // stored in the first case while equal to the geometrical one in the
543 // second case (as it was a single scattering step anyway).
544 fIsEndedUpOnBoundary = false;
545 if (geomStepLength < fTheZPathLenght) {
546 // reached boundary (or very last step) otherwise
547 fIsEndedUpOnBoundary = true;
548 fTheZPathLenght = geomStepLength;
549 fTheTrueStepLenght = geomStepLength;
550 }
551 return fTheTrueStepLenght;
552}

◆ CrossSectionPerVolume()

G4double G4GoudsmitSaundersonMscModel::CrossSectionPerVolume ( const G4Material * mat,
const G4ParticleDefinition * ,
G4double kineticEnergy,
G4double cutEnergy = 0.0,
G4double maxEnergy = DBL_MAX )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 267 of file G4GoudsmitSaundersonMscModel.cc.

271 {
272 // use Moliere's screening (with Mott-corretion if it was requested)
273 kineticEnergy = std::max(kineticEnergy, 10.*CLHEP::eV);
274 // // total mometum square, beta2
275 const G4double pt2 = kineticEnergy*(kineticEnergy + 2.0*electron_mass_c2);
276 const G4double beta2 = pt2/(pt2 + electron_mass_c2*electron_mass_c2);
277 const G4int matindx = static_cast<G4int>(mat->GetIndex());
278 // get the Mott-correcton factors if Mott-correcton was requested by the user
279 fMCtoScrA = 1.0;
280 fMCtoQ1 = 1.0;
281 fMCtoG2PerG1 = 1.0;
282 G4double scpCor = 1.0; // keep this for consistency (this method is not used in normal runs)
283 if (fIsUseMottCorrection) {
284 fGSTable->GetMottCorrectionFactors(G4Log(kineticEnergy), beta2, matindx, fMCtoScrA, fMCtoQ1, fMCtoG2PerG1);
285 // ! no scattering power correction since the current couple is not set before this interface method is called
286 // scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, kineticEnergy);
287 } else if (fIsUsePWACorrection) {
288 fPWACorrection->GetPWACorrectionFactors(G4Log(kineticEnergy), beta2, matindx, fMCtoScrA, fMCtoQ1, fMCtoG2PerG1);
289 // scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, kineticEnergy);
290 }
291 // screening parameter:
292 // - if Mott-corretioncorrection: the Screened-Rutherford times Mott-corretion DCS with this
293 // screening parameter gives back the (elsepa) PWA first transport cross section
294 // - if PWA correction: he Screened-Rutherford DCS with this screening parameter
295 // gives back the (elsepa) PWA first transport cross section
296 const G4double bc = fGSTable->GetMoliereBc(matindx);
297 fScrA = fGSTable->GetMoliereXc2(matindx)/(4.0*pt2*bc)*fMCtoScrA;
298 // elastic mean free path in Geant4 internal lenght units: the neglected (1+screening parameter) term is corrected
299 // (if Mott-corretion: the corrected screening parameter is used for this (1+A) correction + Moliere b_c is also
300 // corrected with the screening parameter correction)
301 fLambda0 = beta2*(1.+fScrA)*fMCtoScrA/(bc*scpCor);
302 // first transport coefficient (if Mott-corretion: the corrected screening parameter is used (it will be fully
303 // consistent with the one used during the pre-computation of the Mott-correted GS angular distributions))
304 fG1 = 2.0*fScrA*((1.0+fScrA)*G4Log(1.0/fScrA+1.0)-1.0);
305 // first transport mean free path
306 fLambda1 = fLambda0/fG1;
307 // return with the macroscopic first transport cross section
308 return fLambda1 > 0.0 ? 1.0/fLambda1 : 1.0E+20;
309}
std::size_t GetIndex() const

◆ GetGSTable()

G4GoudsmitSaundersonTable * G4GoudsmitSaundersonMscModel::GetGSTable ( )
inline

Definition at line 186 of file G4GoudsmitSaundersonMscModel.hh.

186{ return fGSTable; }

Referenced by InitialiseLocal().

◆ GetOptionMottCorrection()

G4bool G4GoudsmitSaundersonMscModel::GetOptionMottCorrection ( ) const
inline

Definition at line 180 of file G4GoudsmitSaundersonMscModel.hh.

180{ return fIsUseMottCorrection; }

Referenced by InitialiseLocal().

◆ GetOptionOptimisation()

G4bool G4GoudsmitSaundersonMscModel::GetOptionOptimisation ( ) const
inline

Definition at line 184 of file G4GoudsmitSaundersonMscModel.hh.

184{ return fIsUseOptimisation; }

◆ GetOptionPWACorrection()

G4bool G4GoudsmitSaundersonMscModel::GetOptionPWACorrection ( ) const
inline

Definition at line 176 of file G4GoudsmitSaundersonMscModel.hh.

176{ return fIsUsePWACorrection; }

Referenced by InitialiseLocal().

◆ GetPWACorrection()

G4GSPWACorrections * G4GoudsmitSaundersonMscModel::GetPWACorrection ( )
inline

Definition at line 188 of file G4GoudsmitSaundersonMscModel.hh.

188{ return fPWACorrection; }

Referenced by InitialiseLocal().

◆ GetTransportMeanFreePath()

G4double G4GoudsmitSaundersonMscModel::GetTransportMeanFreePath ( const G4ParticleDefinition * ,
G4double kineticEnergy )

Definition at line 314 of file G4GoudsmitSaundersonMscModel.cc.

315 {
316 // use Moliere's screening (with Mott-corretion if it was requested)
317 kineticEnergy = std::max(kineticEnergy, 10.0*CLHEP::eV);
318 // total mometum square, beta2
319 const G4double pt2 = kineticEnergy*(kineticEnergy + 2.0*electron_mass_c2);
320 const G4double beta2 = pt2/(pt2 + electron_mass_c2*electron_mass_c2);
321 const G4int matindx = static_cast<G4int>(currentCouple->GetMaterial()->GetIndex());
322 // get the Mott and scattering power correcton factors
323 // (if Mott-correcton was activated)
324 fMCtoScrA = 1.0;
325 fMCtoQ1 = 1.0;
326 fMCtoG2PerG1 = 1.0;
327 G4double scpCor = 1.0;
328 if (fIsUseMottCorrection) {
329 fGSTable->GetMottCorrectionFactors(G4Log(kineticEnergy), beta2, matindx, fMCtoScrA, fMCtoQ1, fMCtoG2PerG1);
330 scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, kineticEnergy);
331 } else if (fIsUsePWACorrection) {
332 fPWACorrection->GetPWACorrectionFactors(G4Log(kineticEnergy), beta2, matindx, fMCtoScrA, fMCtoQ1, fMCtoG2PerG1);
333 // scpCor = fGSTable->ComputeScatteringPowerCorrection(currentCouple, kineticEnergy);
334 }
335 // screening parameter:
336 // - if Mott-corretioncorrection: the Screened-Rutherford times Mott-corretion DCS with this
337 // screening parameter gives back the (elsepa) PWA first transport cross section
338 // - if PWA correction: he Screened-Rutherford DCS with this screening parameter
339 // gives back the (elsepa) PWA first transport cross section
340 const G4double bc = fGSTable->GetMoliereBc(matindx);
341 fScrA = fGSTable->GetMoliereXc2(matindx)/(4.0*pt2*bc)*fMCtoScrA;
342 // elastic mean free path in Geant4 internal lenght units: the neglected (1+screening parameter) term is corrected
343 // (if Mott-corretion: the corrected screening parameter is used for this (1+A) correction + Moliere b_c is also
344 // corrected with the screening parameter correction)
345 fLambda0 = beta2*(1.+fScrA)*fMCtoScrA/(bc*scpCor);
346 // first transport coefficient (if Mott-corretion: the corrected screening parameter is used (it will be fully
347 // consistent with the one used during the pre-computation of the Mott-correted GS angular distributions))
348 fG1 = 2.0*fScrA*((1.0+fScrA)*G4Log(1.0/fScrA+1.0)-1.0);
349 // first transport mean free path
350 fLambda1 = fLambda0/fG1;
351
352 return fLambda1;
353}

Referenced by ComputeTruePathLengthLimit(), and SampleMSC().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 209 of file G4GoudsmitSaundersonMscModel.cc.

209 {
210 SetParticle(p);
212 // -create GoudsmitSaundersonTable and init its Mott-correction member if
213 // Mott-correction was required
214 if (IsMaster()) {
215 // get the Mott-correction flag from EmParameters
216 if (G4EmParameters::Instance()->UseMottCorrection()) {
217 fIsUseMottCorrection = true;
218 }
219 // Mott-correction includes other way of PWA x-section corrections so deactivate it even if it was true
220 // when Mott-correction is activated by the user
221 if (fIsUseMottCorrection) {
222 fIsUsePWACorrection = false;
223 }
224 // clear GS-table
225 if (fGSTable) {
226 delete fGSTable;
227 fGSTable = nullptr;
228 }
229 // clear PWA corrections table if any
230 if (fPWACorrection) {
231 delete fPWACorrection;
232 fPWACorrection = nullptr;
233 }
234 // create GS-table
235 G4bool isElectron = true;
236 if (p->GetPDGCharge()>0.) {
237 isElectron = false;
238 }
239 fGSTable = new G4GoudsmitSaundersonTable(isElectron);
240 // G4GSTable will be initialised:
241 // - Screened-Rutherford DCS based GS angular distributions will be loaded only if they are not there yet
242 // - Mott-correction will be initialised if Mott-correction was requested to be used
243 fGSTable->SetOptionMottCorrection(fIsUseMottCorrection);
244 // - set PWA correction (correction to integrated quantites from Dirac-PWA)
245 fGSTable->SetOptionPWACorrection(fIsUsePWACorrection);
246 // init
247 fGSTable->Initialise(LowEnergyLimit(),HighEnergyLimit());
248 // create PWA corrections table if it was requested (and not deactivated because active Mott-correction)
249 if (fIsUsePWACorrection) {
250 fPWACorrection = new G4GSPWACorrections(isElectron);
251 fPWACorrection->Initialise();
252 }
253 }
254 fParticleChange = GetParticleChangeForMSC(p);
255}
bool G4bool
Definition G4Types.hh:86
static G4EmParameters * Instance()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
void InitialiseParameters(const G4ParticleDefinition *)
G4bool isElectron(G4int ityp)

◆ InitialiseLocal()

void G4GoudsmitSaundersonMscModel::InitialiseLocal ( const G4ParticleDefinition * p,
G4VEmModel * masterModel )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 258 of file G4GoudsmitSaundersonMscModel.cc.

258 {
259 fGSTable = static_cast<G4GoudsmitSaundersonMscModel*>(masterModel)->GetGSTable();
260 fIsUseMottCorrection = static_cast<G4GoudsmitSaundersonMscModel*>(masterModel)->GetOptionMottCorrection();
261 fIsUsePWACorrection = static_cast<G4GoudsmitSaundersonMscModel*>(masterModel)->GetOptionPWACorrection();
262 fPWACorrection = static_cast<G4GoudsmitSaundersonMscModel*>(masterModel)->GetPWACorrection();
263}
G4GoudsmitSaundersonMscModel(const G4String &nam="GoudsmitSaunderson")
G4GoudsmitSaundersonTable * GetGSTable()

◆ operator=()

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

◆ SampleMSC()

void G4GoudsmitSaundersonMscModel::SampleMSC ( )

Definition at line 581 of file G4GoudsmitSaundersonMscModel.cc.

581 {
582 fIsNoScatteringInMSC = false;
583 // Energy loss correction: calculate effective energy and step length
584 // (Eqs.(77) in my notes)
585 //
586 // `currentKinEnergy` is the pre-step point kinetic energy `E_0`
587 // (mean) energy loss (estimate only assuming a step length = `fTheTrueStepLenght`)
588 const G4double eloss = currentKinEnergy - GetEnergy(particle, currentRange - fTheTrueStepLenght, currentCouple);
589 const G4double midEkin = currentKinEnergy - 0.5*eloss; // mid-step energy `\tilde{E} := E_0 - 0.5*\Delta E`
590 const G4double tau = midEkin/electron_mass_c2; // `\tau := \tilde{E}/(mc^2)`
591 const G4double tau2 = tau*tau;
592 const G4double eps0 = eloss/currentKinEnergy; // energy loss fraction to pre-step energy: `\epsilon := \Delta E/E_0`
593 const G4double epsm = eloss/midEkin; // energy loss fraction to the mid-step energy: `\epsilon := \Delta E/\tilde{E}`
594 const G4double epsm2 = epsm*epsm;
595 const G4double effEner = midEkin*(1.0 - epsm2*(6.0 + 10.0*tau + 5.0*tau2)/(24.0*tau2 + 72.0*tau + 48.0));
596 const G4double dum = 1.0/((tau + 1.0)*(tau + 2.0));
597 const G4double effStep = fTheTrueStepLenght*(1.0 - 0.166666*epsm2*dum*dum*(4.0 + tau*(6.0 + tau*(7.0 + tau*(4.0 + tau)))));
598 // compute elastic mfp, first transport mfp, screening parameter, and G1 at this `E_eff`
599 // (with Mott-correction if it was requested by the user)
600 fLambda1 = GetTransportMeanFreePath(particle, effEner);
601 // s/lambda_el i.e. mean number elastic scattering along the step
602 // (using the effective step length and energy)
603 const G4double lambdan = fLambda0 > 0.0 ? effStep/fLambda0 : 0.0;
604 // do nothing in case of very small mean elastic scattering
605 if (lambdan <= 1.0E-12) {
606 fTheZPathLenght = fTheTrueStepLenght;
607 fIsNoScatteringInMSC = true;
608 return;
609 }
610 // first moment: 2.* lambdan *scrA*((1.+scrA)*log(1.+1./scrA)-1.);
611 // (which is equal to Q1 = s*G1/\lambda = s/\lambda_1)
612 G4double Qn1 = lambdan *fG1;
613 // sample scattering angles: divide the step into two and combine at the end
614 // NOTE: this shouldn't happen in normal case as we limit the step such that
615 // `s_max <= 0.5*\lambda_1` leading to `Q1_max = s/|lambda_1 = 0.5`.
616 // But we can see `Q1 > 0.5` in some corner case (e.g. due to the energy
617 // loss along the step as we assumed `\lambda_1` to be constant).
618 // Anyway, our first GS angular distribution table covers the [0.001,0.99]
619 // which is already well above the validity of the MSC theory as well as
620 // the LLCA algorithm that gives accurate post-step positions only up to
621 // Q1 < ~ 0.5. This is why we have that `s_max = 0.5*\lambda_1` step limit.
622 // (while we have a second GS angular distribution table up to Q1=7.99
623 // that is a smooth convergence to isotropic angular distributions)
624 G4double cosTheta1 = 1.0, sinTheta1 = 0.0, cosTheta2 = 1.0, sinTheta2 = 0.0;
625 if (0.5*Qn1 < 7.0) {
626 // (expected to have Qn1 that are usually < ~0.5)
627 // sample 2 scattering cost1, sint1, cost2 and sint2 for half of the step
628 const G4double lekin = G4Log(effEner);
629 const G4double pt2 = effEner*(effEner + 2.0*CLHEP::electron_mass_c2);
630 const G4double beta2 = pt2/(pt2 + CLHEP::electron_mass_c2*CLHEP::electron_mass_c2);
631 // backup GS angular dtr pointer (kinetic energy and delta index in case of Mott-correction)
632 // if the first was an msc sampling (the same will be used if the second is also an msc step)
633 G4GoudsmitSaundersonTable::GSMSCAngularDtr *gsDtr = nullptr;
634 G4int mcEkinIdx = -1;
635 G4int mcDeltIdx = -1;
636 G4double transfPar = 0.0;
637 G4bool isMsc = fGSTable->Sampling(0.5*lambdan, 0.5*Qn1, fScrA, cosTheta1, sinTheta1, lekin, beta2,
638 currentMaterialIndex, &gsDtr, mcEkinIdx, mcDeltIdx, transfPar,
639 true);
640 fGSTable->Sampling(0.5*lambdan, 0.5*Qn1, fScrA, cosTheta2, sinTheta2, lekin, beta2,
641 currentMaterialIndex, &gsDtr, mcEkinIdx, mcDeltIdx, transfPar, !isMsc);
642 if (cosTheta1 + cosTheta2 == 2.0) { // no scattering happened
643 fTheZPathLenght = fTheTrueStepLenght;
644 fIsNoScatteringInMSC = true;
645 return;
646 }
647 } else {
648 // corner case e.g. last step or "range-rejection" like opt.--> isotropic scattering
649 cosTheta1 = 1.-2.*G4UniformRand();
650 sinTheta1 = std::sqrt((1.-cosTheta1)*(1.+cosTheta1));
651 cosTheta2 = 1.-2.*G4UniformRand();
652 sinTheta2 = std::sqrt((1.-cosTheta2)*(1.+cosTheta2));
653 }
654 // sample 2 azimuthal angles
655 const G4double phi1 = CLHEP::twopi*G4UniformRand();
656 const G4double sinPhi1 = std::sin(phi1);
657 const G4double cosPhi1 = std::cos(phi1);
658 const G4double phi2 = CLHEP::twopi*G4UniformRand();
659 const G4double sinPhi2 = std::sin(phi2);
660 const G4double cosPhi2 = std::cos(phi2);
661 // compute final direction after the two scatterings (from the initial 0,0,1 dir)
662 const G4double u2 = sinTheta2*cosPhi2;
663 const G4double v2 = sinTheta2*sinPhi2;
664 const G4double tm = cosTheta1*u2 + sinTheta1*cosTheta2;
665 const G4double uf = tm*cosPhi1 - v2*sinPhi1;
666 const G4double vf = tm*sinPhi1 + v2*cosPhi1;
667 const G4double wf = cosTheta1*cosTheta2 - sinTheta1*u2;
668 // set the new direction (still in the scattering frame that is 0,0,1)
669 fTheNewDirection.set(uf, vf, wf);
670 //
671 // "range-rejection" like optimisation --> no dispalcement
672 if (fIsSimplified) {
673 return;
674 }
675 //
676 // == Compute final position (using the accurate LLCA)
677 // [NIMB 142(1998) Eq.76-80 with energy loss correction from Med.Phys. 27 (2000))
678 // apply correction to Q1
679 Qn1 *= fMCtoQ1;
680 // energy loss correction (gamma, delta, eta with alpha1 and alpha2 for eloss cor.)
681 // note: quantities are based on the SR DCS, gamma and delta with (Mott, screening, DPWA) corrections
682 // while the energy loss correction (alpha1, alpha2) is pure SR DCS based.
683 // compute alpha1, delta + alpha2 and gamma for energy loss correction
684 const G4double tp2 = tau + 2.0;
685 const G4double tp1 = tau + 1.0;
686 const G4double loga = G4Log(1.0 + 1.0/fScrA);
687 const G4double lgf1 = loga*(1.0 + fScrA) - 1.0;
688 const G4double lgf2 = loga*(1.0 + 2.0*fScrA) - 2.0;
689 // energy loss correction factors \alpha_1 and 1.0-alpha1
690 const G4double alpha1 = ((2.0 + tau*tp2)/tp1 - tp1/lgf1)*epsm/tp2;
691 // \alpha2 (plus the -0.25.. contrib from higer order)
692 const G4double alpha2 = 0.4082483*(eps0*tp1/(tp2*lgf1*lgf2) - 0.25*alpha1*alpha1);
693 // then gamma and delta (with energy loss correction: \delta -> \delta + \alpha_2)
694 const G4double gamma = (6.0*fScrA*lgf2*(1.0 + fScrA)/fG1)*fMCtoG2PerG1;
695 // calculate \delta --> \delta + alpha_2
696 const G4double delta = 0.9082483 - (0.1020621 - 0.0263747*gamma)*Qn1 + alpha2;
697 //
698 // ready to calculate the final position:
699 // sample \eta from p(eta)=2*eta i.e. P(eta) = eta_square ;-> P(eta) = rand --> eta = sqrt(rand)
700 const G4double eta = std::sqrt(G4UniformRand());
701 const G4double eta0 = 0.5*(1 - eta)*(1.0 + alpha1); // \eta_0 --> \eta_0(1 + \alpha_1)
702 const G4double eta1 = eta*delta; // \eta_1
703 const G4double eta2 = eta*(1.0 - delta); // \eta_2
704 // calculate the post-step point: xf, yf, zf are the final position coordinates
705 // divided by the (true) step length `s` and as given in NIMB 142 (1998)
706 // with the energy loss correction derived in Med.Phys. 27 (2000)
707 const G4double a1 = 0.5*(1 - eta)*(1.0 - alpha1);
708 const G4double w1v2 = cosTheta1*v2;
709 const G4double dum1 = eta1*sinTheta1;
710 const G4double xf = dum1*cosPhi1 + eta2*(cosPhi1*u2 - sinPhi1*w1v2) + a1*uf;
711 const G4double yf = dum1*sinPhi1 + eta2*(sinPhi1*u2 + cosPhi1*w1v2) + a1*vf;
712 const G4double zf = eta0 + eta1*cosTheta1 + eta2*cosTheta2 + a1*wf;
713 // final position relative to the pre-step point in the scattering frame
714 // (rx, ry, rz) = (xf, yf, zf)*s
715 const G4double rx = xf*fTheTrueStepLenght;
716 const G4double ry = yf*fTheTrueStepLenght;
717 const G4double rz = zf*fTheTrueStepLenght;
718 // calculate the transport distance (with protection: tr_distance <= true_step_length)
719 const G4double transportDistance = std::min(std::sqrt(rx*rx + ry*ry + rz*rz), fTheTrueStepLenght);
720 // set the z-path length i.e. the geometrical length to be the transport
721 // distance the displacement such that the final post-step point will
722 // eventually be at (rx, ry, rz) relative to the pre-step point (with rotations)
723 fTheZPathLenght = transportDistance;
724 fTheDisplacementVector.set(rx, ry, rz - fTheZPathLenght);
725}
const G4double alpha2
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)

Referenced by ComputeTruePathLengthLimit().

◆ SampleScattering()

G4ThreeVector & G4GoudsmitSaundersonMscModel::SampleScattering ( const G4ThreeVector & oldDirection,
G4double safety )
overridevirtual

Implements G4VMscModel.

Definition at line 555 of file G4GoudsmitSaundersonMscModel.cc.

555 {
556 // do nothing on the boundary (reached that in single scattering mode)
557 if (!fIsEndedUpOnBoundary) {
558 // "range-rejection" like optimisation
559 if (fIsSimplified && !fIsNoScatteringInMSC) {
560 fTheNewDirection.rotateUz(oldDirection);
561 fParticleChange->ProposeMomentumDirection(fTheNewDirection);
562 return fTheDisplacementVector;
563 }
564 // if single scattering mode (i.e. within skin) and it's actually happened
565 if (fIsSingleScattering) {
566 fTheNewDirection.rotateUz(oldDirection);
567 fParticleChange->ProposeMomentumDirection(fTheNewDirection);
568 return fTheDisplacementVector;
569 }
570 // if multiple scattering mode (i.e. out of skin) and it's actually happened
571 if (fIsMultipleScattering && !fIsNoScatteringInMSC) {
572 fTheNewDirection.rotateUz(oldDirection);
573 fTheDisplacementVector.rotateUz(oldDirection);
574 fParticleChange->ProposeMomentumDirection(fTheNewDirection);
575 }
576 }
577 // on the boundary: reached in single scattering (or optimisation) -> displacement is (0,0,0))
578 return fTheDisplacementVector;
579}

◆ SetOptionMottCorrection()

void G4GoudsmitSaundersonMscModel::SetOptionMottCorrection ( G4bool opt)
inline

Definition at line 178 of file G4GoudsmitSaundersonMscModel.hh.

178{ fIsUseMottCorrection = opt; }

◆ SetOptionOptimisation()

void G4GoudsmitSaundersonMscModel::SetOptionOptimisation ( G4bool opt)
inline

Definition at line 182 of file G4GoudsmitSaundersonMscModel.hh.

182{ fIsUseOptimisation = opt; }

◆ SetOptionPWACorrection()

void G4GoudsmitSaundersonMscModel::SetOptionPWACorrection ( G4bool opt)
inline

Definition at line 174 of file G4GoudsmitSaundersonMscModel.hh.

174{ fIsUsePWACorrection = opt; }

◆ StartTracking()

void G4GoudsmitSaundersonMscModel::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 392 of file G4GoudsmitSaundersonMscModel.cc.

392 {
393 SetParticle(track->GetDynamicParticle()->GetDefinition());
394}
G4ParticleDefinition * GetDefinition() const

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