Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eMultipleScattering Class Reference

#include <G4eMultipleScattering.hh>

+ Inheritance diagram for G4eMultipleScattering:

Public Member Functions

 G4eMultipleScattering (const G4String &processName="msc")
 
 ~G4eMultipleScattering () override
 
G4bool IsApplicable (const G4ParticleDefinition &p) final
 
void ProcessDescription (std::ostream &) const override
 
G4eMultipleScatteringoperator= (const G4eMultipleScattering &right)=delete
 
 G4eMultipleScattering (const G4eMultipleScattering &)=delete
 
- Public Member Functions inherited from G4VMultipleScattering
 G4VMultipleScattering (const G4String &name="msc", G4ProcessType type=fElectromagnetic)
 
 ~G4VMultipleScattering () override
 
void ProcessDescription (std::ostream &outFile) const override
 
void StreamInfo (std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
void StartTracking (G4Track *) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4double ContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
 
 G4VMultipleScattering (G4VMultipleScattering &)=delete
 
G4VMultipleScatteringoperator= (const G4VMultipleScattering &right)=delete
 
G4VEmModelSelectModel (G4double kinEnergy, size_t idx)
 
void AddEmModel (G4int order, G4VMscModel *, const G4Region *region=nullptr)
 
void SetEmModel (G4VMscModel *, G4int idx=0)
 
G4VMscModelEmModel (size_t index=0) const
 
G4int NumberOfModels () const
 
G4VMscModelGetModelByIndex (G4int idx, G4bool ver=false) const
 
G4bool LateralDisplasmentFlag () const
 
G4double Skin () const
 
G4double RangeFactor () const
 
G4double GeomFactor () const
 
G4double PolarAngleLimit () const
 
G4bool UseBaseMaterial () 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 ()
 
G4VContinuousDiscreteProcessoperator= (const G4VContinuousDiscreteProcess &)=delete
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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 const G4VProcessGetCreatorProcess () const
 
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

void StreamProcessInfo (std::ostream &outFile) const override
 
void InitialiseProcess (const G4ParticleDefinition *) override
 
- Protected Member Functions inherited from G4VMultipleScattering
G4double GetMeanFreePath (const G4Track &track, G4double, G4ForceCondition *condition) override
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VMultipleScattering
G4ParticleChangeForMSC fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 59 of file G4eMultipleScattering.hh.

Constructor & Destructor Documentation

◆ G4eMultipleScattering() [1/2]

G4eMultipleScattering::G4eMultipleScattering ( const G4String & processName = "msc")
explicit

Definition at line 54 of file G4eMultipleScattering.cc.

55 : G4VMultipleScattering(processName)
56{
57 isInitialized = false;
58}
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)

◆ ~G4eMultipleScattering()

G4eMultipleScattering::~G4eMultipleScattering ( )
overridedefault

◆ G4eMultipleScattering() [2/2]

G4eMultipleScattering::G4eMultipleScattering ( const G4eMultipleScattering & )
delete

Member Function Documentation

◆ InitialiseProcess()

void G4eMultipleScattering::InitialiseProcess ( const G4ParticleDefinition * )
overrideprotectedvirtual

Implements G4VMultipleScattering.

Definition at line 73 of file G4eMultipleScattering.cc.

74{
75 if(isInitialized) { return; }
76 if(nullptr == EmModel(0)) { SetEmModel(new G4UrbanMscModel()); }
77 AddEmModel(1, EmModel(0));
78 if(nullptr != EmModel(1)) { AddEmModel(1, EmModel(1)); }
79 isInitialized = true;
80}
void AddEmModel(G4int order, G4VMscModel *, const G4Region *region=nullptr)
void SetEmModel(G4VMscModel *, G4int idx=0)
G4VMscModel * EmModel(size_t index=0) const

◆ IsApplicable()

G4bool G4eMultipleScattering::IsApplicable ( const G4ParticleDefinition & p)
finalvirtual

Reimplemented from G4VProcess.

Definition at line 66 of file G4eMultipleScattering.cc.

67{
68 return (p.GetPDGCharge() != 0.0);
69}

◆ operator=()

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

◆ ProcessDescription()

void G4eMultipleScattering::ProcessDescription ( std::ostream & out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 97 of file G4eMultipleScattering.cc.

98{
99 out <<
100 " Multiple scattering. Simulates combined effects of elastic scattering\n"<<
101 " at the end of the step, to save computing time. May be combined with\n"<<
102 " Coulomb scattering in a 'mixed' scattering algorithm.";
104}
void ProcessDescription(std::ostream &outFile) const override

◆ StreamProcessInfo()

void G4eMultipleScattering::StreamProcessInfo ( std::ostream & outFile) const
overrideprotectedvirtual

Reimplemented from G4VMultipleScattering.

Definition at line 84 of file G4eMultipleScattering.cc.

85{
86 out << " RangeFactor= " << RangeFactor()
87 << ", stepLimType: " << StepLimitType()
88 << ", latDisp: " << LateralDisplasmentFlag();
90 out << ", skin= " << Skin() << ", geomFactor= " << GeomFactor();
91 }
92 out << G4endl;
93}
@ fUseDistanceToBoundary
#define G4endl
Definition G4ios.hh:67
G4MscStepLimitType StepLimitType() const

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