Geant4  10.02.p03
G4VMultipleScattering Class Referenceabstract

#include <G4VMultipleScattering.hh>

Inheritance diagram for G4VMultipleScattering:
Collaboration diagram for G4VMultipleScattering:

Public Member Functions

 G4VMultipleScattering (const G4String &name="msc", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VMultipleScattering ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void PrintInfoDefinition ()
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
 
void StartTracking (G4Track *)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
G4double ContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
 
G4VEmModelSelectModel (G4double kinEnergy, size_t idx)
 
void AddEmModel (G4int order, G4VEmModel *, const G4Region *region=0)
 
void SetEmModel (G4VMscModel *, G4int index=1)
 
G4VMscModelEmModel (G4int index=1) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
void SetIonisation (G4VEnergyLossProcess *)
 
G4bool LateralDisplasmentFlag () const
 
void SetLateralDisplasmentFlag (G4bool val)
 
G4double Skin () const
 
void SetSkin (G4double val)
 
G4double RangeFactor () const
 
void SetRangeFactor (G4double val)
 
G4double GeomFactor () const
 
G4double PolarAngleLimit () const
 
G4MscStepLimitType StepLimitType () const
 
void SetStepLimitType (G4MscStepLimitType val)
 
G4double LowestKinEnergy () const
 
void SetLowestKinEnergy (G4double val)
 
const G4ParticleDefinitionFirstParticle () const
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
G4double GetMeanFreePath (const G4Track &track, G4double, G4ForceCondition *condition)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ParticleChangeForMSC fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Private Member Functions

 G4VMultipleScattering (G4VMultipleScattering &)
 
G4VMultipleScatteringoperator= (const G4VMultipleScattering &right)
 

Private Attributes

G4EmModelManagermodelManager
 
G4LossTableManageremManager
 
G4EmParameterstheParameters
 
G4SafetyHelpersafetyHelper
 
std::vector< G4VMscModel * > mscModels
 
G4int numberOfModels
 
const G4ParticleDefinitionfirstParticle
 
const G4ParticleDefinitioncurrParticle
 
G4MscStepLimitType stepLimit
 
G4double facrange
 
G4double lowestKinEnergy
 
G4bool latDisplacement
 
G4bool isIon
 
G4bool fDispBeyondSafety
 
G4VMscModelcurrentModel
 
G4VEnergyLossProcessfIonisation
 
G4double physStepLimit
 
G4double tPathLength
 
G4double gPathLength
 
G4ThreeVector fNewPosition
 
G4ThreeVector fNewDirection
 
G4bool fPositionChanged
 
G4bool isActive
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 92 of file G4VMultipleScattering.hh.

Constructor & Destructor Documentation

◆ G4VMultipleScattering() [1/2]

G4VMultipleScattering::G4VMultipleScattering ( const G4String name = "msc",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 95 of file G4VMultipleScattering.cc.

96  :
98  numberOfModels(0),
99  firstParticle(nullptr),
100  currParticle(nullptr),
102  facrange(0.04),
103  latDisplacement(true),
104  isIon(false),
105  fDispBeyondSafety(false),
106  fNewPosition(0.,0.,0.),
107  fNewDirection(0.,0.,1.)
108 {
110  SetVerboseLevel(1);
112  if("ionmsc" == name) { firstParticle = G4GenericIon::GenericIon(); }
113 
115 
117  fIonisation = nullptr;
118 
120  safetyHelper = nullptr;
121  fPositionChanged = false;
122  isActive = false;
123 
124  currentModel = nullptr;
127  emManager->Register(this);
128 }
static G4LossTableManager * Instance()
G4EmModelManager * modelManager
G4VEnergyLossProcess * fIonisation
const G4ParticleDefinition * currParticle
G4ParticleChangeForMSC fParticleChange
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void Register(G4VEnergyLossProcess *p)
G4LossTableManager * emManager
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
const G4ParticleDefinition * firstParticle
static G4EmParameters * Instance()
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4MscStepLimitType stepLimit
static const double eV
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
Here is the call graph for this function:

◆ ~G4VMultipleScattering()

G4VMultipleScattering::~G4VMultipleScattering ( )
virtual

Definition at line 132 of file G4VMultipleScattering.cc.

133 {
134  /*
135  if(1 < verboseLevel) {
136  G4cout << "G4VMultipleScattering destruct " << GetProcessName()
137  << G4endl;
138  }
139  */
140  delete modelManager;
141  emManager->DeRegister(this);
142 }
G4EmModelManager * modelManager
void DeRegister(G4VEnergyLossProcess *p)
G4LossTableManager * emManager
Here is the call graph for this function:

◆ G4VMultipleScattering() [2/2]

G4VMultipleScattering::G4VMultipleScattering ( G4VMultipleScattering )
private

Member Function Documentation

◆ AddEmModel()

void G4VMultipleScattering::AddEmModel ( G4int  order,
G4VEmModel p,
const G4Region region = 0 
)

Definition at line 146 of file G4VMultipleScattering.cc.

148 {
149  G4VEmFluctuationModel* fm = nullptr;
150  modelManager->AddEmModel(order, p, fm, region);
151  if(p) { p->SetParticleChange(pParticleChange); }
152 }
G4EmModelManager * modelManager
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AlongStepDoIt()

G4VParticleChange * G4VMultipleScattering::AlongStepDoIt ( const G4Track &  track,
const G4Step &  step 
)
virtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 476 of file G4VMultipleScattering.cc.

477 {
478  fParticleChange.ProposeMomentumDirection(
479  step.GetPostStepPoint()->GetMomentumDirection());
480  fNewPosition = step.GetPostStepPoint()->GetPosition();
481  fParticleChange.ProposePosition(fNewPosition);
482  fPositionChanged = false;
483 
484  G4double geomLength = step.GetStepLength();
485 
486  // very small step - no msc
487  if(!isActive) {
488  tPathLength = geomLength;
489 
490  // sample msc
491  } else {
492  G4double range =
493  currentModel->GetRange(currParticle,track.GetKineticEnergy(),
494  track.GetMaterialCutsCouple());
495 
497 
498  // protection against wrong t->g->t conversion
499  /*
500  if(currParticle->GetPDGMass() > 0.9*GeV)
501  G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
502  << geomLength
503  << " tPathLength= " << tPathLength
504  << " physStepLimit= " << physStepLimit
505  << " dr= " << range - tPathLength
506  << " ekin= " << track.GetKineticEnergy() << G4endl;
507  */
509 
510  // do not sample scattering at the last or at a small step
511  if(tPathLength < range && tPathLength > geomMin) {
512 
514  step.GetPostStepPoint()->GetMomentumDirection(),minSafety);
515 
516  G4double r2 = displacement.mag2();
517  //G4cout << " R= " << sqrt(r2) << " Rmin= " << sqrt(minDisplacement2)
518  // << " flag= " << fDispBeyondSafety << G4endl;
519  if(r2 > minDisplacement2) {
520 
521  fPositionChanged = true;
522  const G4double sFact = 0.99;
523  G4double postSafety =
524  sFact*(step.GetPreStepPoint()->GetSafety() - geomLength);
525  //G4cout<<" R= "<<sqrt(r2)<<" postSafety= "<<postSafety<<G4endl;
526 
527  // far away from geometry boundary
528  if(postSafety > 0.0 && r2 <= postSafety*postSafety) {
529  fNewPosition += displacement;
530 
531  //near the boundary
532  } else {
533  G4double dispR = std::sqrt(r2);
534  postSafety =
535  sFact*safetyHelper->ComputeSafety(fNewPosition, dispR);
536 
537  // displaced point is definitely within the volume
538  //G4cout<<" R= "<<dispR<<" postSafety= "<<postSafety<<G4endl;
539  if(dispR < postSafety) {
540  fNewPosition += displacement;
541 
542  // optional extra mechanism is applied only if a particle
543  // is stopped by the boundary
544  } else if(fDispBeyondSafety && 0.0 == postSafety) {
545  fNewPosition += displacement;
546  G4double maxshift =
547  std::min(2.0*dispR, geomLength*(physStepLimit/tPathLength - 1.0));
548  G4double dist = 0.0;
549  G4double safety = postSafety + dispR;
550  fNewDirection = *(fParticleChange.GetMomentumDirection());
551  /*
552  G4cout << "##MSC before Recheck maxshift= " << maxshift
553  << " postsafety= " << postSafety
554  << " Ekin= " << track.GetKineticEnergy()
555  << " " << track.GetDefinition()->GetParticleName()
556  << G4endl;
557  */
558  // check if it is possible to shift to the boundary
559  // and the shift is not large
561  fNewDirection, maxshift, &dist, &safety)
562  && std::abs(dist) < maxshift) {
563  /*
564  G4cout << "##MSC after Recheck dist= " << dist
565  << " postsafety= " << postSafety
566  << " t= " << tPathLength
567  << " g= " << geomLength
568  << " p= " << physStepLimit
569  << G4endl;
570  */
571  // shift is positive
572  if(dist >= 0.0) {
573  tPathLength *= (1.0 + dist/geomLength);
574  fNewPosition += dist*fNewDirection;
575 
576  // shift is negative cannot be larger than geomLength
577  } else {
578  maxshift = std::min(maxshift, geomLength);
579  if(0.0 < maxshift + dist) {
580  G4ThreeVector postpoint = step.GetPostStepPoint()->GetPosition();
582  G4double R2 = (postpoint - point).mag2();
583  G4double newdist = dist;
584  // check not more than 10 extra boundaries
585  for(G4int i=0; i<10; ++i) {
586  dist = 0.0;
588  point, fNewDirection, maxshift, &dist, &safety)
589  && std::abs(newdist + dist) < maxshift) {
590  point += dist*fNewDirection;
591  G4double R2new = (postpoint - point).mag2();
592  //G4cout << "Backward i= " << i << " dist= " << dist
593  // << " R2= " << R2new << G4endl;
594  if(dist >= 0.0 || R2new > R2) { break; }
595  R2 = R2new;
596  fNewPosition = point;
597  newdist += dist;
598  } else {
599  break;
600  }
601  }
602  tPathLength *= (1.0 + newdist/geomLength);
603  // shift on boundary is not possible for negative disp
604  } else {
605  fNewPosition += displacement*(postSafety/dispR - 1.0);
606  }
607  }
608  // shift on boundary is not possible for any disp
609  } else {
610  fNewPosition += displacement*(postSafety/dispR - 1.0);
611  }
612  // reduced displacement
613  } else if(postSafety > geomMin) {
614  fNewPosition += displacement*(postSafety/dispR);
615 
616  // very small postSafety
617  } else {
618  fPositionChanged = false;
619  }
620  }
621  //safetyHelper->ReLocateWithinVolume(fNewPosition);
622  }
623  }
624  }
625  fParticleChange.ProposeTrueStepLength(tPathLength);
626  //fParticleChange.ProposePosition(fNewPosition);
627  return &fParticleChange;
628 }
static const G4double geomMin
int G4int
Definition: G4Types.hh:78
static const G4double minDisplacement2
double mag2() const
const G4ParticleDefinition * currParticle
G4ParticleChangeForMSC fParticleChange
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:160
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:295
G4bool RecheckDistanceToCurrentBoundary(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=0) const
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
static const G4double minSafety
double G4double
Definition: G4Types.hh:76
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
Definition: G4VMscModel.cc:139
Here is the call graph for this function:

