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

#include <G4DynamicParticleMSC.hh>

Inheritance diagram for G4DynamicParticleMSC:

Public Member Functions

 G4DynamicParticleMSC ()
 ~G4DynamicParticleMSC () override
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition) override
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
void ProcessDescription (std::ostream &) const override
G4DynamicParticleMSCoperator= (const G4DynamicParticleMSC &right)=delete
 G4DynamicParticleMSC (const G4DynamicParticleMSC &)=delete
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)
virtual G4bool IsApplicable (const G4ParticleDefinition &)
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
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 StartTracking (G4Track *)
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 &)

Additional Inherited Members

Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
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 ()
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 58 of file G4DynamicParticleMSC.hh.

Constructor & Destructor Documentation

◆ G4DynamicParticleMSC() [1/2]

G4DynamicParticleMSC::G4DynamicParticleMSC ( )

Definition at line 60 of file G4DynamicParticleMSC.cc.

61 : G4VContinuousDiscreteProcess("dynPartMSC")
62{
65 lManager = G4LossTableManager::Instance();
66 lManager->Register(this);
67}
@ fDynamicMultipleScattering
static G4LossTableManager * Instance()
G4VContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
void SetVerboseLevel(G4int value)
void SetProcessSubType(G4int)

Referenced by G4DynamicParticleMSC(), and operator=().

◆ ~G4DynamicParticleMSC()

G4DynamicParticleMSC::~G4DynamicParticleMSC ( )
override

Definition at line 71 of file G4DynamicParticleMSC.cc.

72{
73 lManager->DeRegister(this);
74}

◆ G4DynamicParticleMSC() [2/2]

G4DynamicParticleMSC::G4DynamicParticleMSC ( const G4DynamicParticleMSC & )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4DynamicParticleMSC::AlongStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 114 of file G4DynamicParticleMSC.cc.

116{
117 fParticleChange.InitialiseMSC(track, step);
118
119 // no energy loss
120 if (fCharge == 0.0) { return &fParticleChange; }
121
122 G4double geomLength = step.GetStepLength();
123 G4double y = geomLength/fMaterial->GetRadlen();
124 G4double theta0 = c_highland*std::abs(fCharge)*std::sqrt(y)*
125 (1.0 + 0.038*G4Log(y*fCharge*fCharge/(fBeta*fBeta)))/fBeta;
126
127 if (theta0 < 0.001) { return &fParticleChange; }
128 G4double cost = 1.0;
130 if (theta0 < 1.0) {
131 G4double theta2 = theta0*theta0;
132 cost -= theta2*G4Log(1.0 + r*(G4Exp(2.0/theta2) - 1.0));
133 } else {
134 cost -= 2.0*r;
135 }
136 G4double phi = CLHEP::twopi*G4UniformRand();
137 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
138 fNewDir.set(sint*std::cos(phi), sint*std::sin(phi), cost);
139 fNewDir.rotateUz(step.GetPostStepPoint()->GetMomentumDirection());
140
141 fParticleChange.ProposeMomentumDirection(fNewDir);
142 fParticleChange.ProposeTrueStepLength(geomLength);
143 return &fParticleChange;
144}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4double G4Log(G4double x)
Definition G4Log.hh:169
double G4double
Definition G4Types.hh:83
#define G4UniformRand()
Definition Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const

◆ AlongStepGetPhysicalInteractionLength()

G4double G4DynamicParticleMSC::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety,
G4GPILSelection * selection )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 91 of file G4DynamicParticleMSC.cc.

94{
95 *selection = CandidateForSelection;
96 PreStepInitialisation(track);
97
98 // no step limit for the time being
99 return DBL_MAX;
100}
@ CandidateForSelection
#define DBL_MAX
Definition templates.hh:62

Referenced by GetContinuousStepLimit().

◆ GetContinuousStepLimit()

G4double G4DynamicParticleMSC::GetContinuousStepLimit ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimalStep,
G4double & currentSafety )
overridevirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 157 of file G4DynamicParticleMSC.cc.

161{
163 G4double x = AlongStepGetPhysicalInteractionLength(track, previousStepSize,
164 currentMinimalStep,
165 currentSafety, &selection);
166 return x;
167}
G4GPILSelection
@ NotCandidateForSelection
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override

◆ GetMeanFreePath()

G4double G4DynamicParticleMSC::GetMeanFreePath ( const G4Track & ,
G4double ,
G4ForceCondition * condition )
overridevirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 148 of file G4DynamicParticleMSC.cc.

150{
151 *condition = Forced;
152 return DBL_MAX;
153}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ operator=()

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

◆ PostStepGetPhysicalInteractionLength()

G4double G4DynamicParticleMSC::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 104 of file G4DynamicParticleMSC.cc.

107{
109 return DBL_MAX;
110}
@ NotForced

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 171 of file G4DynamicParticleMSC.cc.

172{
173 out << "G4DynamicParticleMSC: no delta rays" << G4endl;
174}
#define G4endl
Definition G4ios.hh:67

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