Geant4 10.7.0
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 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

G4ParticleChangeForMSCGetParticleChangeForMSC (const G4ParticleDefinition *p=nullptr)
 
G4double ConvertTrueToGeom (G4double &tLength, G4double &gLength)
 
- 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
 
G4double facgeom
 
G4double facsafety
 
G4double skin
 
G4double dtrl
 
G4double lambdalimit
 
G4double geomMin
 
G4double geomMax
 
G4ThreeVector fDisplacement
 
G4MscStepLimitType steppingAlgorithm
 
G4bool samplez
 
G4bool latDisplasment
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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 safetyHelper(nullptr),
62 ionisation(nullptr),
63 facrange(0.04),
64 facgeom(2.5),
65 facsafety(0.6),
66 skin(1.0),
67 dtrl(0.05),
68 lambdalimit(1.*CLHEP::mm),
69 geomMin(1.e-6*CLHEP::mm),
70 geomMax(1.e50*CLHEP::mm),
71 fDisplacement(0.,0.,0.),
73 samplez(false),
74 latDisplasment(true)
75{
76 dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g;
77 localrange = DBL_MAX;
78 localtkin = 0.0;
79 currentPart = nullptr;
81}
@ fUseSafety
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:743
G4double dtrl
Definition: G4VMscModel.hh:197
G4double facrange
Definition: G4VMscModel.hh:193
G4bool samplez
Definition: G4VMscModel.hh:205
G4double skin
Definition: G4VMscModel.hh:196
G4double geomMin
Definition: G4VMscModel.hh:199
G4double geomMax
Definition: G4VMscModel.hh:200
G4double lambdalimit
Definition: G4VMscModel.hh:198
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:203
G4bool latDisplasment
Definition: G4VMscModel.hh:206
G4double facsafety
Definition: G4VMscModel.hh:195
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:202
G4double facgeom
Definition: G4VMscModel.hh:194
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4VMscModel()

G4VMscModel::~G4VMscModel ( )
override

Definition at line 85 of file G4VMscModel.cc.

86{}

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

290{
291 return safetyHelper->CheckNextStep(
292 track.GetStep()->GetPreStepPoint()->GetPosition(),
293 track.GetMomentumDirection(),
294 limit, presafety);
295}
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 G4LowEWentzelVIModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::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()

virtual G4double G4VMscModel::ComputeTruePathLengthLimit ( const G4Track track,
G4double stepLimit 
)
pure virtual

◆ ComputeTrueStepLength()

virtual G4double G4VMscModel::ComputeTrueStepLength ( G4double  geomPathLength)
pure virtual

◆ ConvertTrueToGeom()

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

◆ DumpParameters()

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

Definition at line 160 of file G4VMscModel.cc.

161{
162 G4String alg = "UseSafety";
163 if (steppingAlgorithm == fUseDistanceToBoundary) alg = "DistanceToBoundary";
164 else if (steppingAlgorithm == fMinimal) alg = "Minimal";
165 else if (steppingAlgorithm == fUseSafetyPlus) alg = "SafetyPlus";
166
167 out << std::setw(22) << "StepLim=" << alg << " Rfact=" << facrange
168 << " Gfact=" << facgeom << " Sfact=" << facsafety << " DispFlag:" << latDisplasment
169 << " Skin=" << skin << " Llimit=" << lambdalimit << G4endl;
170}
@ fMinimal
@ fUseSafetyPlus
@ fUseDistanceToBoundary
#define G4endl
Definition: G4ios.hh:57

Referenced by G4EmModelManager::DumpModelList().

◆ GetDEDX() [1/2]

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

Definition at line 300 of file G4VMscModel.hh.

302{
303 G4double x;
304 if (ionisation) {
305 x = ionisation->GetDEDX(kinEnergy, couple);
306 } else {
307 const G4double q = part->GetPDGCharge()*inveplus;
308 x = dedx*q*q;
309 }
310 return x;
311}
double G4double
Definition: G4Types.hh:83
G4double GetPDGCharge() const
G4double inveplus
Definition: G4VEmModel.hh:446
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 
)
inline

