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

#include <G4ParticleChangeForLoss.hh>

+ Inheritance diagram for G4ParticleChangeForLoss:

Public Member Functions

 G4ParticleChangeForLoss ()
 
virtual ~G4ParticleChangeForLoss ()
 
G4StepUpdateStepForAlongStep (G4Step *Step)
 
G4StepUpdateStepForPostStep (G4Step *Step)
 
void InitializeForAlongStep (const G4Track &)
 
void InitializeForPostStep (const G4Track &)
 
G4double GetProposedCharge () const
 
void SetProposedCharge (G4double theCharge)
 
G4double GetCharge () const
 
void ProposeCharge (G4double finalCharge)
 
G4double GetProposedKineticEnergy () const
 
void SetProposedKineticEnergy (G4double proposedKinEnergy)
 
const G4ThreeVectorGetProposedMomentumDirection () const
 
void SetProposedMomentumDirection (const G4ThreeVector &dir)
 
const G4ThreeVectorGetMomentumDirection () const
 
void ProposeMomentumDirection (G4double Px, G4double Py, G4double Pz)
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetProposedPolarization () const
 
void ProposePolarization (const G4ThreeVector &dir)
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
const G4TrackGetCurrentTrack () const
 
void SetLowEnergyLimit (G4double elimit)
 
virtual void DumpInfo () const
 
virtual G4bool CheckIt (const G4Track &)
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()
 
G4bool operator== (const G4VParticleChange &right) const
 
G4bool operator!= (const G4VParticleChange &right) const
 
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)
 
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
 

Protected Member Functions

 G4ParticleChangeForLoss (const G4ParticleChangeForLoss &right)
 
G4ParticleChangeForLossoperator= (const G4ParticleChangeForLoss &right)
 
- Protected Member Functions inherited from G4VParticleChange
 G4VParticleChange (const G4VParticleChange &right)
 
G4VParticleChangeoperator= (const G4VParticleChange &right)
 
G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeTrueStepLength (const G4Track &)
 
void InitializeLocalEnergyDeposit (const G4Track &)
 
void InitializeSteppingControl (const G4Track &)
 
void InitializeParentWeight (const G4Track &)
 