◆ AlongStepGetPhysicalInteractionLength()

G4double G4VMultipleScattering::AlongStepGetPhysicalInteractionLength ( const G4Track &  track,
G4double  previousStepSize,
G4double  currentMinimalStep,
G4double currentSafety,
G4GPILSelection *  selection 
)
virtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 409 of file G4VMultipleScattering.cc.

415 {
416  // get Step limit proposed by the process
417  *selection = NotCandidateForSelection;
418  physStepLimit = gPathLength = tPathLength = currentMinimalStep;
419 
420  G4double ekin = track.GetKineticEnergy();
421  /*
422  G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
423  << " " << currParticle->GetParticleName()
424  << " currMod " << currentModel
425  << G4endl;
426  */
427  // isIon flag is used only to select a model
428  if(isIon) {
429  ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass();
430  }
431 
432  // select new model
433  if(1 < numberOfModels) {
434  currentModel = static_cast<G4VMscModel*>(
435  SelectModel(ekin,track.GetMaterialCutsCouple()->GetIndex()));
436  }
437  // msc is active is model is active, energy above the limit,
438  // and step size is above the limit;
439  // if it is active msc may limit the step
441  && ekin >= lowestKinEnergy) {
442  isActive = true;
443  tPathLength =
445  if (tPathLength < physStepLimit) {
446  *selection = CandidateForSelection;
447  }
448  } else { isActive = false; }
449 
450  //if(currParticle->GetPDGMass() > GeV)
451  /*
452  G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
453  << " " << currParticle->GetParticleName()
454  << " gPathLength= " << gPathLength
455  << " tPathLength= " << tPathLength
456  << " currentMinimalStep= " << currentMinimalStep
457  << " isActive " << isActive << G4endl;
458  */
459  return gPathLength;
460 }
static const G4double geomMin
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:146
float proton_mass_c2
Definition: hepunit.py:275
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:753
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildPhysicsTable()

