Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VMscModel Class Reference

#include <G4VMscModel.hh>

Inheritance diagram for G4VMscModel:
Collaboration diagram for G4VMscModel:

Public Member Functions

 G4VMscModel (const G4String &nam)
 
virtual ~G4VMscModel ()
 
virtual G4double ComputeTruePathLengthLimit (const G4Track &track, G4double &stepLimit)
 
virtual G4double ComputeGeomPathLength (G4double truePathLength)
 
virtual G4double ComputeTrueStepLength (G4double geomPathLength)
 
virtual G4ThreeVectorSampleScattering (const G4ThreeVector &, G4double safety)
 
virtual 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 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 GetRange (const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetEnergy (const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
 
G4double GetTransportMeanFreePath (const G4ParticleDefinition *part, G4double kinEnergy)
 
- 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 *, 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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, 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 *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

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 std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 67 of file G4VMscModel.hh.

Constructor & Destructor Documentation

G4VMscModel::G4VMscModel ( const G4String nam)
explicit

Definition at line 60 of file G4VMscModel.cc.

60  :
61  G4VEmModel(nam),
62  safetyHelper(nullptr),
63  ionisation(nullptr),
64  facrange(0.04),
65  facgeom(2.5),
66  facsafety(0.3),
67  skin(1.0),
68  dtrl(0.05),
69  lambdalimit(mm),
70  geomMin(1.e-6*CLHEP::mm),
71  geomMax(1.e50*CLHEP::mm),
72  fDisplacement(0.,0.,0.),
74  samplez(false),
75  latDisplasment(true)
76 {
77  dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g;
78  localrange = DBL_MAX;
79  localtkin = 0.0;
80  currentPart = 0;
81 }
G4double facgeom
Definition: G4VMscModel.hh:176
static constexpr double mm
Definition: G4SIunits.hh:115
G4double dtrl
Definition: G4VMscModel.hh:179
G4double facrange
Definition: G4VMscModel.hh:175
static constexpr double cm2
G4double lambdalimit
Definition: G4VMscModel.hh:180
G4bool samplez
Definition: G4VMscModel.hh:187
G4double skin
Definition: G4VMscModel.hh:178
G4bool latDisplasment
Definition: G4VMscModel.hh:188
G4double geomMax
Definition: G4VMscModel.hh:182
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double g
G4double geomMin
Definition: G4VMscModel.hh:181
static constexpr double mm
Definition: SystemOfUnits.h:95
static constexpr double MeV
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:184
G4double facsafety
Definition: G4VMscModel.hh:177
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:185
#define DBL_MAX
Definition: templates.hh:83
G4VMscModel::~G4VMscModel ( )
virtual

Definition at line 85 of file G4VMscModel.cc.

86 {}

Member Function Documentation

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

Definition at line 255 of file G4VMscModel.hh.

258 {
259  return safetyHelper->CheckNextStep(
260  track.GetStep()->GetPreStepPoint()->GetPosition(),
261  track.GetMomentumDirection(),
262  limit, presafety);
263 }
const G4Step * GetStep() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VMscModel::ComputeGeomPathLength ( G4double  truePathLength)
virtual

Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanAdjointMscModel, G4UrbanMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 156 of file G4VMscModel.cc.

157 {
158  return truePathLength;
159 }

Here is the caller graph for this function:

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

Definition at line 237 of file G4VMscModel.hh.

239 {
240  return safetyHelper->ComputeSafety(position, limit);
241 }
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanAdjointMscModel, G4UrbanMscModel, G4WentzelVIModel, G4WentzelVIRelModel, and G4LowEWentzelVIModel.

Definition at line 149 of file G4VMscModel.cc.

150 {
151  return DBL_MAX;
152 }
#define DBL_MAX
Definition: templates.hh:83

Here is the caller graph for this function:

G4double G4VMscModel::ComputeTrueStepLength ( G4double  geomPathLength)
virtual

Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanAdjointMscModel, G4UrbanMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 163 of file G4VMscModel.cc.

164 {
165  return geomPathLength;
166 }

Here is the caller graph for this function:

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

Definition at line 245 of file G4VMscModel.hh.

247 {
248  glength = ComputeGeomPathLength(tlength);
249  // should return true length
250  return tlength;
251 }
virtual G4double ComputeGeomPathLength(G4double truePathLength)
Definition: G4VMscModel.cc:156

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 268 of file G4VMscModel.hh.

270 {
271  G4double x;
272  if(ionisation) { x = ionisation->GetDEDX(kinEnergy, couple); }
273  else {
274  G4double q = part->GetPDGCharge()*inveplus;
275  x = dedx*q*q;
276  }
277  return x;
278 }
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
static const G4double inveplus
Definition: G4VEmModel.hh:427
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 303 of file G4VMscModel.hh.

305 {
306  G4double e;
307  //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
308  // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin
309  // << G4endl;
310  if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
311  else {
312  e = localtkin;
313  if(localrange > range) {
314  G4double q = part->GetPDGCharge()*inveplus;
315  e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
316  }
317  }
318  return e;
319 }
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4double GetDensity() const
Definition: G4Material.hh:180
static const G4double inveplus
Definition: G4VEmModel.hh:427
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4VEnergyLossProcess * G4VMscModel::GetIonisation ( ) const
inline

Definition at line 323 of file G4VMscModel.hh.

324 {
325  return ionisation;
326 }
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) {
96  ->GetSafetyHelper();
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) {
128  xSectionTable = builder->BuildTableForModel(xSectionTable, this, p,
129  emin, emax, true);
130  }
131  }
132  theDensityFactor = builder->GetDensityFactors();
133  theDensityIdx = builder->GetCoupleIndexes();
134  }
135  }
136  return change;
137 }
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:654
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:647
size_t idxTable
Definition: G4VEmModel.hh:426
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4SafetyHelper * GetSafetyHelper() const
const std::vector< G4double > * GetDensityFactors()
void InitialiseHelper()
G4double MaxKinEnergy() const
static G4LossTableManager * Instance()
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:689
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4String & GetParticleName() const
G4LossTableBuilder * GetTableBuilder()
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:425
G4double MinKinEnergy() const
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:424
static G4TransportationManager * GetTransportationManager()
static const G4double emax
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4EmParameters * Instance()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:422
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:423
double G4double
Definition: G4Types.hh:76
void clearAndDestroy()
const std::vector< G4int > * GetCoupleIndexes()

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 283 of file G4VMscModel.hh.

