Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleChangeForMSC Class Referencefinal

#include <G4ParticleChangeForMSC.hh>

+ Inheritance diagram for G4ParticleChangeForMSC:

Public Member Functions

 G4ParticleChangeForMSC ()
 
 ~G4ParticleChangeForMSC () override=default
 
 G4ParticleChangeForMSC (const G4ParticleChangeForMSC &right)=delete
 
G4ParticleChangeForMSCoperator= (const G4ParticleChangeForMSC &right)=delete
 
G4StepUpdateStepForAlongStep (G4Step *step) final
 
void InitialiseMSC (const G4Track &, const G4Step &step)
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetProposedMomentumDirection () const
 
void ProposePosition (const G4ThreeVector &finalPosition)
 
const G4ThreeVectorGetProposedPosition () const
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()=default
 
 G4VParticleChange (const G4VParticleChange &right)=delete
 
G4VParticleChangeoperator= (const G4VParticleChange &right)=delete
 
virtual G4StepUpdateStepForAtRest (G4Step *Step)
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
virtual G4StepUpdateStepForPostStep (G4Step *Step)
 
virtual void Initialize (const G4Track &)
 
G4double GetTrueStepLength () const
 
void ProposeTrueStepLength (G4double truePathLength)
 
G4double GetLocalEnergyDeposit () const
 
void ProposeLocalEnergyDeposit (G4double anEnergyPart)
 
G4double GetNonIonizingEnergyDeposit () const
 
void ProposeNonIonizingEnergyDeposit (G4double anEnergyPart)
 
G4TrackStatus GetTrackStatus () const
 
void ProposeTrackStatus (G4TrackStatus status)
 
const G4TrackGetCurrentTrack () const
 
G4SteppingControl GetSteppingControl () const
 
void ProposeSteppingControl (G4SteppingControl StepControlFlag)
 
G4bool GetFirstStepInVolume () const
 
G4bool GetLastStepInVolume () const
 
void ProposeFirstStepInVolume (G4bool flag)
 
void ProposeLastStepInVolume (G4bool flag)
 
void Clear ()
 
void SetNumberOfSecondaries (G4int totSecondaries)
 
G4int GetNumberOfSecondaries () const
 
G4TrackGetSecondary (G4int anIndex) const
 
void AddSecondary (G4Track *aSecondary)
 
G4double GetWeight () const
 
G4double GetParentWeight () const
 
void ProposeWeight (G4double finalWeight)
 
void ProposeParentWeight (G4double finalWeight)
 
void SetSecondaryWeightByProcess (G4bool)
 
G4bool IsSecondaryWeightSetByProcess () const
 
void SetParentWeightByProcess (G4bool)
 
G4bool IsParentWeightSetByProcess () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int vLevel)
 
G4int GetVerboseLevel () const
 
virtual G4bool CheckIt (const G4Track &)
 
void ClearDebugFlag ()
 
void SetDebugFlag ()
 
G4bool GetDebugFlag () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VParticleChange
G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeLocalEnergyDeposit ()
 
void InitializeSteppingControl ()
 
void InitializeParentWeight (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries ()
 
void InitializeFromStep (const G4Step *)
 
G4double ComputeBeta (G4double kinEnergy)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 
- Protected Attributes inherited from G4VParticleChange
const G4TracktheCurrentTrack = nullptr
 
std::vector< G4Track * > theListOfSecondaries
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4double theTrueStepLength = 0.0
 
G4double theParentWeight = 1.0
 
G4double theParentGlobalTime = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4int theSizeOftheListOfSecondaries = 0
 
G4int verboseLevel = 1
 
G4int nError = 0
 
G4bool theFirstStepInVolume = false
 
G4bool theLastStepInVolume = false
 
G4bool isParentWeightProposed = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool debugFlag = false
 
- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 
static const G4int maxError = 10
 

Detailed Description

Definition at line 44 of file G4ParticleChangeForMSC.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForMSC() [1/2]

G4ParticleChangeForMSC::G4ParticleChangeForMSC ( )

Definition at line 40 of file G4ParticleChangeForMSC.cc.

41{
42 // Disable flag that is enabled in G4VParticleChange if G4VERBOSE.
43 debugFlag = false;
44}

◆ ~G4ParticleChangeForMSC()

G4ParticleChangeForMSC::~G4ParticleChangeForMSC ( )
overridedefault

◆ G4ParticleChangeForMSC() [2/2]

G4ParticleChangeForMSC::G4ParticleChangeForMSC ( const G4ParticleChangeForMSC right)
delete

Member Function Documentation

◆ GetProposedMomentumDirection()

const G4ThreeVector * G4ParticleChangeForMSC::GetProposedMomentumDirection ( ) const
inline

Definition at line 85 of file G4ParticleChangeForMSC.hh.

86{
87 return &theMomentumDirection;
88}

◆ GetProposedPosition()

const G4ThreeVector & G4ParticleChangeForMSC::GetProposedPosition ( ) const
inline

Definition at line 95 of file G4ParticleChangeForMSC.hh.

96{
97 return thePosition;
98}

Referenced by G4VMultipleScattering::AlongStepDoIt().

◆ InitialiseMSC()

void G4ParticleChangeForMSC::InitialiseMSC ( const G4Track track,
const G4Step step 
)
inline

Definition at line 101 of file G4ParticleChangeForMSC.hh.

102{
104 auto poststep = step.GetPostStepPoint();
105 thePosition = poststep->GetPosition();
106 theMomentumDirection = poststep->GetMomentumDirection();
107 theCurrentTrack = &track;
108}
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
G4TrackStatus theStatusChange
const G4Track * theCurrentTrack

Referenced by G4VMultipleScattering::AlongStepDoIt().

◆ operator=()

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

◆ ProposeMomentumDirection()

void G4ParticleChangeForMSC::ProposeMomentumDirection ( const G4ThreeVector Pfinal)
inline

◆ ProposePosition()

void G4ParticleChangeForMSC::ProposePosition ( const G4ThreeVector finalPosition)
inline

Definition at line 90 of file G4ParticleChangeForMSC.hh.

91{
92 thePosition = P;
93}

Referenced by G4VMultipleScattering::AlongStepDoIt().

◆ UpdateStepForAlongStep()

G4Step * G4ParticleChangeForMSC::UpdateStepForAlongStep ( G4Step step)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 47 of file G4ParticleChangeForMSC.cc.

48{
49 // update the G4Step specific attributes
50 pStep->SetStepLength(theTrueStepLength);
51
52 // multiple scattering calculates the final state of the particle
53 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
54
55 // update momentum direction
56 pPostStepPoint->SetMomentumDirection(theMomentumDirection);
57
58 // update position
59 pPostStepPoint->SetPosition(thePosition);
60 return pStep;
61}
void SetPosition(const G4ThreeVector &aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)

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