void G4VMultipleScattering::BuildPhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4VProcess.

Definition at line 286 of file G4VMultipleScattering.cc.

287 {
288  G4String num = part.GetParticleName();
289  if(1 < verboseLevel) {
290  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
291  << GetProcessName()
292  << " and particle " << num
293  << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
294  << G4endl;
295  }
296  G4bool master = true;
297  const G4VMultipleScattering* masterProcess =
298  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
299  if(masterProcess && masterProcess != this) { master = false; }
300 
301  if(firstParticle == &part) {
302  /*
303  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
304  << GetProcessName()
305  << " and particle " << num
306  << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
307  << " " << this
308  << G4endl;
309  */
311 
312  if(!master) {
313  // initialisation of models
314  G4bool printing = true;
316  /*
317  G4cout << "### G4VMultipleScattering::SlaveBuildPhysicsTable() for "
318  << GetProcessName()
319  << " and particle " << num
320  << " Nmod= " << numberOfModels << " " << this
321  << G4endl;
322  */
323  for(G4int i=0; i<numberOfModels; ++i) {
324  G4VMscModel* msc =
325  static_cast<G4VMscModel*>(GetModelByIndex(i, printing));
326  G4VMscModel* msc0=
327  static_cast<G4VMscModel*>(masterProcess->GetModelByIndex(i,printing));
328  msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false);
329  msc->InitialiseLocal(firstParticle, msc0);
330  }
331  }
332  }
333 
334  // explicitly defined printout by particle name
335  if(1 < verboseLevel ||
336  (0 < verboseLevel && (num == "e-" ||
337  num == "e+" || num == "mu+" ||
338  num == "mu-" || num == "proton"||
339  num == "pi+" || num == "pi-" ||
340  num == "kaon+" || num == "kaon-" ||
341  num == "alpha" || num == "anti_proton" ||
342  num == "GenericIon")))
343  {
344  G4cout << G4endl << GetProcessName()
345  << ": for " << num
346  << " SubType= " << GetProcessSubType()
347  << G4endl;
348  PrintInfo();
350  }
351 
352  if(1 < verboseLevel) {
353  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
354  << GetProcessName()
355  << " and particle " << num
356  << G4endl;
357  }
358 }
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
static G4LossTableManager * Instance()
G4EmModelManager * modelManager
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:219
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:826
int G4int
Definition: G4Types.hh:78
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:418
G4LossTableManager * emManager
void DumpModelList(G4int verb)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
const G4ParticleDefinition * firstParticle
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
#define G4endl
Definition: G4ios.hh:61
virtual void PrintInfo()=0
G4int NumberOfModels() const
Here is the call graph for this function:

