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

#include <G4ParticleChangeForLoss.hh>

+ Inheritance diagram for G4ParticleChangeForLoss:

Public Member Functions

 G4ParticleChangeForLoss ()
 
 ~G4ParticleChangeForLoss () override=default
 
 G4ParticleChangeForLoss (const G4ParticleChangeForLoss &right)=delete
 
G4ParticleChangeForLossoperator= (const G4ParticleChangeForLoss &right)=delete
 
G4StepUpdateStepForAlongStep (G4Step *step) final
 
G4StepUpdateStepForPostStep (G4Step *step) final
 
void InitializeForAlongStep (const G4Track &)
 
void InitializeForPostStep (const G4Track &)
 
G4double GetProposedCharge () const
 
void SetProposedCharge (G4double theCharge)
 
G4double GetProposedKineticEnergy () const
 
void SetProposedKineticEnergy (G4double proposedKinEnergy)
 
const G4ThreeVectorGetProposedMomentumDirection () const
 
void SetProposedMomentumDirection (const G4ThreeVector &dir)
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetProposedPolarization () const
 
void ProposePolarization (const G4ThreeVector &dir)
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
void DumpInfo () const final
 
- 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 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
 
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 42 of file G4ParticleChangeForLoss.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForLoss() [1/2]

G4ParticleChangeForLoss::G4ParticleChangeForLoss ( )

Definition at line 37 of file G4ParticleChangeForLoss.cc.

38{
39 // Disable flag that is enabled in G4VParticleChange if G4VERBOSE.
40 debugFlag = false;
41}

◆ ~G4ParticleChangeForLoss()

G4ParticleChangeForLoss::~G4ParticleChangeForLoss ( )
overridedefault

◆ G4ParticleChangeForLoss() [2/2]

G4ParticleChangeForLoss::G4ParticleChangeForLoss ( const G4ParticleChangeForLoss & right)
delete

Member Function Documentation

◆ DumpInfo()

void G4ParticleChangeForLoss::DumpInfo ( ) const
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 44 of file G4ParticleChangeForLoss.cc.

45{
46 // use base-class DumpInfo
48
49 G4long oldprc = G4cout.precision(8);
50 G4cout << " -----------------------------------------------" << G4endl;
51 G4cout << " G4ParticleChangeForLoss proposes: " << G4endl;
52 G4cout << " Charge (eplus) : " << std::setw(20)
53 << currentCharge / eplus << G4endl;
54 G4cout << " Kinetic Energy (MeV): " << std::setw(20)
55 << proposedKinEnergy / MeV << G4endl;
56 G4cout << " Momentum Direct - x : " << std::setw(20)
57 << proposedMomentumDirection.x() << G4endl;
58 G4cout << " Momentum Direct - y : " << std::setw(20)
59 << proposedMomentumDirection.y() << G4endl;
60 G4cout << " Momentum Direct - z : " << std::setw(20)
61 << proposedMomentumDirection.z() << G4endl;
62 G4cout.precision(oldprc);
63}
long G4long
Definition G4Types.hh:87
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
virtual void DumpInfo() const

◆ GetProposedCharge()

G4double G4ParticleChangeForLoss::GetProposedCharge ( ) const
inline

Definition at line 113 of file G4ParticleChangeForLoss.hh.

114{
115 return currentCharge;
116}

◆ GetProposedKineticEnergy()

G4double G4ParticleChangeForLoss::GetProposedKineticEnergy ( ) const
inline

Definition at line 102 of file G4ParticleChangeForLoss.hh.

103{
104 return proposedKinEnergy;
105}

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), and G4VEnergyLossProcess::PostStepDoIt().

◆ GetProposedMomentumDirection()

const G4ThreeVector & G4ParticleChangeForLoss::GetProposedMomentumDirection ( ) const
inline

Definition at line 126 of file G4ParticleChangeForLoss.hh.

127{
128 return proposedMomentumDirection;
129}

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), and G4PolarizedBremsstrahlungModel::SampleSecondaries().

◆ GetProposedPolarization()

const G4ThreeVector & G4ParticleChangeForLoss::GetProposedPolarization ( ) const
inline

Definition at line 145 of file G4ParticleChangeForLoss.hh.

146{
147 return proposedPolarization;
148}

◆ InitializeForAlongStep()

void G4ParticleChangeForLoss::InitializeForAlongStep ( const G4Track & track)
inline

Definition at line 164 of file G4ParticleChangeForLoss.hh.

165{
170 proposedKinEnergy = track.GetKineticEnergy();
171 currentCharge = track.GetDynamicParticle()->GetCharge();
172}
G4double GetCharge() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
void InitializeLocalEnergyDeposit()
void InitializeStatusChange(const G4Track &)
void InitializeSecondaries()
void InitializeParentWeight(const G4Track &)

