Geant4  10.02.p03
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)
 
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, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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=0)
 
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

G4ParticleChangeForMSC * GetParticleChangeForMSC (const G4ParticleDefinition *p=0)
 
G4double ConvertTrueToGeom (G4double &tLength, G4double &gLength)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
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
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4VMscModeloperator= (const G4VMscModel &right)
 
 G4VMscModel (const G4VMscModel &)
 

Private Attributes

G4SafetyHelpersafetyHelper
 
G4VEnergyLossProcessionisation
 
const G4ParticleDefinitioncurrentPart
 
G4LossTableManagerman
 
G4double dedx
 
G4double localtkin
 
G4double localrange
 

Additional Inherited Members

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

Detailed Description

Definition at line 68 of file G4VMscModel.hh.

Constructor & Destructor Documentation

◆ G4VMscModel() [1/2]

G4VMscModel::G4VMscModel ( const G4String nam)

Definition at line 58 of file G4VMscModel.cc.

58  :
59  G4VEmModel(nam),
60  safetyHelper(nullptr),
61  ionisation(nullptr),
62  facrange(0.04),
63  facgeom(2.5),
64  facsafety(0.3),
65  skin(1.0),
66  dtrl(0.05),
67  lambdalimit(mm),
68  geomMin(1.e-6*CLHEP::mm),
69  geomMax(1.e50*CLHEP::mm),
70  fDisplacement(0.,0.,0.),
72  samplez(false),
73  latDisplasment(true)
74 {
77  localtkin = 0.0;
79  currentPart = 0;
80 }
G4double facgeom
Definition: G4VMscModel.hh:178
static G4LossTableManager * Instance()
G4double dtrl
Definition: G4VMscModel.hh:181
G4double facrange
Definition: G4VMscModel.hh:177
G4double lambdalimit
Definition: G4VMscModel.hh:182
G4bool samplez
Definition: G4VMscModel.hh:189
G4double skin
Definition: G4VMscModel.hh:180
G4double dedx
Definition: G4VMscModel.hh:171
G4bool latDisplasment
Definition: G4VMscModel.hh:190
G4double geomMax
Definition: G4VMscModel.hh:184
G4double localtkin
Definition: G4VMscModel.hh:172
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4double geomMin
Definition: G4VMscModel.hh:183
static const double g
static const double cm2
Definition: SystemOfUnits.h:99
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:186
G4double localrange
Definition: G4VMscModel.hh:173
static const double mm
Definition: SystemOfUnits.h:94
const G4ParticleDefinition * currentPart
Definition: G4VMscModel.hh:168
G4double facsafety
Definition: G4VMscModel.hh:179
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:167
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:187
#define DBL_MAX
Definition: templates.hh:83
G4SafetyHelper * safetyHelper
Definition: G4VMscModel.hh:166
static const double mm
Definition: G4SIunits.hh:114
static const double MeV
G4LossTableManager * man
Definition: G4VMscModel.hh:169
Here is the call graph for this function:

◆ ~G4VMscModel()

G4VMscModel::~G4VMscModel ( )
virtual

Definition at line 84 of file G4VMscModel.cc.

85 {}

◆ G4VMscModel() [2/2]

G4VMscModel::G4VMscModel ( const G4VMscModel )
private

Member Function Documentation

◆ ComputeGeomLimit()

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

Definition at line 257 of file G4VMscModel.hh.

260 {
261  /*
262  G4double res = geomMax;
263  if(track.GetVolume() != safetyHelper->GetWorldVolume()) {
264  G4double res = safetyHelper->CheckNextStep(
265  track.GetStep()->GetPreStepPoint()->GetPosition(),
266  track.GetMomentumDirection(),
267  limit, presafety);
268  }
269  return res;
270  */
271  return safetyHelper->CheckNextStep(
272  track.GetStep()->GetPreStepPoint()->GetPosition(),
273  track.GetMomentumDirection(),
274  limit, presafety);
275 }
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
G4SafetyHelper * safetyHelper
Definition: G4VMscModel.hh:166
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeGeomPathLength()