◆ ContinuousStepLimit()

G4double G4VMultipleScattering::ContinuousStepLimit ( const G4Track &  track,
G4double  previousStepSize,
G4double  currentMinimalStep,
G4double currentSafety 
)

Definition at line 663 of file G4VMultipleScattering.cc.

668 {
669  return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
670  currentSafety);
671 }
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
Here is the call graph for this function:

◆ EmModel()

G4VMscModel * G4VMultipleScattering::EmModel ( G4int  index = 1) const

Definition at line 165 of file G4VMultipleScattering.cc.

166 {
167  G4VMscModel* p = nullptr;
168  if(index >= 0 && index < G4int(mscModels.size())) { p = mscModels[index]; }
169  return p;
170 }
Int_t index
int G4int
Definition: G4Types.hh:78
std::vector< G4VMscModel * > mscModels
Here is the caller graph for this function:

◆ FirstParticle()

const G4ParticleDefinition * G4VMultipleScattering::FirstParticle ( ) const
inline

Definition at line 399 of file G4VMultipleScattering.hh.

400 {
401  return firstParticle;
402 }
const G4ParticleDefinition * firstParticle

◆ GeomFactor()

G4double G4VMultipleScattering::GeomFactor ( ) const
inline

Definition at line 355 of file G4VMultipleScattering.hh.

356 {
357  return theParameters->MscGeomFactor();
358 }
G4double MscGeomFactor() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetContinuousStepLimit()

G4double G4VMultipleScattering::GetContinuousStepLimit ( const G4Track &  track,
G4double  previousStepSize,
G4double  currentMinimalStep,
G4double currentSafety 
)
protectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 647 of file G4VMultipleScattering.cc.

652 {
653  G4GPILSelection selection = NotCandidateForSelection;
655  currentMinimalStep,
656  currentSafety,
657  &selection);
658  return x;
659 }
double G4double
Definition: G4Types.hh:76
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanFreePath()

G4double G4VMultipleScattering::GetMeanFreePath ( const G4Track &  track,
G4double  ,
G4ForceCondition *  condition 
)
protectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 675 of file G4VMultipleScattering.cc.