Referenced by G4NuclearStopping::AlongStepDoIt(), G4VEnergyLossProcess::AlongStepDoIt(), and InitializeForPostStep().

◆ InitializeForPostStep()

void G4ParticleChangeForLoss::InitializeForPostStep ( const G4Track & track)
inline

Definition at line 175 of file G4ParticleChangeForLoss.hh.

176{
178 proposedMomentumDirection = track.GetMomentumDirection();
179 proposedPolarization = track.GetPolarization();
180}
void InitializeForAlongStep(const G4Track &)
const G4ThreeVector & GetMomentumDirection() const
const G4ThreeVector & GetPolarization() const

Referenced by G4VEnergyLossProcess::PostStepDoIt().

◆ operator=()

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

◆ ProposeMomentumDirection()

void G4ParticleChangeForLoss::ProposeMomentumDirection ( const G4ThreeVector & Pfinal)
inline

◆ ProposePolarization() [1/2]

void G4ParticleChangeForLoss::ProposePolarization ( const G4ThreeVector & dir)
inline

Definition at line 151 of file G4ParticleChangeForLoss.hh.

152{
153 proposedPolarization = dir;
154}

Referenced by G4PolarizedBremsstrahlungModel::SampleSecondaries(), and G4PolarizedIonisationModel::SampleSecondaries().

◆ ProposePolarization() [2/2]

void G4ParticleChangeForLoss::ProposePolarization ( G4double Px,
G4double Py,
G4double Pz )
inline

Definition at line 156 of file G4ParticleChangeForLoss.hh.

159{
160 proposedPolarization.set(Px, Py, Pz);
161}
void set(double x, double y, double z)

◆ SetProposedCharge()

void G4ParticleChangeForLoss::SetProposedCharge ( G4double theCharge)
inline

Definition at line 119 of file G4ParticleChangeForLoss.hh.

120{
121 currentCharge = theCharge;
122}

Referenced by G4NuclearStopping::AlongStepDoIt(), and G4VEnergyLossProcess::AlongStepDoIt().

◆ SetProposedKineticEnergy()

◆ SetProposedMomentumDirection()

◆ UpdateStepForAlongStep()

G4Step * G4ParticleChangeForLoss::UpdateStepForAlongStep ( G4Step * step)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 66 of file G4ParticleChangeForLoss.cc.

67{
68 const G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
69 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
70
71 // accumulate change of the kinetic energy
72 G4double preKinEnergy = pPreStepPoint->GetKineticEnergy();
73 G4double kinEnergy =
74 pPostStepPoint->GetKineticEnergy() + (proposedKinEnergy - preKinEnergy);
75
76 pPostStepPoint->SetCharge(currentCharge);
77
78 // calculate velocity
79 if(kinEnergy > 0.0)
80 {
81 pPostStepPoint->SetKineticEnergy(kinEnergy);
82
83 // assuming that mass>0, zero mass particles do not have energy loss
84 pPostStepPoint->SetVelocity(CLHEP::c_light*ComputeBeta(kinEnergy));
85 }
86 else
87 {
88 pPostStepPoint->SetKineticEnergy(0.0);
89 pPostStepPoint->SetVelocity(0.0);
90 }
91
93 {
94 pPostStepPoint->SetWeight(theParentWeight);
95 }
96
97 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
98 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
99 return pStep;
100}
double G4double
Definition G4Types.hh:83
void SetKineticEnergy(const G4double aValue)
void SetWeight(G4double aValue)
void SetCharge(G4double value)
void SetVelocity(G4double v)
G4double GetKineticEnergy() const
G4double ComputeBeta(G4double kinEnergy)
G4double theNonIonizingEnergyDeposit

◆ UpdateStepForPostStep()

G4Step * G4ParticleChangeForLoss::UpdateStepForPostStep ( G4Step * step)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 103 of file G4ParticleChangeForLoss.cc.

104{
105 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
106
107 pPostStepPoint->SetCharge(currentCharge);
108 pPostStepPoint->SetMomentumDirection(proposedMomentumDirection);
109 if(proposedKinEnergy > 0.0)
110 {
111 pPostStepPoint->SetKineticEnergy(proposedKinEnergy);
112
113 // assuming that mass>0, zero mass particles do not have energy loss
114 pPostStepPoint->SetVelocity(CLHEP::c_light*ComputeBeta(proposedKinEnergy));
115 }
116 else
117 {
118 pPostStepPoint->SetKineticEnergy(0.0);
119 pPostStepPoint->SetVelocity(0.0);
120 }
121 pPostStepPoint->SetPolarization(proposedPolarization);
122
124 {
125 pPostStepPoint->SetWeight(theParentWeight);
126 }
127
128 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
129 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
130 return pStep;
131}
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)

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