G4double G4VMscModel::ComputeGeomPathLength ( G4double  truePathLength)
virtual

Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 153 of file G4VMscModel.cc.

154 {
155  return truePathLength;
156 }
Here is the caller graph for this function:

◆ ComputeSafety()

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

Definition at line 239 of file G4VMscModel.hh.

241 {
242  return safetyHelper->ComputeSafety(position, limit);
243 }
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
G4SafetyHelper * safetyHelper
Definition: G4VMscModel.hh:166
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeTruePathLengthLimit()

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

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

Definition at line 146 of file G4VMscModel.cc.

147 {
148  return DBL_MAX;
149 }
#define DBL_MAX
Definition: templates.hh:83
Here is the caller graph for this function:

◆ ComputeTrueStepLength()

G4double G4VMscModel::ComputeTrueStepLength ( G4double  geomPathLength)
virtual

Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 160 of file G4VMscModel.cc.

161 {
162  return geomPathLength;
163 }
Here is the caller graph for this function:

◆ ConvertTrueToGeom()

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

Definition at line 247 of file G4VMscModel.hh.

249 {
250  glength = ComputeGeomPathLength(tlength);
251  // should return true length
252  return tlength;
253 }
virtual G4double ComputeGeomPathLength(G4double truePathLength)
Definition: G4VMscModel.cc:153
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetDEDX()

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

Definition at line 280 of file G4VMscModel.hh.

282 {
283  G4double x;
284  if(ionisation) { x = ionisation->GetDEDX(kinEnergy, couple); }
285  else {
286  G4double q = part->GetPDGCharge()*inveplus;
287  x = dedx*q*q;
288  }
289  return x;
290 }
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double dedx
Definition: G4VMscModel.hh:171
static const G4double inveplus
Definition: G4VEmModel.hh:427
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:167
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:

◆ GetEnergy()

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

Definition at line 315 of file G4VMscModel.hh.