677 {
678  *condition = Forced;
679  return DBL_MAX;
680 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83

◆ GetModelByIndex()

G4VEmModel * G4VMultipleScattering::GetModelByIndex ( G4int  idx = 0,
G4bool  ver = false 
) const

Definition at line 175 of file G4VMultipleScattering.cc.

176 {
177  return modelManager->GetModel(idx, ver);
178 }
G4EmModelManager * modelManager
G4VEmModel * GetModel(G4int, G4bool ver=false)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseProcess()

virtual void G4VMultipleScattering::InitialiseProcess ( const G4ParticleDefinition )
protectedpure virtual

Implemented in G4AdjointhMultipleScattering, G4hMultipleScattering, G4MuMultipleScattering, and G4eMultipleScattering.

Here is the caller graph for this function:

◆ IsApplicable()

virtual G4bool G4VMultipleScattering::IsApplicable ( const G4ParticleDefinition p)
pure virtual

◆ LateralDisplasmentFlag()

G4bool G4VMultipleScattering::LateralDisplasmentFlag ( ) const
inline

Definition at line 310 of file G4VMultipleScattering.hh.

311 {
312  return latDisplacement;
313 }
Here is the caller graph for this function:

◆ LowestKinEnergy()

G4double G4VMultipleScattering::LowestKinEnergy ( ) const
inline

Definition at line 385 of file G4VMultipleScattering.hh.

386 {
387  return lowestKinEnergy;
388 }

◆ operator=()

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

◆ PolarAngleLimit()

G4double G4VMultipleScattering::PolarAngleLimit ( ) const
inline

Definition at line 362 of file G4VMultipleScattering.hh.

363 {
364  return theParameters->MscThetaLimit();
365 }
G4double MscThetaLimit() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PostStepDoIt()

G4VParticleChange * G4VMultipleScattering::PostStepDoIt ( const G4Track &  track,
const G4Step &   
)
virtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 633 of file G4VMultipleScattering.cc.

634 {
635  fParticleChange.Initialize(track);
636 
637  if(fPositionChanged) {
639  fParticleChange.ProposePosition(fNewPosition);
640  }
641 
642  return &fParticleChange;
643 }
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
G4ParticleChangeForMSC fParticleChange
Here is the call graph for this function:

◆ PostStepGetPhysicalInteractionLength()

G4double G4VMultipleScattering::PostStepGetPhysicalInteractionLength ( const G4Track &  ,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
virtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 465 of file G4VMultipleScattering.cc.

467 {
468  *condition = Forced;
469  //*condition = NotForced;
470  return DBL_MAX;
471 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83

◆ PreparePhysicsTable()

void G4VMultipleScattering::PreparePhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4VProcess.

Definition at line 183 of file G4VMultipleScattering.cc.

184 {
185  G4bool master = true;
186  const G4VMultipleScattering* masterProc =
187  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
188  if(masterProc && masterProc != this) { master = false; }
189 
190  if(!firstParticle) { firstParticle = &part; }
191  if(part.GetParticleType() == "nucleus") {
193  latDisplacement = false;
194  facrange = 0.2;
195  G4String pname = part.GetParticleName();
196  if(pname != "deuteron" && pname != "triton" &&
197  pname != "alpha+" && pname != "helium" &&
198  pname != "alpha" && pname != "He3" &&
199  pname != "hydrogen") {
200 
201  const G4ParticleDefinition* theGenericIon =
203  if(&part == theGenericIon) { isIon = true; }
204 
205  if(theGenericIon && firstParticle != theGenericIon) {
206  G4ProcessManager* pm = theGenericIon->GetProcessManager();
208  size_t n = v->size();
209  for(size_t j=0; j<n; ++j) {
210  if((*v)[j] == this) {
211  firstParticle = theGenericIon;
212  isIon = true;
213  break;
214  }
215  }
216  }
217  }
218  }
219 
220  emManager->PreparePhysicsTable(&part, this, master);
221  currParticle = 0;
222 
223  if(1 < verboseLevel) {
224  G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
225  << GetProcessName()
226  << " and particle " << part.GetParticleName()
227  << " local particle " << firstParticle->GetParticleName()
228  << " isIon= " << isIon
229  << G4endl;
230  }
231 
232  if(firstParticle == &part) {
233 
234  // initialise process
236 
237  // heavy particles and not ions
238  if(!isIon) {
239  if(part.GetPDGMass() > MeV) {
243  } else {
247  }
248  if(latDisplacement) {
250  }
251  }
252  if(master) { SetVerboseLevel(theParameters->Verbose()); }
254 
255  // initialisation of models
257  for(G4int i=0; i<numberOfModels; ++i) {
258  G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
259  msc->SetIonisation(0, firstParticle);
260  msc->SetMasterThread(master);
261  if(0 == i) { currentModel = msc; }
264  msc->SetSkin(theParameters->MscSkin());
265  msc->SetRangeFactor(facrange);
268  G4double emax =
270  msc->SetHighEnergyLimit(emax);
271  }
272 
274  10.0, verboseLevel);
275 
276  if(!safetyHelper) {
278  ->GetSafetyHelper();
280  }
281  }
282 }
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void InitialiseHelper()
G4EmModelManager * modelManager
G4int verboseLevel
Definition: G4VProcess.hh:368
G4MscStepLimitType MscStepLimitType() const
G4VEmModel * GetModel(G4int, G4bool ver=false)
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:197
const G4String & GetParticleType() const
G4ProcessManager * GetProcessManager() const
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:225
int G4int
Definition: G4Types.hh:78
G4MscStepLimitType MscMuHadStepLimitType() const
G4int Verbose() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
G4bool MuHadLateralDisplacement() const
const G4ParticleDefinition * currParticle
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
Char_t n[5]
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double MscRangeFactor() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
bool G4bool
Definition: G4Types.hh:79
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:211
TString part[npart]
G4SafetyHelper * GetSafetyHelper() const
G4int size() const
G4LossTableManager * emManager
string pname
Definition: eplot.py:33
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
G4double MscThetaLimit() const
static G4TransportationManager * GetTransportationManager()
static const G4double emax
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:711
G4bool LatDisplacementBeyondSafety() const
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:342
static G4ParticleTable * GetParticleTable()
G4double MaxKinEnergy() const
G4double MscSkin() const
const G4ParticleDefinition * firstParticle
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4MscStepLimitType stepLimit
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:218
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4int WorkerVerbose() const
void SetSkin(G4double)
Definition: G4VMscModel.hh:204
G4double MscMuHadRangeFactor() const
G4int NumberOfModels() const
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:760
G4bool LateralDisplacement() const
G4double MscGeomFactor() const
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
Here is the call graph for this function:

◆ PrintInfo()

virtual void G4VMultipleScattering::PrintInfo ( )
pure virtual

Implemented in G4AdjointhMultipleScattering, G4hMultipleScattering, G4MuMultipleScattering, and G4eMultipleScattering.

Here is the caller graph for this function:

◆ PrintInfoDefinition()

void G4VMultipleScattering::PrintInfoDefinition ( )

Definition at line 362 of file G4VMultipleScattering.cc.

363 {
364  if (0 < verboseLevel) {
365  G4cout << G4endl << GetProcessName()
366  << ": for " << firstParticle->GetParticleName()
367  << " SubType= " << GetProcessSubType()
368  << G4endl;
369  PrintInfo();
371  }
372 }
G4EmModelManager * modelManager
G4int verboseLevel
Definition: G4VProcess.hh:368
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
void DumpModelList(G4int verb)
const G4ParticleDefinition * firstParticle
#define G4endl
Definition: G4ios.hh:61
virtual void PrintInfo()=0
Here is the call graph for this function:

◆ ProcessDescription()

void G4VMultipleScattering::ProcessDescription ( std::ostream &  outFile) const
virtual

Definition at line 748 of file G4VMultipleScattering.cc.

749 {
750  outFile << "Multiple scattering process <" << GetProcessName()
751  << ">" << G4endl;
752 }
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ RangeFactor()

G4double G4VMultipleScattering::RangeFactor ( ) const
inline

Definition at line 338 of file G4VMultipleScattering.hh.

339 {
340  return facrange;
341 }
Here is the caller graph for this function:

◆ RetrievePhysicsTable()

G4bool G4VMultipleScattering::RetrievePhysicsTable ( const G4ParticleDefinition ,
const G4String directory,
G4bool  ascii 
)
virtual

Reimplemented from G4VProcess.

Definition at line 729 of file G4VMultipleScattering.cc.

732 {
733  return true;
734 }

◆ SelectModel()

G4VEmModel * G4VMultipleScattering::SelectModel ( G4double  kinEnergy,
size_t  idx 
)
inline

Definition at line 303 of file G4VMultipleScattering.hh.

304 {
305  return modelManager->SelectModel(kinEnergy, coupleIndex);
306 }
G4EmModelManager * modelManager
G4VEmModel * SelectModel(G4double &energy, size_t &index)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetEmModel()

void G4VMultipleScattering::SetEmModel ( G4VMscModel p,
G4int  index = 1 
)

Definition at line 156 of file G4VMultipleScattering.cc.

157 {
158  G4int n = mscModels.size();
159  if(index >= n) { for(G4int i=n; i<=index; ++i) { mscModels.push_back(0); } }
160  mscModels[index] = p;
161 }
Int_t index
int G4int
Definition: G4Types.hh:78
Char_t n[5]
std::vector< G4VMscModel * > mscModels
Here is the caller graph for this function:

◆ SetIonisation()

void G4VMultipleScattering::SetIonisation ( G4VEnergyLossProcess p)

Definition at line 738 of file G4VMultipleScattering.cc.

739 {
740  for(G4int i=0; i<numberOfModels; ++i) {
741  G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i, true));
742  msc->SetIonisation(p, firstParticle);
743  }
744 }
int G4int
Definition: G4Types.hh:78
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:342
const G4ParticleDefinition * firstParticle
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
Here is the call graph for this function:

