Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VMscModel Class Referenceabstract

#include <G4VMscModel.hh>

+ Inheritance diagram for G4VMscModel:

Public Member Functions

 G4VMscModel (const G4String &nam)
 
 ~G4VMscModel () override
 
virtual G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &stepLimit)=0
 
virtual G4double ComputeGeomPathLength (G4double truePathLength)=0
 
virtual G4double ComputeTrueStepLength (G4double geomPathLength)=0
 
virtual G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety)=0
 
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 Initialise (const G4ParticleDefinition *, const G4DataVector &)=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 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)
 
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)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

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

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 66 of file G4VMscModel.hh.

Constructor & Destructor Documentation

◆ G4VMscModel() [1/2]

G4VMscModel::G4VMscModel ( const G4String & nam)
explicit

Definition at line 59 of file G4VMscModel.cc.

59 :
60 G4VEmModel(nam),
61 lambdalimit(1.*CLHEP::mm),
62 geomMin(1.e-6*CLHEP::mm),
63 geomMax(1.e50*CLHEP::mm),
64 fDisplacement(0.,0.,0.),
66{
67 dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g;
68}
@ fUseSafety
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
G4double geomMin
G4double geomMax
G4double lambdalimit
G4MscStepLimitType steppingAlgorithm
G4ThreeVector fDisplacement

◆ ~G4VMscModel()

G4VMscModel::~G4VMscModel ( )
overridedefault

◆ G4VMscModel() [2/2]

G4VMscModel::G4VMscModel ( const G4VMscModel & )
delete

Member Function Documentation

◆ ComputeGeomLimit()

G4double G4VMscModel::ComputeGeomLimit ( const G4Track & track,
G4double & presafety,
G4double limit )
inline

Definition at line 296 of file G4VMscModel.hh.

299{
300 return safetyHelper->CheckNextStep(
301 track.GetStep()->GetPreStepPoint()->GetPosition(),
302 track.GetMomentumDirection(),
303 limit, presafety);
304}
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
const G4Step * GetStep() const

Referenced by G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), and G4WentzelVIModel::ComputeTruePathLengthLimit().

◆ ComputeGeomPathLength()

virtual G4double G4VMscModel::ComputeGeomPathLength ( G4double truePathLength)
pure virtual

◆ ComputeSafety()

G4double G4VMscModel::ComputeSafety ( const G4ThreeVector & position,
G4double limit = DBL_MAX )
inline

◆ ComputeTruePathLengthLimit()

◆ ComputeTrueStepLength()

◆ ConvertTrueToGeom()

G4double G4VMscModel::ConvertTrueToGeom ( G4double & tLength,
G4double & gLength )
inlineprotected

◆ DumpParameters()

void G4VMscModel::DumpParameters ( std::ostream & out) const

Definition at line 136 of file G4VMscModel.cc.

137{
138 G4String alg = "UseSafety";
139 if (steppingAlgorithm == fUseDistanceToBoundary) alg = "DistanceToBoundary";
140 else if (steppingAlgorithm == fMinimal) alg = "Minimal";
141 else if (steppingAlgorithm == fUseSafetyPlus) alg = "SafetyPlus";
142
143 out << std::setw(18) << "StepLim=" << alg << " Rfact=" << facrange
144 << " Gfact=" << facgeom << " Sfact=" << facsafety << " DispFlag:" << latDisplasment
145 << " Skin=" << skin << " Llim=" << lambdalimit/CLHEP::mm << " mm" << G4endl;
146}
@ fUseSafetyPlus
@ fUseDistanceToBoundary
#define G4endl
Definition G4ios.hh:67
G4double facrange
G4double skin
G4bool latDisplasment
G4double facsafety
G4double facgeom

Referenced by G4EmModelManager::DumpModelList().

◆ GetDEDX() [1/2]

G4double G4VMscModel::GetDEDX ( const G4ParticleDefinition * part,
G4double kineticEnergy,
const G4MaterialCutsCouple * couple )