void InitializeParentGlobalTime (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries (const G4Track &)
 
void InitializeStepInVolumeFlags (const G4Track &)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VParticleChange
G4TrackFastVectortheListOfSecondaries = nullptr
 
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
 
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
 

Detailed Description

Definition at line 44 of file G4ParticleChangeForLoss.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForLoss() [1/2]

G4ParticleChangeForLoss::G4ParticleChangeForLoss ( )

Definition at line 40 of file G4ParticleChangeForLoss.cc.

42 , lowEnergyLimit(1.0 * eV)
43{
45 debugFlag = false;
46}
@ NormalCondition
G4SteppingControl theSteppingControlFlag

◆ ~G4ParticleChangeForLoss()

G4ParticleChangeForLoss::~G4ParticleChangeForLoss ( )
virtual

Definition at line 49 of file G4ParticleChangeForLoss.cc.

50{
51}

◆ G4ParticleChangeForLoss() [2/2]

G4ParticleChangeForLoss::G4ParticleChangeForLoss ( const G4ParticleChangeForLoss right)
protected

Definition at line 54 of file G4ParticleChangeForLoss.cc.

56 : G4VParticleChange(right)
57{
58 currentTrack = right.currentTrack;
59 proposedKinEnergy = right.proposedKinEnergy;
60 lowEnergyLimit = right.lowEnergyLimit;
61 currentCharge = right.currentCharge;
62 proposedMomentumDirection = right.proposedMomentumDirection;
63}

Member Function Documentation

◆ CheckIt()

G4bool G4ParticleChangeForLoss::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 124 of file G4ParticleChangeForLoss.cc.

125{
126 G4bool itsOK = true;
127 G4bool exitWithError = false;
128
129 G4double accuracy;
130
131 // Energy should not be lager than initial value
132 accuracy = (proposedKinEnergy - aTrack.GetKineticEnergy()) / MeV;
133 if(accuracy > accuracyForWarning)
134 {
135 itsOK = false;
136 exitWithError = (accuracy > accuracyForException);
137#ifdef G4VERBOSE
138 G4cout << "G4ParticleChangeForLoss::CheckIt: ";
139 G4cout << "KinEnergy become larger than the initial value!"
140 << " Difference: " << accuracy << "[MeV] " << G4endl;
142 << " E=" << aTrack.GetKineticEnergy() / MeV
143 << " pos=" << aTrack.GetPosition().x() / m << ", "
144 << aTrack.GetPosition().y() / m << ", "
145 << aTrack.GetPosition().z() / m << G4endl;
146#endif
147 }
148
149 // dump out information of this particle change
150#ifdef G4VERBOSE
151 if(!itsOK) { DumpInfo(); }
152#endif
153
154 // exit with error
155 if(exitWithError)
156 {
157 G4Exception("G4ParticleChangeForLoss::CheckIt()", "TRACK004",
158 EventMustBeAborted, "energy was illegal");
159 }
160
161 // correction
162 if(!itsOK)
163 {
164 proposedKinEnergy = aTrack.GetKineticEnergy();
165 }
166
167 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
168 return itsOK;
169}
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
const G4String & GetParticleName() const
const G4ThreeVector & GetPosition() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
virtual G4bool CheckIt(const G4Track &)
static const G4double accuracyForException
static const G4double accuracyForWarning

◆ DumpInfo()

void G4ParticleChangeForLoss::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 104 of file G4ParticleChangeForLoss.cc.

105{
106 // use base-class DumpInfo
108
109 G4int oldprc = G4cout.precision(3);
110 G4cout << " Charge (eplus) : " << std::setw(20)
111 << currentCharge / eplus << G4endl;
112 G4cout << " Kinetic Energy (MeV): " << std::setw(20)
113 << proposedKinEnergy / MeV << G4endl;
114 G4cout << " Momentum Direct - x : " << std::setw(20)
115 << proposedMomentumDirection.x() << G4endl;
116 G4cout << " Momentum Direct - y : " << std::setw(20)
117 << proposedMomentumDirection.y() << G4endl;
118 G4cout << " Momentum Direct - z : " << std::setw(20)
119 << proposedMomentumDirection.z() << G4endl;
120 G4cout.precision(oldprc);
121}
int G4int
Definition: G4Types.hh:85
virtual void DumpInfo() const

Referenced by CheckIt().

◆ GetCharge()

G4double G4ParticleChangeForLoss::GetCharge ( ) const
inline

Definition at line 147 of file G4ParticleChangeForLoss.hh.

148{
149 return currentCharge;
150}

◆ GetCurrentTrack()

const G4Track * G4ParticleChangeForLoss::GetCurrentTrack ( ) const
inline

Definition at line 201 of file G4ParticleChangeForLoss.hh.

202{
203 return currentTrack;
204}

Referenced by G4PolarizedMollerBhabhaModel::SampleSecondaries().

◆ GetMomentumDirection()

const G4ThreeVector & G4ParticleChangeForLoss::GetMomentumDirection ( ) const
inline

Definition at line 172 of file G4ParticleChangeForLoss.hh.

173{
174 return proposedMomentumDirection;
175}

◆ GetProposedCharge()

G4double G4ParticleChangeForLoss::GetProposedCharge ( ) const
inline

Definition at line 141 of file G4ParticleChangeForLoss.hh.

142{
143 return currentCharge;
144}

◆ GetProposedKineticEnergy()

G4double G4ParticleChangeForLoss::GetProposedKineticEnergy ( ) const
inline

Definition at line 129 of file G4ParticleChangeForLoss.hh.

130{
131 return proposedKinEnergy;
132}

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

◆ GetProposedMomentumDirection()

const G4ThreeVector & G4ParticleChangeForLoss::GetProposedMomentumDirection ( ) const
inline

Definition at line 166 of file G4ParticleChangeForLoss.hh.

167{
168 return proposedMomentumDirection;
169}

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

◆ GetProposedPolarization()

const G4ThreeVector & G4ParticleChangeForLoss::GetProposedPolarization ( ) const
inline

Definition at line 207 of file G4ParticleChangeForLoss.hh.

208{
209 return proposedPolarization;
210}

◆ InitializeForAlongStep()

void G4ParticleChangeForLoss::InitializeForAlongStep ( const G4Track track)
inline

Definition at line 229 of file G4ParticleChangeForLoss.hh.

230{
235 theParentWeight = track.GetWeight();
236 // isParentWeightProposed = false;
237 proposedKinEnergy = track.GetKineticEnergy();
238 currentCharge = track.GetDynamicParticle()->GetCharge();
239}
G4double GetCharge() const
G4TrackStatus GetTrackStatus() const
G4double GetWeight() const
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit
void InitializeSecondaries(const G4Track &)

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

◆ InitializeForPostStep()

void G4ParticleChangeForLoss::InitializeForPostStep ( const G4Track track)
inline

Definition at line 242 of file G4ParticleChangeForLoss.hh.

243{
248 theParentWeight = track.GetWeight();
249 // isParentWeightProposed = false;
250 proposedKinEnergy = track.GetKineticEnergy();
251 currentCharge = track.GetDynamicParticle()->GetCharge();
252 proposedMomentumDirection = track.GetMomentumDirection();
253 proposedPolarization = track.GetPolarization();
254 currentTrack = &track;
255}
const G4ThreeVector & GetMomentumDirection() const
const G4ThreeVector & GetPolarization() const

Referenced by G4VEnergyLossProcess::PostStepDoIt().

◆ operator=()

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

Definition at line 66 of file G4ParticleChangeForLoss.cc.

68{
69 if(this != &right)
70 {
72 {
73 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
74 {
75 if((*theListOfSecondaries)[index])
76 delete(*theListOfSecondaries)[index];
77 }
78 }
82 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
83 {
84 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index]));
85 theListOfSecondaries->SetElement(index, newTrack);
86 }
87
94
95 currentTrack = right.currentTrack;
96 proposedKinEnergy = right.proposedKinEnergy;
97 currentCharge = right.currentCharge;
98 proposedMomentumDirection = right.proposedMomentumDirection;
99 }
100 return *this;
101}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
G4TrackFastVector * theListOfSecondaries