Definition at line 314 of file G4VMscModel.hh.

316{
317 G4double x;
318 if (ionisation) {
319 x = ionisation->GetDEDX(kinEnergy, couple, logKinEnergy);
320 } else {
321 const G4double q = part->GetPDGCharge()*inveplus;
322 x = dedx*q*q;
323 }
324 return x;
325}

◆ GetEnergy()

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

Definition at line 368 of file G4VMscModel.hh.

370{
371 G4double e;
372 //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
373 // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin
374 // << G4endl;
375 if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
376 else {
377 e = localtkin;
378 if(localrange > range) {
379 G4double q = part->GetPDGCharge()*inveplus;
380 e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
381 }
382 }
383 return e;
384}
const G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:178
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)

Referenced by 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 388 of file G4VMscModel.hh.

389{
390 return ionisation;
391}

◆ GetParticleChangeForMSC()

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

Definition at line 91 of file G4VMscModel.cc.

92{
93 // recomputed for each new run
94 if(!safetyHelper) {
97 safetyHelper->InitialiseHelper();
98 }
99 G4ParticleChangeForMSC* change = nullptr;
100 if (pParticleChange) {
101 change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
102 } else {
103 change = new G4ParticleChangeForMSC();
104 }
105 if(p) {
106
107 // table is never built for GenericIon
108 if(p->GetParticleName() == "GenericIon") {
109 if(xSectionTable) {
111 delete xSectionTable;
112 xSectionTable = nullptr;
113 }
114
115 // table is always built for low mass particles
116 } else if(p->GetPDGMass() < 4.5*GeV || ForceBuildTableFlag()) {
117
119 idxTable = 0;
120 G4LossTableBuilder* builder =
122 if(IsMaster()) {
125 emin = std::max(emin, param->MinKinEnergy());
126 emax = std::min(emax, param->MaxKinEnergy());
127 if(emin < emax) {
129 emin, emax, true);
130 }
131 }
132 }
133 }
134 return change;
135}
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
void clearAndDestroy()
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:440
size_t idxTable
Definition: G4VEmModel.hh:444
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:439
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:659
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:666
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:701

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 
)
inline

Definition at line 330 of file G4VMscModel.hh.

332{
333 //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
334 // << ionisation << " " << part->GetParticleName()
335 // << G4endl;
336 localtkin = kinEnergy;
337 if (ionisation) {
338 localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
339 } else {
340 const G4double q = part->GetPDGCharge()*inveplus;
341 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
342 }
343 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
344 return localrange;
345}
G4double GetRangeForLoss(G4double kineticEnergy, const G4MaterialCutsCouple *)

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

◆ GetRange() [2/2]

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

Definition at line 348 of file G4VMscModel.hh.

350{
351 //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
352 // << ionisation << " " << part->GetParticleName()
353 // << G4endl;
354 localtkin = kinEnergy;
355 if (ionisation) {
356 localrange = ionisation->GetRangeForLoss(kinEnergy, couple, logKinEnergy);
357 } else {
358 const G4double q = part->GetPDGCharge()*inveplus;
359 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
360 }
361 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
362 return localrange;
363}

◆ GetTransportMeanFreePath() [1/2]

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

Definition at line 405 of file G4VMscModel.hh.

407{
408 G4double x;
409 if (xSectionTable) {
410 const G4int idx = CurrentCouple()->GetIndex();
411 x = (*xSectionTable)[idx]->Value(ekin, idxTable)/(ekin*ekin);
412 } else {
413 x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
414 0.0, DBL_MAX);
415 }
416 return (x > 0.0) ? 1.0/x : DBL_MAX;
417}
int G4int
Definition: G4Types.hh:85
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:254
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:480

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

◆ GetTransportMeanFreePath() [2/2]

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

Definition at line 420 of file G4VMscModel.hh.