Definition at line 158 of file G4VMscModel.cc.

160{
161 G4double x;
162 if (nullptr != ionisation) {
163 x = ionisation->GetDEDX(kinEnergy, couple);
164 } else {
165 const G4double q = part->GetPDGCharge()*inveplus;
166 x = dedx*q*q;
167 }
168 return x;
169}
double G4double
Definition G4Types.hh:83
G4double inveplus
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)

Referenced by G4UrbanAdjointMscModel::SampleScattering(), and G4UrbanMscModel::SampleScattering().

◆ GetDEDX() [2/2]

G4double G4VMscModel::GetDEDX ( const G4ParticleDefinition * part,
G4double kineticEnergy,
const G4MaterialCutsCouple * couple,
G4double logKineticEnergy )

Definition at line 173 of file G4VMscModel.cc.

175{
176 G4double x;
177 if (nullptr != ionisation) {
178 x = ionisation->GetDEDX(kinEnergy, couple, logKinEnergy);
179 } else {
180 const G4double q = part->GetPDGCharge()*inveplus;
181 x = dedx*q*q;
182 }
183 return x;
184}

◆ GetEnergy()

G4double G4VMscModel::GetEnergy ( const G4ParticleDefinition * part,
G4double range,
const G4MaterialCutsCouple * couple )

Definition at line 223 of file G4VMscModel.cc.

225{
226 G4double e;
227 //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
228 // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin << G4endl;
229 if(nullptr != ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
230 else {
231 e = localtkin;
232 if(localrange > range) {
233 G4double q = part->GetPDGCharge()*inveplus;
234 e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
235 }
236 }
237 return e;
238}
const G4Material * GetMaterial() const
G4double GetDensity() const
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)

Referenced by G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4UrbanAdjointMscModel::ComputeGeomPathLength(), G4UrbanMscModel::ComputeGeomPathLength(), G4WentzelVIModel::ComputeGeomPathLength(), G4WentzelVIModel::ComputeTrueStepLength(), G4GoudsmitSaundersonMscModel::SampleMSC(), G4UrbanAdjointMscModel::SampleScattering(), and G4UrbanMscModel::SampleScattering().

◆ GetIonisation()

G4VEnergyLossProcess * G4VMscModel::GetIonisation ( ) const
inline

Definition at line 308 of file G4VMscModel.hh.

309{
310 return ionisation;
311}

◆ GetParticleChangeForMSC()

G4ParticleChangeForMSC * G4VMscModel::GetParticleChangeForMSC ( const G4ParticleDefinition * p = nullptr)
protected

Definition at line 77 of file G4VMscModel.cc.

78{
79 // recomputed for each new run
80 if(nullptr == safetyHelper) {
83 safetyHelper->InitialiseHelper();
84 }
85 G4ParticleChangeForMSC* change = nullptr;
86 if (nullptr != pParticleChange) {
87 change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
88 } else {
89 change = new G4ParticleChangeForMSC();
90 }
91 if(IsMaster() && nullptr != p) {
92
93 // table is always built for low mass particles
94 if(p->GetParticleName() != "GenericIon" &&
95 (p->GetPDGMass() < CLHEP::GeV || ForceBuildTableFlag()) ) {
96
98 G4LossTableBuilder* builder =
102 emin = std::max(emin, param->MinKinEnergy());
103 emax = std::min(emax, param->MaxKinEnergy());
104 if(emin < emax) {
106 emin, emax, useSpline);
107 }
108 }
109 }
110 return change;
111}
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
const G4String & GetParticleName() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4PhysicsTable * xSectionTable
G4double LowEnergyLimit() const
G4bool IsMaster() const
G4double HighEnergyLimit() const
G4VParticleChange * pParticleChange
G4double HighEnergyActivationLimit() const
G4double LowEnergyActivationLimit() const
G4bool ForceBuildTableFlag() const