◆ SetLateralDisplasmentFlag()

void G4VMultipleScattering::SetLateralDisplasmentFlag ( G4bool  val)
inline

Definition at line 317 of file G4VMultipleScattering.hh.

318 {
319  latDisplacement = val;
320 }

◆ SetLowestKinEnergy()

void G4VMultipleScattering::SetLowestKinEnergy ( G4double  val)
inline

Definition at line 392 of file G4VMultipleScattering.hh.

393 {
394  lowestKinEnergy = val;
395 }

◆ SetRangeFactor()

void G4VMultipleScattering::SetRangeFactor ( G4double  val)
inline

Definition at line 345 of file G4VMultipleScattering.hh.

346 {
347  if(val > 0.0 && val < 1.0) {
348  facrange = val;
350  }
351 }
void SetMscRangeFactor(G4double val)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetSkin()

void G4VMultipleScattering::SetSkin ( G4double  val)
inline

Definition at line 331 of file G4VMultipleScattering.hh.

332 {
334 }
void SetMscSkin(G4double val)
Here is the call graph for this function:

◆ SetStepLimitType()

void G4VMultipleScattering::SetStepLimitType ( G4MscStepLimitType  val)
inline

Definition at line 376 of file G4VMultipleScattering.hh.

377 {
378  stepLimit = val;
379  if(val == fMinimal) { SetRangeFactor(0.2); }
381 }
void SetMscStepLimitType(G4MscStepLimitType val)
G4MscStepLimitType stepLimit
void SetRangeFactor(G4double val)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Skin()