317 {
318  G4double e;
319  //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
320  // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin
321  // << G4endl;
322  if(ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
323  else {
324  e = localtkin;
325  if(localrange > range) {
326  G4double q = part->GetPDGCharge()*inveplus;
327  e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
328  }
329  }
330  return e;
331 }
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
const G4Material * GetMaterial() const
G4double dedx
Definition: G4VMscModel.hh:171
G4double GetDensity() const
Definition: G4Material.hh:180
G4double localtkin
Definition: G4VMscModel.hh:172
static const G4double inveplus
Definition: G4VEmModel.hh:427
G4double localrange
Definition: G4VMscModel.hh:173
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:167
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:

◆ GetIonisation()

G4VEnergyLossProcess * G4VMscModel::GetIonisation ( ) const
inline

Definition at line 335 of file G4VMscModel.hh.

336 {
337  return ionisation;
338 }
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:167

◆ GetParticleChangeForMSC()

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

Definition at line 90 of file G4VMscModel.cc.

91 {
92  // recomputed for each new run
93  if(!safetyHelper) {
95  ->GetSafetyHelper();
97  }
98  G4ParticleChangeForMSC* change = nullptr;
99  if (pParticleChange) {
100  change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
101  } else {
102  change = new G4ParticleChangeForMSC();
103  }
104  if(p) {
105 
106  // table is never built for GenericIon
107  if(p->GetParticleName() == "GenericIon") {
108  if(xSectionTable) {
110  delete xSectionTable;
111  xSectionTable = nullptr;
112  }
113 
114  // table is always built for low mass particles
115  } else if(p->GetPDGMass() < 4.5*GeV || ForceBuildTableFlag()) {
116 
117  idxTable = 0;
118  G4LossTableBuilder* builder = man->GetTableBuilder();
119  if(IsMaster()) {
122  emin = std::max(emin, man->MinKinEnergy());
123  emax = std::min(emax, man->MaxKinEnergy());
124  if(emin < emax) {
125  xSectionTable = builder->BuildTableForModel(xSectionTable, this, p,
126  emin, emax, true);
127  }
128  }
129  theDensityFactor = builder->GetDensityFactors();
130  theDensityIdx = builder->GetCoupleIndexes();
131  }
132  }
133  return change;
134 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
size_t idxTable
Definition: G4VEmModel.hh:426
const std::vector< G4double > * GetDensityFactors()
void InitialiseHelper()
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:690
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4LossTableBuilder * GetTableBuilder()
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:425
const G4String & GetParticleName() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4SafetyHelper * GetSafetyHelper() const
static const double GeV
Definition: G4SIunits.hh:214
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:655
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:424
static G4TransportationManager * GetTransportationManager()
static const G4double emax
G4double MaxKinEnergy() const
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:422
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:423
double G4double
Definition: G4Types.hh:76
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:648
G4SafetyHelper * safetyHelper
Definition: G4VMscModel.hh:166
G4double MinKinEnergy() const
void clearAndDestroy()
G4LossTableManager * man
Definition: G4VMscModel.hh:169
const std::vector< G4int > * GetCoupleIndexes()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetRange()

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

Definition at line 295 of file G4VMscModel.hh.

297 {
298  //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
299  // << ionisation << " " << part->GetParticleName()
300  // << G4endl;
301  localtkin = kinEnergy;
302  if(ionisation) {
303  localrange = ionisation->GetRangeForLoss(kinEnergy, couple);
304  } else {
305  G4double q = part->GetPDGCharge()*inveplus;
306  localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
307  }
308  //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
309  return localrange;
310 }
const G4Material * GetMaterial() const
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4double dedx
Definition: G4VMscModel.hh:171
G4double GetDensity() const
Definition: G4Material.hh:180
G4double localtkin
Definition: G4VMscModel.hh:172
static const G4double inveplus
Definition: G4VEmModel.hh:427
G4double localrange
Definition: G4VMscModel.hh:173
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:167
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:

◆ GetTransportMeanFreePath()

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

Definition at line 352 of file G4VMscModel.hh.

354 {
355  G4double x;
356  if(xSectionTable) {
357  G4int idx = CurrentCouple()->GetIndex();
358  x = (*xSectionTable)[(*theDensityIdx)[idx]]->Value(ekin, idxTable)
359  *(*theDensityFactor)[idx]/(ekin*ekin);
360  } else {
361  x = CrossSectionPerVolume(CurrentCouple()->GetMaterial(), part, ekin,
362  0.0, DBL_MAX);
363  }
364  if(0.0 >= x) { x = DBL_MAX; }
365  else { x = 1.0/x; }
366  return x;
367 }
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:258
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
int G4int
Definition: G4Types.hh:78
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:

◆ operator=()

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

◆ SampleScattering()

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

Reimplemented in G4GoudsmitSaundersonMscModel, G4UrbanMscModel, G4WentzelVIModel, and G4WentzelVIRelModel.

Definition at line 139 of file G4VMscModel.cc.

140 {
141  return fDisplacement;
142 }
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:186
Here is the caller graph for this function:

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4DummyModel.

Definition at line 167 of file G4VMscModel.cc.

171 {}

◆ SetGeomFactor()

void G4VMscModel::SetGeomFactor ( G4double  val)
inline

Definition at line 218 of file G4VMscModel.hh.

219 {
220  if(!IsLocked()) { facgeom = val; }
221 }
G4double facgeom
Definition: G4VMscModel.hh:178
G4bool IsLocked() const
Definition: G4VEmModel.hh:833
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetIonisation()

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

Definition at line 342 of file G4VMscModel.hh.

344 {
345  ionisation = p;
346  currentPart = part;
347 }
TString part[npart]
const G4ParticleDefinition * currentPart
Definition: G4VMscModel.hh:168
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:167
Here is the caller graph for this function:

◆ SetLateralDisplasmentFlag()

void G4VMscModel::SetLateralDisplasmentFlag ( G4bool  val)
inline

Definition at line 197 of file G4VMscModel.hh.

198 {
199  if(!IsLocked()) { latDisplasment = val; }
200 }
G4bool latDisplasment
Definition: G4VMscModel.hh:190
G4bool IsLocked() const
Definition: G4VEmModel.hh:833
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetRangeFactor()

void G4VMscModel::SetRangeFactor ( G4double  val)
inline

Definition at line 211 of file G4VMscModel.hh.

212 {
213  if(!IsLocked()) { facrange = val; }
214 }
G4double facrange
Definition: G4VMscModel.hh:177
G4bool IsLocked() const
Definition: G4VEmModel.hh:833
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetSampleZ()

void G4VMscModel::SetSampleZ ( G4bool  val)
inline

Definition at line 232 of file G4VMscModel.hh.

233 {
234  if(!IsLocked()) { samplez = val; }
235 }
G4bool samplez
Definition: G4VMscModel.hh:189
G4bool IsLocked() const
Definition: G4VEmModel.hh:833
Here is the call graph for this function:

◆ SetSkin()

void G4VMscModel::SetSkin ( G4double  val)
inline

Definition at line 204 of file G4VMscModel.hh.

205 {
206  if(!IsLocked()) { skin = val; }
207 }
G4double skin
Definition: G4VMscModel.hh:180
G4bool IsLocked() const
Definition: G4VEmModel.hh:833
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetStepLimitType()

void G4VMscModel::SetStepLimitType ( G4MscStepLimitType  val)
inline

Definition at line 225 of file G4VMscModel.hh.

226 {
227  if(!IsLocked()) { steppingAlgorithm = val; }
228 }
G4bool IsLocked() const
Definition: G4VEmModel.hh:833
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:187
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ currentPart

const G4ParticleDefinition* G4VMscModel::currentPart
private

Definition at line 168 of file G4VMscModel.hh.

◆ dedx

G4double G4VMscModel::dedx
private

Definition at line 171 of file G4VMscModel.hh.

◆ dtrl

G4double G4VMscModel::dtrl
protected

Definition at line 181 of file G4VMscModel.hh.

◆ facgeom

G4double G4VMscModel::facgeom
protected

Definition at line 178 of file G4VMscModel.hh.

◆ facrange

G4double G4VMscModel::facrange
protected

Definition at line 177 of file G4VMscModel.hh.

◆ facsafety

G4double G4VMscModel::facsafety
protected

Definition at line 179 of file G4VMscModel.hh.

◆ fDisplacement

G4ThreeVector G4VMscModel::fDisplacement
protected

Definition at line 186 of file G4VMscModel.hh.

◆ geomMax

G4double G4VMscModel::geomMax
protected

Definition at line 184 of file G4VMscModel.hh.

◆ geomMin

G4double G4VMscModel::geomMin
protected

Definition at line 183 of file G4VMscModel.hh.

◆ ionisation

G4VEnergyLossProcess* G4VMscModel::ionisation
private

Definition at line 167 of file G4VMscModel.hh.

◆ lambdalimit

G4double G4VMscModel::lambdalimit
protected

Definition at line 182 of file G4VMscModel.hh.

◆ latDisplasment

G4bool G4VMscModel::latDisplasment
protected

Definition at line 190 of file G4VMscModel.hh.

◆ localrange

G4double G4VMscModel::localrange
private

Definition at line 173 of file G4VMscModel.hh.

◆ localtkin

G4double G4VMscModel::localtkin
private

Definition at line 172 of file G4VMscModel.hh.

◆ man

G4LossTableManager* G4VMscModel::man
private

Definition at line 169 of file G4VMscModel.hh.

◆ safetyHelper

G4SafetyHelper* G4VMscModel::safetyHelper
private

Definition at line 166 of file G4VMscModel.hh.

◆ samplez

G4bool G4VMscModel::samplez
protected

Definition at line 189 of file G4VMscModel.hh.

◆ skin

G4double G4VMscModel::skin
protected

Definition at line 180 of file G4VMscModel.hh.

◆ steppingAlgorithm

G4MscStepLimitType G4VMscModel::steppingAlgorithm
protected

Definition at line 187 of file G4VMscModel.hh.


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