Referenced by G4GoudsmitSaundersonMscModel::Initialise(), G4UrbanAdjointMscModel::Initialise(), G4UrbanMscModel::Initialise(), and G4WentzelVIModel::Initialise().

◆ GetRange() [1/2]

G4double G4VMscModel::GetRange ( const G4ParticleDefinition * part,
G4double kineticEnergy,
const G4MaterialCutsCouple * couple )

Definition at line 188 of file G4VMscModel.cc.

190{
191 // << ionisation << " " << part->GetParticleName() << G4endl;
192 localtkin = kinEnergy;
193 if (nullptr != ionisation) {
194 localrange = ionisation->GetRange(kinEnergy, couple);
195 } else {
196 const G4double q = part->GetPDGCharge()*inveplus;
197 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
198 }
199 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
200 return localrange;
201}
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)

Referenced by G4VMultipleScattering::AlongStepDoIt(), G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), and G4WentzelVIModel::ComputeTruePathLengthLimit().

◆ GetRange() [2/2]

G4double G4VMscModel::GetRange ( const G4ParticleDefinition * part,
G4double kineticEnergy,
const G4MaterialCutsCouple * couple,
G4double logKineticEnergy )

Definition at line 205 of file G4VMscModel.cc.

207{
208 //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
209 // << ionisation << " " << part->GetParticleName() << G4endl;
210 localtkin = kinEnergy;
211 if (nullptr != ionisation) {
212 localrange = ionisation->GetRange(kinEnergy, couple, logKinEnergy);
213 } else {
214 const G4double q = part->GetPDGCharge()*inveplus;
215 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
216 }
217 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
218 return localrange;
219}

◆ GetTransportMeanFreePath() [1/2]

G4double G4VMscModel::GetTransportMeanFreePath ( const G4ParticleDefinition * part,
G4double kinEnergy )
inline

Definition at line 325 of file G4VMscModel.hh.

327{
328 G4double x;
329 if (nullptr != xSectionTable) {
330 x = pFactor*(*xSectionTable)[basedCoupleIndex]->Value(ekin)/(ekin*ekin);
331 } else {
332 x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
333 }
334 return (x > 0.0) ? 1.0/x : DBL_MAX;
335}
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
std::size_t basedCoupleIndex
const G4Material * pBaseMaterial
G4double pFactor
#define DBL_MAX
Definition templates.hh:62

Referenced by G4UrbanAdjointMscModel::ComputeGeomPathLength(), G4UrbanMscModel::ComputeGeomPathLength(), G4WentzelVIModel::ComputeGeomPathLength(), G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), and G4WentzelVIModel::ComputeTrueStepLength().

◆ GetTransportMeanFreePath() [2/2]

G4double G4VMscModel::GetTransportMeanFreePath ( const G4ParticleDefinition * part,
G4double kinEnergy,
G4double logKinEnergy )
inline

Definition at line 340 of file G4VMscModel.hh.

342{
343 G4double x;
344 if (nullptr != xSectionTable) {
345 x = pFactor*(*xSectionTable)[basedCoupleIndex]->LogVectorValue(ekin, logekin)/(ekin*ekin);
346 } else {
347 x = pFactor*CrossSectionPerVolume(pBaseMaterial, part, ekin, 0.0, DBL_MAX);
348 }
349 return (x > 0.0) ? 1.0/x : DBL_MAX;
350}

◆ InitialiseParameters()

void G4VMscModel::InitialiseParameters ( const G4ParticleDefinition * part)

Definition at line 115 of file G4VMscModel.cc.

