#include <G4VMscModel.hh>
Inheritance diagram for G4VMscModel:
Definition at line 68 of file G4VMscModel.hh.
G4VMscModel::G4VMscModel | ( | const G4String & | nam | ) |
Definition at line 58 of file G4VMscModel.cc.
References DBL_MAX, and G4LossTableManager::Instance().
00058 : 00059 G4VEmModel(nam), 00060 safetyHelper(0), 00061 ionisation(0), 00062 facrange(0.04), 00063 facgeom(2.5), 00064 facsafety(0.3), 00065 skin(1.0), 00066 dtrl(0.05), 00067 lambdalimit(mm), 00068 geomMin(1.e-6*CLHEP::mm), 00069 geomMax(1.e50*CLHEP::mm), 00070 steppingAlgorithm(fUseSafety), 00071 samplez(false), 00072 latDisplasment(true) 00073 { 00074 dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g; 00075 localrange = DBL_MAX; 00076 localtkin = 0.0; 00077 man = G4LossTableManager::Instance(); 00078 currentPart = 0; 00079 }
G4VMscModel::~G4VMscModel | ( | ) | [virtual] |
G4double G4VMscModel::ComputeGeomLimit | ( | const G4Track & | , | |
G4double & | presafety, | |||
G4double | limit | |||
) | [inline] |
Definition at line 256 of file G4VMscModel.hh.
References G4SafetyHelper::CheckNextStep(), geomMax, G4Track::GetMomentumDirection(), G4StepPoint::GetPosition(), G4Step::GetPreStepPoint(), G4Track::GetStep(), G4Track::GetVolume(), and G4SafetyHelper::GetWorldVolume().
Referenced by G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), and G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit().
00259 { 00260 G4double res = geomMax; 00261 if(track.GetVolume() != safetyHelper->GetWorldVolume()) { 00262 res = safetyHelper->CheckNextStep( 00263 track.GetStep()->GetPreStepPoint()->GetPosition(), 00264 track.GetMomentumDirection(), 00265 limit, presafety); 00266 } 00267 return res; 00268 }
Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanMscModel90, G4UrbanMscModel92, G4UrbanMscModel93, G4UrbanMscModel95, G4UrbanMscModel96, G4WentzelVIModel, and G4WentzelVIRelModel.
Definition at line 145 of file G4VMscModel.cc.
Referenced by ConvertTrueToGeom().
G4double G4VMscModel::ComputeSafety | ( | const G4ThreeVector & | position, | |
G4double | limit | |||
) | [inline] |
Definition at line 238 of file G4VMscModel.hh.
References G4SafetyHelper::ComputeSafety().
Referenced by G4VMultipleScattering::AlongStepDoIt(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), and G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit().
00240 { 00241 return safetyHelper->ComputeSafety(position); 00242 }
G4double G4VMscModel::ComputeTruePathLengthLimit | ( | const G4Track & | track, | |
G4double & | stepLimit | |||
) | [virtual] |
Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanMscModel90, G4UrbanMscModel92, G4UrbanMscModel93, G4UrbanMscModel95, G4UrbanMscModel96, G4WentzelVIModel, and G4WentzelVIRelModel.
Definition at line 138 of file G4VMscModel.cc.
References DBL_MAX.
Referenced by G4VMultipleScattering::AlongStepGetPhysicalInteractionLength().
00139 { 00140 return DBL_MAX; 00141 }
Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanMscModel90, G4UrbanMscModel92, G4UrbanMscModel93, G4UrbanMscModel95, G4UrbanMscModel96, G4WentzelVIModel, and G4WentzelVIRelModel.
Definition at line 152 of file G4VMscModel.cc.
Referenced by G4VMultipleScattering::AlongStepDoIt().
G4double G4VMscModel::ConvertTrueToGeom | ( | G4double & | tLength, | |
G4double & | gLength | |||
) | [inline, protected] |
Definition at line 246 of file G4VMscModel.hh.
References ComputeGeomPathLength().
Referenced by G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), and G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit().
00248 { 00249 glength = ComputeGeomPathLength(tlength); 00250 // should return true length 00251 return tlength; 00252 }
G4double G4VMscModel::GetDEDX | ( | const G4ParticleDefinition * | part, | |
G4double | kineticEnergy, | |||
const G4MaterialCutsCouple * | couple | |||
) | [inline] |
Definition at line 273 of file G4VMscModel.hh.
References G4VEnergyLossProcess::GetDEDX(), and G4ParticleDefinition::GetPDGCharge().
Referenced by G4UrbanMscModel96::SampleScattering(), G4UrbanMscModel95::SampleScattering(), G4UrbanMscModel93::SampleScattering(), G4UrbanMscModel92::SampleScattering(), G4UrbanMscModel90::SampleScattering(), and G4GoudsmitSaundersonMscModel::SampleScattering().
00275 { 00276 G4double x; 00277 if(ionisation) { x = ionisation->GetDEDX(kinEnergy, couple); } 00278 else { 00279 G4double q = part->GetPDGCharge()/CLHEP::eplus; 00280 x = dedx*q*q; 00281 } 00282 return x; 00283 }
G4double G4VMscModel::GetEnergy | ( | const G4ParticleDefinition * | part, | |
G4double | range, | |||
const G4MaterialCutsCouple * | couple | |||
) | [inline] |
Definition at line 304 of file G4VMscModel.hh.
References G4Material::GetDensity(), G4VEnergyLossProcess::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), and G4ParticleDefinition::GetPDGCharge().
Referenced by G4WentzelVIRelModel::ComputeGeomPathLength(), G4WentzelVIModel::ComputeGeomPathLength(), G4UrbanMscModel96::ComputeGeomPathLength(), G4UrbanMscModel95::ComputeGeomPathLength(), G4UrbanMscModel93::ComputeGeomPathLength(), G4UrbanMscModel92::ComputeGeomPathLength(), G4UrbanMscModel90::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4WentzelVIRelModel::ComputeTrueStepLength(), G4WentzelVIModel::ComputeTrueStepLength(), G4UrbanMscModel96::SampleScattering(), G4UrbanMscModel95::SampleScattering(), G4UrbanMscModel93::SampleScattering(), G4UrbanMscModel92::SampleScattering(), G4UrbanMscModel90::SampleScattering(), and G4GoudsmitSaundersonMscModel::SampleScattering().
00306 { 00307 G4double e; 00308 if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); } 00309 else { 00310 e = localtkin; 00311 if(localrange > range) { 00312 G4double q = part->GetPDGCharge()/CLHEP::eplus; 00313 e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity(); 00314 } 00315 } 00316 return e; 00317 }
G4VEnergyLossProcess * G4VMscModel::GetIonisation | ( | ) | const [inline] |
G4ParticleChangeForMSC * G4VMscModel::GetParticleChangeForMSC | ( | const G4ParticleDefinition * | p = 0 |
) | [protected] |
Definition at line 89 of file G4VMscModel.cc.
References G4VEmModel::ForceBuildTableFlag(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4TransportationManager::GetSafetyHelper(), G4LossTableManager::GetTableBuilder(), G4TransportationManager::GetTransportationManager(), G4VEmModel::HighEnergyActivationLimit(), G4VEmModel::HighEnergyLimit(), G4SafetyHelper::InitialiseHelper(), G4VEmModel::LowEnergyActivationLimit(), G4VEmModel::LowEnergyLimit(), G4LossTableManager::MaxKinEnergy(), G4LossTableManager::MinKinEnergy(), G4VEmModel::pParticleChange, G4VEmModel::theDensityFactor, G4VEmModel::theDensityIdx, and G4VEmModel::xSectionTable.
Referenced by G4WentzelVIRelModel::Initialise(), G4WentzelVIModel::Initialise(), G4UrbanMscModel96::Initialise(), G4UrbanMscModel95::Initialise(), G4UrbanMscModel93::Initialise(), G4UrbanMscModel92::Initialise(), G4UrbanMscModel90::Initialise(), and G4GoudsmitSaundersonMscModel::Initialise().
00090 { 00091 if(!safetyHelper) { 00092 safetyHelper = G4TransportationManager::GetTransportationManager() 00093 ->GetSafetyHelper(); 00094 safetyHelper->InitialiseHelper(); 00095 } 00096 G4ParticleChangeForMSC* change = 0; 00097 if (pParticleChange) { 00098 change = static_cast<G4ParticleChangeForMSC*>(pParticleChange); 00099 } else { 00100 change = new G4ParticleChangeForMSC(); 00101 } 00102 if(p) { 00103 00104 // table is never built for GenericIon 00105 if(p->GetParticleName() == "GenericIon") { 00106 if(xSectionTable) { 00107 xSectionTable->clearAndDestroy(); 00108 delete xSectionTable; 00109 xSectionTable = 0; 00110 } 00111 00112 // table is always built for low mass particles 00113 } else if(p->GetPDGMass() < 4.5*GeV || ForceBuildTableFlag()) { 00114 G4double emin = std::max(LowEnergyLimit(), LowEnergyActivationLimit()); 00115 G4double emax = std::min(HighEnergyLimit(), HighEnergyActivationLimit()); 00116 emin = std::max(emin, man->MinKinEnergy()); 00117 emax = std::min(emax, man->MaxKinEnergy()); 00118 G4LossTableBuilder* builder = man->GetTableBuilder(); 00119 xSectionTable = builder->BuildTableForModel(xSectionTable, this, p, 00120 emin, emax, true); 00121 theDensityFactor = builder->GetDensityFactors(); 00122 theDensityIdx = builder->GetCoupleIndexes(); 00123 } 00124 } 00125 return change; 00126 }
G4double G4VMscModel::GetRange | ( | const G4ParticleDefinition * | part, | |
G4double | kineticEnergy, | |||
const G4MaterialCutsCouple * | couple | |||
) | [inline] |
Definition at line 288 of file G4VMscModel.hh.
References G4Material::GetDensity(), G4MaterialCutsCouple::GetMaterial(), G4ParticleDefinition::GetPDGCharge(), and G4VEnergyLossProcess::GetRangeForLoss().
Referenced by G4VMultipleScattering::AlongStepDoIt(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), and G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit().
00290 { 00291 if(ionisation) { 00292 localrange = ionisation->GetRangeForLoss(kinEnergy, couple); 00293 } else { 00294 G4double q = part->GetPDGCharge()/CLHEP::eplus; 00295 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity()); 00296 localtkin = kinEnergy; 00297 } 00298 return localrange; 00299 }
G4double G4VMscModel::GetTransportMeanFreePath | ( | const G4ParticleDefinition * | part, | |
G4double | kinEnergy | |||
) | [inline] |
Definition at line 332 of file G4VMscModel.hh.
References G4VEmModel::CrossSectionPerVolume(), G4VEmModel::CurrentCouple(), DBL_MAX, G4MaterialCutsCouple::GetIndex(), and G4VEmModel::xSectionTable.
Referenced by G4WentzelVIRelModel::ComputeGeomPathLength(), G4WentzelVIModel::ComputeGeomPathLength(), G4UrbanMscModel96::ComputeGeomPathLength(), G4UrbanMscModel95::ComputeGeomPathLength(), G4UrbanMscModel93::ComputeGeomPathLength(), G4UrbanMscModel92::ComputeGeomPathLength(), G4UrbanMscModel90::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4WentzelVIRelModel::ComputeTrueStepLength(), and G4WentzelVIModel::ComputeTrueStepLength().
00334 { 00335 G4double x; 00336 if(xSectionTable) { 00337 G4int idx = CurrentCouple()->GetIndex(); 00338 x = (*xSectionTable)[(*theDensityIdx)[idx]]->Value(ekin) 00339 *(*theDensityFactor)[idx]/(ekin*ekin); 00340 } else { 00341 x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin, 00342 0.0, DBL_MAX); 00343 } 00344 if(0.0 >= x) { x = DBL_MAX; } 00345 else { x = 1.0/x; } 00346 return x; 00347 }
G4ThreeVector & G4VMscModel::SampleScattering | ( | const G4ThreeVector & | , | |
G4double | safety | |||
) | [virtual] |
Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanMscModel90, G4UrbanMscModel92, G4UrbanMscModel93, G4UrbanMscModel95, G4UrbanMscModel96, G4WentzelVIModel, and G4WentzelVIRelModel.
Definition at line 131 of file G4VMscModel.cc.
References fDisplacement.
Referenced by G4VMultipleScattering::AlongStepDoIt().
00132 { 00133 return fDisplacement; 00134 }
void G4VMscModel::SampleSecondaries | ( | std::vector< G4DynamicParticle * > * | , | |
const G4MaterialCutsCouple * | , | |||
const G4DynamicParticle * | , | |||
G4double | tmin, | |||
G4double | tmax | |||
) | [virtual] |
Implements G4VEmModel.
Reimplemented in G4UrbanMscModel90, and G4DummyModel.
Definition at line 159 of file G4VMscModel.cc.
void G4VMscModel::SetGeomFactor | ( | G4double | ) | [inline] |
Definition at line 217 of file G4VMscModel.hh.
References facgeom.
Referenced by G4AdjointhMultipleScattering::InitialiseProcess().
00218 { 00219 facgeom = val; 00220 }
void G4VMscModel::SetIonisation | ( | G4VEnergyLossProcess * | , | |
const G4ParticleDefinition * | part | |||
) | [inline] |
Definition at line 324 of file G4VMscModel.hh.
Referenced by G4VMultipleScattering::StartTracking().
void G4VMscModel::SetLateralDisplasmentFlag | ( | G4bool | val | ) | [inline] |
Definition at line 196 of file G4VMscModel.hh.
References latDisplasment.
Referenced by G4AdjointhMultipleScattering::InitialiseProcess().
00197 { 00198 latDisplasment = val; 00199 }
void G4VMscModel::SetRangeFactor | ( | G4double | ) | [inline] |
Definition at line 210 of file G4VMscModel.hh.
References facrange.
Referenced by G4AdjointhMultipleScattering::InitialiseProcess().
00211 { 00212 facrange = val; 00213 }
void G4VMscModel::SetSampleZ | ( | G4bool | ) | [inline] |
Definition at line 231 of file G4VMscModel.hh.
References samplez.
Referenced by G4UrbanMscModel93::G4UrbanMscModel93(), G4UrbanMscModel95::G4UrbanMscModel95(), and G4UrbanMscModel96::G4UrbanMscModel96().
00232 { 00233 samplez = val; 00234 }
void G4VMscModel::SetSkin | ( | G4double | ) | [inline] |
Definition at line 203 of file G4VMscModel.hh.
References skin.
Referenced by G4AdjointhMultipleScattering::InitialiseProcess().
00204 { 00205 skin = val; 00206 }
void G4VMscModel::SetStepLimitType | ( | G4MscStepLimitType | ) | [inline] |
Definition at line 224 of file G4VMscModel.hh.
References steppingAlgorithm.
Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), and G4VMultipleScattering::PreparePhysicsTable().
00225 { 00226 steppingAlgorithm = val; 00227 }
G4double G4VMscModel::dtrl [protected] |
Definition at line 180 of file G4VMscModel.hh.
Referenced by G4UrbanMscModel96::ComputeGeomPathLength(), G4UrbanMscModel95::ComputeGeomPathLength(), G4UrbanMscModel93::ComputeGeomPathLength(), G4UrbanMscModel92::ComputeGeomPathLength(), G4UrbanMscModel90::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4UrbanMscModel96::SampleScattering(), G4UrbanMscModel95::SampleScattering(), G4UrbanMscModel93::SampleScattering(), G4UrbanMscModel92::SampleScattering(), G4UrbanMscModel90::SampleScattering(), and G4GoudsmitSaundersonMscModel::SampleScattering().
G4double G4VMscModel::facgeom [protected] |
Definition at line 177 of file G4VMscModel.hh.
Referenced by G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), and SetGeomFactor().
G4double G4VMscModel::facrange [protected] |
Definition at line 176 of file G4VMscModel.hh.
Referenced by G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), and SetRangeFactor().
G4double G4VMscModel::facsafety [protected] |
Definition at line 178 of file G4VMscModel.hh.
Referenced by G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), and G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit().
G4ThreeVector G4VMscModel::fDisplacement [protected] |
Definition at line 185 of file G4VMscModel.hh.
Referenced by G4WentzelVIRelModel::SampleScattering(), G4WentzelVIModel::SampleScattering(), SampleScattering(), G4UrbanMscModel96::SampleScattering(), G4UrbanMscModel95::SampleScattering(), G4UrbanMscModel93::SampleScattering(), G4UrbanMscModel92::SampleScattering(), G4UrbanMscModel90::SampleScattering(), and G4GoudsmitSaundersonMscModel::SampleScattering().
G4double G4VMscModel::geomMax [protected] |
G4double G4VMscModel::geomMin [protected] |
Definition at line 182 of file G4VMscModel.hh.
G4double G4VMscModel::lambdalimit [protected] |
Definition at line 181 of file G4VMscModel.hh.
Referenced by G4UrbanMscModel90::ComputeTruePathLengthLimit().
G4bool G4VMscModel::latDisplasment [protected] |
Definition at line 189 of file G4VMscModel.hh.
Referenced by G4WentzelVIRelModel::SampleScattering(), G4WentzelVIModel::SampleScattering(), G4UrbanMscModel96::SampleScattering(), G4UrbanMscModel95::SampleScattering(), G4UrbanMscModel93::SampleScattering(), G4UrbanMscModel92::SampleScattering(), G4UrbanMscModel90::SampleScattering(), and SetLateralDisplasmentFlag().
G4bool G4VMscModel::samplez [protected] |
Definition at line 188 of file G4VMscModel.hh.
Referenced by G4UrbanMscModel96::ComputeGeomPathLength(), G4UrbanMscModel95::ComputeGeomPathLength(), G4UrbanMscModel93::ComputeGeomPathLength(), G4UrbanMscModel92::ComputeGeomPathLength(), G4UrbanMscModel90::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::ComputeGeomPathLength(), G4GoudsmitSaundersonMscModel::G4GoudsmitSaundersonMscModel(), and SetSampleZ().
G4double G4VMscModel::skin [protected] |
Definition at line 179 of file G4VMscModel.hh.
Referenced by G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4UrbanMscModel90::G4UrbanMscModel90(), G4UrbanMscModel92::G4UrbanMscModel92(), G4UrbanMscModel93::G4UrbanMscModel93(), G4UrbanMscModel95::G4UrbanMscModel95(), G4UrbanMscModel96::G4UrbanMscModel96(), G4UrbanMscModel96::Initialise(), G4UrbanMscModel95::Initialise(), G4UrbanMscModel93::Initialise(), G4UrbanMscModel92::Initialise(), G4UrbanMscModel90::Initialise(), G4GoudsmitSaundersonMscModel::Initialise(), and SetSkin().
G4MscStepLimitType G4VMscModel::steppingAlgorithm [protected] |
Definition at line 186 of file G4VMscModel.hh.
Referenced by G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4UrbanMscModel96::ComputeTruePathLengthLimit(), G4UrbanMscModel95::ComputeTruePathLengthLimit(), G4UrbanMscModel93::ComputeTruePathLengthLimit(), G4UrbanMscModel92::ComputeTruePathLengthLimit(), G4UrbanMscModel90::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), and SetStepLimitType().