G4double G4VMultipleScattering::Skin ( ) const
inline

Definition at line 324 of file G4VMultipleScattering.hh.

325 {
326  return theParameters->MscSkin();
327 }
G4double MscSkin() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ StartTracking()

void G4VMultipleScattering::StartTracking ( G4Track *  track)
virtual

Reimplemented from G4VProcess.

Definition at line 376 of file G4VMultipleScattering.cc.

377 {
378  G4VEnergyLossProcess* eloss = nullptr;
379  if(track->GetParticleDefinition() != currParticle) {
380  currParticle = track->GetParticleDefinition();
382  eloss = fIonisation;
383  }
384  /*
385  G4cout << "G4VMultipleScattering::StartTracking Nmod= " << numberOfModels
386  << " " << currParticle->GetParticleName()
387  << " E(MeV)= " << track->GetKineticEnergy()
388  << " Ion= " << eloss << " " << fIonisation << " IsMaster= "
389  << G4LossTableManager::Instance()->IsMaster()
390  << G4endl;
391  */
392  // one model
393  if(1 == numberOfModels) {
396 
397  // many models
398  } else {
399  for(G4int i=0; i<numberOfModels; ++i) {
400  G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i,true));
401  msc->StartTracking(track);
402  if(eloss) { msc->SetIonisation(fIonisation, currParticle); }
403  }
404  }
405 }
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:284
int G4int
Definition: G4Types.hh:78
G4VEnergyLossProcess * fIonisation
const G4ParticleDefinition * currParticle
G4LossTableManager * emManager
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:342
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
Here is the call graph for this function:

◆ StepLimitType()

G4MscStepLimitType G4VMultipleScattering::StepLimitType ( ) const
inline

Definition at line 369 of file G4VMultipleScattering.hh.

370 {
371  return stepLimit;
372 }
G4MscStepLimitType stepLimit
Here is the caller graph for this function:

◆ StorePhysicsTable()

G4bool G4VMultipleScattering::StorePhysicsTable ( const G4ParticleDefinition part,
const G4String directory,
G4bool  ascii = false 
)
virtual

Reimplemented from G4VProcess.

Definition at line 685 of file G4VMultipleScattering.cc.