422{
423 G4double x;
424 if (xSectionTable) {
425 const G4int idx = CurrentCouple()->GetIndex();
426 x = (*xSectionTable)[idx]->LogVectorValue(ekin, logekin)/(ekin*ekin);
427 } else {
428 x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
429 0.0, DBL_MAX);
430 }
431 return (x > 0.0) ? 1.0/x : DBL_MAX;
432}

◆ InitialiseParameters()

void G4VMscModel::InitialiseParameters ( const G4ParticleDefinition part)

Definition at line 139 of file G4VMscModel.cc.

140{
141 if(IsLocked()) { return; }
143 if(std::abs(part->GetPDGEncoding()) == 11) {
145 facrange = param->MscRangeFactor();
147 } else {
149 facrange = param->MscMuHadRangeFactor();
151 }
152 skin = param->MscSkin();
153 facgeom = param->MscGeomFactor();
154 facsafety = param->MscSafetyFactor();
155 lambdalimit = param->MscLambdaLimit();
156}
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
Definition: G4VEmModel.hh:867

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

◆ operator=()

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

◆ SampleScattering()

virtual G4ThreeVector & G4VMscModel::SampleScattering ( const G4ThreeVector ,
G4double  safety 
)
pure virtual

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 174 of file G4VMscModel.cc.

178{}

◆ SetGeomFactor()

void G4VMscModel::SetGeomFactor ( G4double  val)
inline

Definition at line 234 of file G4VMscModel.hh.

235{
236 if(!IsLocked()) { facgeom = val; }
237}

◆ SetIonisation()

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

Definition at line 395 of file G4VMscModel.hh.

397{
398 ionisation = p;
399 currentPart = part;
400}

Referenced by G4VMultipleScattering::PreparePhysicsTable(), G4VMultipleScattering::SetIonisation(), and G4VMultipleScattering::StartTracking().

◆ SetLambdaLimit()

void G4VMscModel::SetLambdaLimit ( G4double  val)
inline

Definition at line 241 of file G4VMscModel.hh.

242{
243 if(!IsLocked()) { lambdalimit = val; }
244}

◆ SetLateralDisplasmentFlag()

void G4VMscModel::SetLateralDisplasmentFlag ( G4bool  val)
inline

Definition at line 213 of file G4VMscModel.hh.

214{
215 if(!IsLocked()) { latDisplasment = val; }
216}

◆ SetRangeFactor()

void G4VMscModel::SetRangeFactor ( G4double  val)
inline

Definition at line 227 of file G4VMscModel.hh.

228{
229 if(!IsLocked()) { facrange = val; }
230}

◆ SetSafetyFactor()

void G4VMscModel::SetSafetyFactor ( G4double  val)
inline

Definition at line 248 of file G4VMscModel.hh.

249{
250 if(!IsLocked()) { facsafety = val; }
251}

◆ SetSampleZ()

void G4VMscModel::SetSampleZ ( G4bool  val)
inline

Definition at line 262 of file G4VMscModel.hh.

263{
264 if(!IsLocked()) { samplez = val; }
265}

◆ SetSkin()

void G4VMscModel::SetSkin ( G4double  val)
inline

Definition at line 220 of file G4VMscModel.hh.

221{
222 if(!IsLocked()) { skin = val; }
223}

◆ SetStepLimitType()

void G4VMscModel::SetStepLimitType ( G4MscStepLimitType  val)
inline

Definition at line 255 of file G4VMscModel.hh.

256{
257 if(!IsLocked()) { steppingAlgorithm = val; }
258}

Member Data Documentation

◆ dtrl

◆ facgeom

◆ facrange

◆ facsafety

◆ fDisplacement

◆ geomMax

G4double G4VMscModel::geomMax
protected

Definition at line 200 of file G4VMscModel.hh.

◆ geomMin

G4double G4VMscModel::geomMin
protected

Definition at line 199 of file G4VMscModel.hh.

◆ lambdalimit

G4double G4VMscModel::lambdalimit
protected

◆ latDisplasment

◆ samplez

G4bool G4VMscModel::samplez
protected

Definition at line 205 of file G4VMscModel.hh.

Referenced by SetSampleZ().

◆ skin

◆ steppingAlgorithm


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