285 {
286  //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
287  // << ionisation << " " << part->GetParticleName()
288  // << G4endl;
289  localtkin = kinEnergy;
290  if(ionisation) {
291  localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
292  } else {
293  G4double q = part->GetPDGCharge()*inveplus;
294  localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
295  }
296  //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
297  return localrange;
298 }
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double GetDensity() const
Definition: G4Material.hh:180
static const G4double inveplus
Definition: G4VEmModel.hh:427
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 340 of file G4VMscModel.hh.

342 {
343  G4double x;
344  if(xSectionTable) {
345  G4int idx = CurrentCouple()->GetIndex();
346  x = (*xSectionTable)[(*theDensityIdx)[idx]]->Value(ekin, idxTable)
347  *(*theDensityFactor)[idx]/(ekin*ekin);
348  } else {
349  x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
350  0.0, DBL_MAX);
351  }
352  x = (x > 0.0) ? 1.0/x : DBL_MAX;
353  return x;
354 }
size_t idxTable
Definition: G4VEmModel.hh:426
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:257
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:423
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanAdjointMscModel, G4UrbanMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 142 of file G4VMscModel.cc.

143 {
144  return fDisplacement;
145 }
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:184

Here is the caller graph for this function:

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

Implements G4VEmModel.

Reimplemented in G4DummyModel.

Definition at line 170 of file G4VMscModel.cc.

174 {}
void G4VMscModel::SetGeomFactor ( G4double  val)
inline

Definition at line 216 of file G4VMscModel.hh.

217 {
218  if(!IsLocked()) { facgeom = val; }
219 }
G4double facgeom
Definition: G4VMscModel.hh:176
G4bool IsLocked() const
Definition: G4VEmModel.hh:832

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 330 of file G4VMscModel.hh.

332 {
333  ionisation = p;
334  currentPart = part;
335 }
const char * p
Definition: xmltok.h:285

Here is the caller graph for this function:

void G4VMscModel::SetLateralDisplasmentFlag ( G4bool  val)
inline

Definition at line 195 of file G4VMscModel.hh.

196 {
197  if(!IsLocked()) { latDisplasment = val; }
198 }
G4bool latDisplasment
Definition: G4VMscModel.hh:188
G4bool IsLocked() const
Definition: G4VEmModel.hh:832

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VMscModel::SetRangeFactor ( G4double  val)
inline

Definition at line 209 of file G4VMscModel.hh.

210 {
211  if(!IsLocked()) { facrange = val; }
212 }
G4double facrange
Definition: G4VMscModel.hh:175
G4bool IsLocked() const
Definition: G4VEmModel.hh:832

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VMscModel::SetSampleZ ( G4bool  val)
inline

Definition at line 230 of file G4VMscModel.hh.

231 {
232  if(!IsLocked()) { samplez = val; }
233 }
G4bool samplez
Definition: G4VMscModel.hh:187
G4bool IsLocked() const
Definition: G4VEmModel.hh:832

Here is the call graph for this function:

void G4VMscModel::SetSkin ( G4double  val)
inline

Definition at line 202 of file G4VMscModel.hh.

203 {
204  if(!IsLocked()) { skin = val; }
205 }
G4double skin
Definition: G4VMscModel.hh:178
G4bool IsLocked() const
Definition: G4VEmModel.hh:832

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VMscModel::SetStepLimitType ( G4MscStepLimitType  val)
inline

Definition at line 223 of file G4VMscModel.hh.

224 {
225  if(!IsLocked()) { steppingAlgorithm = val; }
226 }
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:185
G4bool IsLocked() const
Definition: G4VEmModel.hh:832

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

G4double G4VMscModel::dtrl
protected

Definition at line 179 of file G4VMscModel.hh.

G4double G4VMscModel::facgeom
protected

Definition at line 176 of file G4VMscModel.hh.

G4double G4VMscModel::facrange
protected

Definition at line 175 of file G4VMscModel.hh.

G4double G4VMscModel::facsafety
protected

Definition at line 177 of file G4VMscModel.hh.

G4ThreeVector G4VMscModel::fDisplacement
protected

Definition at line 184 of file G4VMscModel.hh.

G4double G4VMscModel::geomMax
protected

Definition at line 182 of file G4VMscModel.hh.

G4double G4VMscModel::geomMin
protected

Definition at line 181 of file G4VMscModel.hh.

G4double G4VMscModel::lambdalimit
protected

Definition at line 180 of file G4VMscModel.hh.

G4bool G4VMscModel::latDisplasment
protected

Definition at line 188 of file G4VMscModel.hh.

G4bool G4VMscModel::samplez
protected

Definition at line 187 of file G4VMscModel.hh.

G4double G4VMscModel::skin
protected

Definition at line 178 of file G4VMscModel.hh.

G4MscStepLimitType G4VMscModel::steppingAlgorithm
protected

Definition at line 185 of file G4VMscModel.hh.


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