688 {
689  G4bool yes = true;
690  if(part != firstParticle) { return yes; }
691  const G4VMultipleScattering* masterProcess =
692  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
693  if(masterProcess && masterProcess != this) { return yes; }
694 
696  static const G4String ss[4] = {"1","2","3","4"};
697  for(G4int i=0; i<nmod; ++i) {
698  G4VEmModel* msc = modelManager->GetModel(i);
699  yes = true;
700  G4PhysicsTable* table = msc->GetCrossSectionTable();
701  if (table) {
702  G4int j = std::min(i,3);
703  G4String name =
704  GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii);
705  yes = table->StorePhysicsTable(name,ascii);
706 
707  if ( yes ) {
708  if ( verboseLevel>0 ) {
709  G4cout << "Physics table are stored for "
710  << part->GetParticleName()
711  << " and process " << GetProcessName()
712  << " with a name <" << name << "> " << G4endl;
713  }
714  } else {
715  G4cout << "Fail to store Physics Table for "
716  << part->GetParticleName()
717  << " and process " << GetProcessName()
718  << " in the directory <" << directory
719  << "> " << G4endl;
720  }
721  }
722  }
723  return yes;
724 }
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
G4EmModelManager * modelManager
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VEmModel * GetModel(G4int, G4bool ver=false)
G4String name
Definition: TRTMaterials.hh:40
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:826
int G4int
Definition: G4Types.hh:78
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:186
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4ParticleDefinition * firstParticle
#define G4endl
Definition: G4ios.hh:61
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
G4int NumberOfModels() const
Here is the call graph for this function:

Member Data Documentation

◆ currentModel

G4VMscModel* G4VMultipleScattering::currentModel
private

Definition at line 285 of file G4VMultipleScattering.hh.

◆ currParticle

const G4ParticleDefinition* G4VMultipleScattering::currParticle
private

Definition at line 265 of file G4VMultipleScattering.hh.

◆ emManager

G4LossTableManager* G4VMultipleScattering::emManager
private

Definition at line 254 of file G4VMultipleScattering.hh.

◆ facrange

G4double G4VMultipleScattering::facrange
private

Definition at line 269 of file G4VMultipleScattering.hh.

◆ fDispBeyondSafety

G4bool G4VMultipleScattering::fDispBeyondSafety
private

Definition at line 274 of file G4VMultipleScattering.hh.

◆ fIonisation

G4VEnergyLossProcess* G4VMultipleScattering::fIonisation
private

Definition at line 286 of file G4VMultipleScattering.hh.

◆ firstParticle

const G4ParticleDefinition* G4VMultipleScattering::firstParticle
private

Definition at line 264 of file G4VMultipleScattering.hh.

◆ fNewDirection

G4ThreeVector G4VMultipleScattering::fNewDirection
private

Definition at line 293 of file G4VMultipleScattering.hh.

◆ fNewPosition

G4ThreeVector G4VMultipleScattering::fNewPosition
private

Definition at line 292 of file G4VMultipleScattering.hh.

◆ fParticleChange

G4ParticleChangeForMSC G4VMultipleScattering::fParticleChange
protected

Definition at line 280 of file G4VMultipleScattering.hh.

◆ fPositionChanged

G4bool G4VMultipleScattering::fPositionChanged
private

Definition at line 294 of file G4VMultipleScattering.hh.

◆ gPathLength

G4double G4VMultipleScattering::gPathLength
private

Definition at line 290 of file G4VMultipleScattering.hh.

◆ isActive

G4bool G4VMultipleScattering::isActive
private

Definition at line 295 of file G4VMultipleScattering.hh.

◆ isIon

G4bool G4VMultipleScattering::isIon
private

Definition at line 273 of file G4VMultipleScattering.hh.

◆ latDisplacement

G4bool G4VMultipleScattering::latDisplacement
private

Definition at line 272 of file G4VMultipleScattering.hh.

◆ lowestKinEnergy

G4double G4VMultipleScattering::lowestKinEnergy
private

Definition at line 270 of file G4VMultipleScattering.hh.

◆ modelManager

G4EmModelManager* G4VMultipleScattering::modelManager
private

Definition at line 253 of file G4VMultipleScattering.hh.

◆ mscModels

std::vector<G4VMscModel*> G4VMultipleScattering::mscModels
private

Definition at line 261 of file G4VMultipleScattering.hh.

◆ numberOfModels

G4int G4VMultipleScattering::numberOfModels
private

Definition at line 262 of file G4VMultipleScattering.hh.

◆ physStepLimit

G4double G4VMultipleScattering::physStepLimit
private

Definition at line 288 of file G4VMultipleScattering.hh.

◆ safetyHelper

G4SafetyHelper* G4VMultipleScattering::safetyHelper
private

Definition at line 259 of file G4VMultipleScattering.hh.

◆ stepLimit

G4MscStepLimitType G4VMultipleScattering::stepLimit
private

Definition at line 267 of file G4VMultipleScattering.hh.

◆ theParameters

G4EmParameters* G4VMultipleScattering::theParameters
private

Definition at line 255 of file G4VMultipleScattering.hh.

◆ tPathLength

G4double G4VMultipleScattering::tPathLength
private

Definition at line 289 of file G4VMultipleScattering.hh.


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