◆ ProposeCharge()

void G4ParticleChangeForLoss::ProposeCharge ( G4double  finalCharge)
inline

Definition at line 159 of file G4ParticleChangeForLoss.hh.

160{
161 currentCharge = theCharge;
162}

◆ ProposeMomentumDirection() [1/2]

void G4ParticleChangeForLoss::ProposeMomentumDirection ( const G4ThreeVector Pfinal)
inline

Definition at line 178 of file G4ParticleChangeForLoss.hh.

179{
180 proposedMomentumDirection = dir;
181}

◆ ProposeMomentumDirection() [2/2]

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

Definition at line 191 of file G4ParticleChangeForLoss.hh.

194{
195 proposedMomentumDirection.setX(Px);
196 proposedMomentumDirection.setY(Py);
197 proposedMomentumDirection.setZ(Pz);
198}
void setY(double)
void setZ(double)
void setX(double)

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), G4LivermoreIonisationModel::SampleSecondaries(), G4PenelopeBremsstrahlungModel::SampleSecondaries(), and G4PenelopeIonisationModel::SampleSecondaries().

◆ ProposePolarization() [1/2]

void G4ParticleChangeForLoss::ProposePolarization ( const G4ThreeVector dir)
inline

Definition at line 213 of file G4ParticleChangeForLoss.hh.

214{
215 proposedPolarization = dir;
216}

Referenced by G4ePolarizedBremsstrahlungModel::SampleSecondaries(), and G4PolarizedMollerBhabhaModel::SampleSecondaries().

◆ ProposePolarization() [2/2]

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

Definition at line 219 of file G4ParticleChangeForLoss.hh.

222{
223 proposedPolarization.setX(Px);
224 proposedPolarization.setY(Py);
225 proposedPolarization.setZ(Pz);
226}

◆ SetLowEnergyLimit()

void G4ParticleChangeForLoss::SetLowEnergyLimit ( G4double  elimit)
inline

Definition at line 258 of file G4ParticleChangeForLoss.hh.

259{
260 lowEnergyLimit = elimit;
261}

◆ SetProposedCharge()

void G4ParticleChangeForLoss::SetProposedCharge ( G4double  theCharge)
inline

Definition at line 153 of file G4ParticleChangeForLoss.hh.

154{
155 currentCharge = theCharge;
156}

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

◆ SetProposedKineticEnergy()

◆ SetProposedMomentumDirection()

◆ UpdateStepForAlongStep()

G4Step * G4ParticleChangeForLoss::UpdateStepForAlongStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 172 of file G4ParticleChangeForLoss.cc.

173{
174 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
175 G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
176 G4Track* pTrack = pStep->GetTrack();
177
178 // accumulate change of the kinetic energy
179 G4double preKinEnergy = pPreStepPoint->GetKineticEnergy();
180 G4double kinEnergy =
181 pPostStepPoint->GetKineticEnergy() + (proposedKinEnergy - preKinEnergy);
182
183 // update kinetic energy and charge
184 if(kinEnergy < lowEnergyLimit)
185 {
186 theLocalEnergyDeposit += kinEnergy;
187 kinEnergy = 0.0;
188 pPostStepPoint->SetVelocity(0.0);
189 }
190 else
191 {
192 pPostStepPoint->SetCharge(currentCharge);
193 // calculate velocity
194 pTrack->SetKineticEnergy(kinEnergy);
195 pPostStepPoint->SetVelocity(pTrack->CalculateVelocity());
196 pTrack->SetKineticEnergy(preKinEnergy);
197 }
198 pPostStepPoint->SetKineticEnergy(kinEnergy);
199
201 {
202 pPostStepPoint->SetWeight(theParentWeight);
203 }
204
205 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
206 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
207 return pStep;
208}
void SetKineticEnergy(const G4double aValue)
void SetWeight(G4double aValue)
void SetCharge(G4double value)
void SetVelocity(G4double v)
G4double GetKineticEnergy() const
G4double CalculateVelocity() const
void SetKineticEnergy(const G4double aValue)

◆ UpdateStepForPostStep()

G4Step * G4ParticleChangeForLoss::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 211 of file G4ParticleChangeForLoss.cc.

212{
213 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
214 G4Track* pTrack = pStep->GetTrack();
215
216 pPostStepPoint->SetCharge(currentCharge);
217 pPostStepPoint->SetMomentumDirection(proposedMomentumDirection);
218 pPostStepPoint->SetKineticEnergy(proposedKinEnergy);
219 pTrack->SetKineticEnergy(proposedKinEnergy);
220 if(proposedKinEnergy > 0.0)
221 {
222 pPostStepPoint->SetVelocity(pTrack->CalculateVelocity());
223 }
224 else
225 {
226 pPostStepPoint->SetVelocity(0.0);
227 }
228 pPostStepPoint->SetPolarization(proposedPolarization);
229
231 {
232 pPostStepPoint->SetWeight(theParentWeight);
233 }
234
235 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
236 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
237 return pStep;
238}
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)

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