116{
117 if(IsLocked()) { return; }
119 if(std::abs(part->GetPDGEncoding()) == 11) {
121 facrange = param->MscRangeFactor();
123 } else {
125 facrange = param->MscMuHadRangeFactor();
127 }
128 skin = param->MscSkin();
129 facgeom = param->MscGeomFactor();
130 facsafety = param->MscSafetyFactor();
131 lambdalimit = param->MscLambdaLimit();
132}
G4double MscMuHadRangeFactor() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4double MscSafetyFactor() const
G4MscStepLimitType MscStepLimitType() const
G4double MscGeomFactor() const
G4bool LateralDisplacement() const
G4bool MuHadLateralDisplacement() const
G4double MscLambdaLimit() const
G4double MscRangeFactor() const
G4double MscSkin() const
G4bool IsLocked() const

Referenced by G4GoudsmitSaundersonMscModel::Initialise(), G4UrbanMscModel::Initialise(), and G4WentzelVIModel::Initialise().

◆ operator=()

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

◆ SampleScattering()

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 150 of file G4VMscModel.cc.

154{}

◆ SetGeomFactor()

void G4VMscModel::SetGeomFactor ( G4double val)
inline

Definition at line 243 of file G4VMscModel.hh.

244{
245 if(!IsLocked()) { facgeom = val; }
246}

◆ SetIonisation()

void G4VMscModel::SetIonisation ( G4VEnergyLossProcess * p,
const G4ParticleDefinition * part )
inline

Definition at line 315 of file G4VMscModel.hh.

317{
318 ionisation = p;
319 currentPart = part;
320}

Referenced by G4EmTableUtil::PrepareMscProcess(), and G4VMultipleScattering::StartTracking().

◆ SetLambdaLimit()

void G4VMscModel::SetLambdaLimit ( G4double val)
inline

Definition at line 250 of file G4VMscModel.hh.

251{
252 if(!IsLocked()) { lambdalimit = val; }
253}

◆ SetLateralDisplasmentFlag()

void G4VMscModel::SetLateralDisplasmentFlag ( G4bool val)
inline

Definition at line 222 of file G4VMscModel.hh.

223{
224 if(!IsLocked()) { latDisplasment = val; }
225}

◆ SetRangeFactor()

void G4VMscModel::SetRangeFactor ( G4double val)
inline

Definition at line 236 of file G4VMscModel.hh.

237{
238 if(!IsLocked()) { facrange = val; }
239}

◆ SetSafetyFactor()

void G4VMscModel::SetSafetyFactor ( G4double val)
inline

Definition at line 257 of file G4VMscModel.hh.

258{
259 if(!IsLocked()) { facsafety = val; }
260}

◆ SetSampleZ()

void G4VMscModel::SetSampleZ ( G4bool val)
inline

Definition at line 271 of file G4VMscModel.hh.

272{
273 if(!IsLocked()) { samplez = val; }
274}

◆ SetSkin()

void G4VMscModel::SetSkin ( G4double val)
inline

Definition at line 229 of file G4VMscModel.hh.

230{
231 if(!IsLocked()) { skin = val; }
232}

◆ SetStepLimitType()

void G4VMscModel::SetStepLimitType ( G4MscStepLimitType val)
inline

Definition at line 264 of file G4VMscModel.hh.

265{
266 if(!IsLocked()) { steppingAlgorithm = val; }
267}

◆ SetUseSplineForMSC()

void G4VMscModel::SetUseSplineForMSC ( G4bool val)
inlineprotected

Definition at line 354 of file G4VMscModel.hh.

355{
356 useSpline = val;
357}

Member Data Documentation

◆ dtrl

◆ facgeom

◆ facrange

◆ facsafety

◆ fDisplacement

◆ geomMax

G4double G4VMscModel::geomMax
protected

Definition at line 206 of file G4VMscModel.hh.

◆ geomMin

G4double G4VMscModel::geomMin
protected

Definition at line 205 of file G4VMscModel.hh.

◆ lambdalimit

G4double G4VMscModel::lambdalimit
protected

◆ latDisplasment

◆ samplez

G4bool G4VMscModel::samplez = false
protected

Definition at line 211 of file G4VMscModel.hh.

Referenced by SetSampleZ().

◆ skin

◆ steppingAlgorithm


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