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

#include <G4ParticleChange.hh>

+ Inheritance diagram for G4ParticleChange:

Public Member Functions

 G4ParticleChange ()
 
 ~G4ParticleChange () override=default
 
 G4ParticleChange (const G4ParticleChange &right)=delete
 
G4ParticleChangeoperator= (const G4ParticleChange &right)=delete
 
G4StepUpdateStepForAlongStep (G4Step *Step) override
 
G4StepUpdateStepForAtRest (G4Step *Step) override
 
G4StepUpdateStepForPostStep (G4Step *Step) override
 
void Initialize (const G4Track &) override
 
const G4ThreeVectorGetMomentumDirection () const
 
void ProposeMomentumDirection (G4double Px, G4double Py, G4double Pz)
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetPolarization () const
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
void ProposePolarization (const G4ThreeVector &finalPoralization)
 
G4double GetEnergy () const
 
void ProposeEnergy (G4double finalEnergy)
 
G4double GetVelocity () const
 
void ProposeVelocity (G4double finalVelocity)
 
G4double GetProperTime () const
 
void ProposeProperTime (G4double finalProperTime)
 
const G4ThreeVectorGetPosition () const
 
void ProposePosition (G4double x, G4double y, G4double z)
 
void ProposePosition (const G4ThreeVector &finalPosition)
 
void ProposeGlobalTime (G4double t)
 
void ProposeLocalTime (G4double t)
 
G4double GetGlobalTime (G4double timeDelay=0.0) const
 
G4double GetLocalTime (G4double timeDelay=0.0) const
 
G4double GetMass () const
 
void ProposeMass (G4double finalMass)
 
G4double GetCharge () const
 
void ProposeCharge (G4double finalCharge)
 
G4double GetMagneticMoment () const
 
void ProposeMagneticMoment (G4double finalMagneticMoment)
 
G4ThreeVector GetGlobalPosition (const G4ThreeVector &displacement) const
 
G4ThreeVector CalcMomentum (G4double energy, G4ThreeVector direction, G4double mass) const
 
void AddSecondary (G4Track *aSecondary)
 
void AddSecondary (G4DynamicParticle *aSecondary, G4bool IsGoodForTracking=false)
 
void AddSecondary (G4DynamicParticle *aSecondary, G4ThreeVector position, G4bool IsGoodForTracking=false)
 
void AddSecondary (G4DynamicParticle *aSecondary, G4double time, G4bool IsGoodForTracking=false)
 
void DumpInfo () const override
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()=default
 
 G4VParticleChange (const G4VParticleChange &right)=delete
 
G4VParticleChangeoperator= (const G4VParticleChange &right)=delete
 
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
 

Protected Member Functions

G4StepUpdateStepInfo (G4Step *Step)
 
- 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

G4ThreeVector theMomentumDirectionChange
 
G4ThreeVector thePolarizationChange
 
G4double theEnergyChange = 0.0
 
G4double theVelocityChange = 0.0
 
G4bool isVelocityChanged = false
 
G4ThreeVector thePositionChange
 
G4double theGlobalTime0 = 0.0
 
G4double theLocalTime0 = 0.0
 
G4double theTimeChange = 0.0
 
G4double theProperTimeChange = 0.0
 
G4double theMassChange = 0.0
 
G4double theChargeChange = 0.0
 
G4double theMagneticMomentChange = 0.0
 
- 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
 

Additional Inherited Members

- 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 55 of file G4ParticleChange.hh.

Constructor & Destructor Documentation

◆ G4ParticleChange() [1/2]

G4ParticleChange::G4ParticleChange ( )

Definition at line 40 of file G4ParticleChange.cc.

41{}

◆ ~G4ParticleChange()

G4ParticleChange::~G4ParticleChange ( )
overridedefault

◆ G4ParticleChange() [2/2]

G4ParticleChange::G4ParticleChange ( const G4ParticleChange & right)
delete

Member Function Documentation

◆ AddSecondary() [1/4]

void G4ParticleChange::AddSecondary ( G4DynamicParticle * aSecondary,
G4bool IsGoodForTracking = false )

Definition at line 44 of file G4ParticleChange.cc.

46{
47 // create track
48 G4Track* aTrack = new G4Track(aParticle, GetGlobalTime(), thePositionChange);
49
50 // set IsGoodGorTrackingFlag
51 if(IsGoodForTracking)
52 aTrack->SetGoodForTrackingFlag();
53
54 // touchable handle is copied to keep the pointer
56
57 // add a secondary
59}
G4ThreeVector thePositionChange
G4double GetGlobalTime(G4double timeDelay=0.0) const
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
void SetGoodForTrackingFlag(G4bool value=true)
void AddSecondary(G4Track *aSecondary)
const G4Track * theCurrentTrack

◆ AddSecondary() [2/4]

void G4ParticleChange::AddSecondary ( G4DynamicParticle * aSecondary,
G4double time,
G4bool IsGoodForTracking = false )

Definition at line 81 of file G4ParticleChange.cc.

83{
84 // create track
85 G4Track* aTrack = new G4Track(aParticle, newTime, thePositionChange);
86
87 // set IsGoodGorTrackingFlag
88 if(IsGoodForTracking)
89 aTrack->SetGoodForTrackingFlag();
90
91 // touchable handle is copied to keep the pointer
93
94 // add a secondary
96}

◆ AddSecondary() [3/4]

void G4ParticleChange::AddSecondary ( G4DynamicParticle * aSecondary,
G4ThreeVector position,
G4bool IsGoodForTracking = false )

Definition at line 62 of file G4ParticleChange.cc.

65{
66 // create track
67 G4Track* aTrack = new G4Track(aParticle, GetGlobalTime(), newPosition);
68
69 // set IsGoodGorTrackingFlag
70 if(IsGoodForTracking)
71 aTrack->SetGoodForTrackingFlag();
72
73 // touchable is a temporary object, so you cannot keep the pointer
74 aTrack->SetTouchableHandle(nullptr);
75
76 // add a secondary
78}

◆ AddSecondary() [4/4]

void G4ParticleChange::AddSecondary ( G4Track * aSecondary)

Definition at line 100 of file G4ParticleChange.cc.

101{
102 // add a secondary
104}

Referenced by G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4HadronicProcess::FillResult(), G4BOptnCloning::GenerateBiasingFinalState(), G4AdjointForcedInteractionForGamma::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4NuVacOscProcess::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4XrayReflection::PostStepDoIt(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), and G4AdjointPhotoElectricModel::SampleSecondaries().

◆ CalcMomentum()

G4ThreeVector G4ParticleChange::CalcMomentum ( G4double energy,
G4ThreeVector direction,
G4double mass ) const
inline

Referenced by UpdateStepForAlongStep().

◆ DumpInfo()

void G4ParticleChange::DumpInfo ( ) const
overridevirtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 337 of file G4ParticleChange.cc.

338{
339 // use base-class DumpInfo
341
342 G4long oldprc = G4cout.precision(8);
343
344 G4cout << " Mass (GeV) : " << std::setw(20) << theMassChange / GeV
345 << G4endl;
346 G4cout << " Charge (eplus) : " << std::setw(20)
347 << theChargeChange / eplus << G4endl;
348 G4cout << " MagneticMoment : " << std::setw(20)
350 G4cout << " = : " << std::setw(20)
352 * 2. * theMassChange / c_squared / eplus / hbar_Planck
353 << "*[e hbar]/[2 m]" << G4endl;
354 G4cout << " Position - x (mm) : " << std::setw(20)
355 << thePositionChange.x() / mm << G4endl;
356 G4cout << " Position - y (mm) : " << std::setw(20)
357 << thePositionChange.y() / mm << G4endl;
358 G4cout << " Position - z (mm) : " << std::setw(20)
359 << thePositionChange.z() / mm << G4endl;
360 G4cout << " Time (ns) : " << std::setw(20)
361 << theTimeChange / ns << G4endl;
362 G4cout << " Proper Time (ns) : " << std::setw(20)
363 << theProperTimeChange / ns << G4endl;
364 G4cout << " Momentum Direct - x : " << std::setw(20)
366 G4cout << " Momentum Direct - y : " << std::setw(20)
368 G4cout << " Momentum Direct - z : " << std::setw(20)
370 G4cout << " Kinetic Energy (MeV): " << std::setw(20)
371 << theEnergyChange / MeV << G4endl;
372 G4cout << " Velocity (/c) : " << std::setw(20)
373 << theVelocityChange / c_light << G4endl;
374 G4cout << " Polarization - x : " << std::setw(20)
376 G4cout << " Polarization - y : " << std::setw(20)
378 G4cout << " Polarization - z : " << std::setw(20)
380 G4cout.precision(oldprc);
381}
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
G4ThreeVector theMomentumDirectionChange
G4ThreeVector thePolarizationChange
G4double theMagneticMomentChange
virtual void DumpInfo() const

Referenced by G4ParticleChangeForTransport::DumpInfo().

◆ GetCharge()

G4double G4ParticleChange::GetCharge ( ) const
inline

◆ GetEnergy()

◆ GetGlobalPosition()

G4ThreeVector G4ParticleChange::GetGlobalPosition ( const G4ThreeVector & displacement) const
inline

◆ GetGlobalTime()

G4double G4ParticleChange::GetGlobalTime ( G4double timeDelay = 0.0) const
inline

Referenced by AddSecondary(), and AddSecondary().

◆ GetLocalTime()

G4double G4ParticleChange::GetLocalTime ( G4double timeDelay = 0.0) const
inline

◆ GetMagneticMoment()

G4double G4ParticleChange::GetMagneticMoment ( ) const
inline

◆ GetMass()

G4double G4ParticleChange::GetMass ( ) const
inline

◆ GetMomentumDirection()

◆ GetPolarization()

const G4ThreeVector * G4ParticleChange::GetPolarization ( ) const
inline

◆ GetPosition()

const G4ThreeVector * G4ParticleChange::GetPosition ( ) const
inline

◆ GetProperTime()

G4double G4ParticleChange::GetProperTime ( ) const
inline

◆ GetVelocity()

G4double G4ParticleChange::GetVelocity ( ) const
inline

◆ Initialize()

void G4ParticleChange::Initialize ( const G4Track & track)
overridevirtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 107 of file G4ParticleChange.cc.

108{
109 // use base class's method at first
111
112 // set Energy/Momentum etc. equal to those of the parent particle
113 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
114 theEnergyChange = pParticle->GetKineticEnergy();
116 isVelocityChanged = false;
119 theProperTimeChange = pParticle->GetProperTime();
120
121 // set mass/charge/MagneticMoment of DynamicParticle
122 theMassChange = pParticle->GetMass();
123 theChargeChange = pParticle->GetCharge();
125
126 // set Position equal to those of the parent track
128
129 // set TimeChange equal to local time of the parent track
131
132 // set initial global time of the parent track
134}
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetProperTime() const
G4double GetMagneticMoment() const
const G4ThreeVector & GetPolarization() const
G4double GetVelocity() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void Initialize(const G4Track &)

Referenced by G4AdjointAlongStepWeightCorrection::AlongStepDoIt(), G4AdjointForcedInteractionForGamma::AlongStepDoIt(), G4ContinuousGainOfEnergy::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4BOptnForceFreeFlight::ApplyFinalStateBiasing(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4DNAElectronHoleRecombination::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4DNAMolecularDissociation::DecayIt(), G4BOptnCloning::GenerateBiasingFinalState(), G4AdjointForcedInteractionForGamma::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Channeling::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4ImportanceProcess::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4NeutronGeneralProcess::PostStepDoIt(), G4NuVacOscProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4PhononDownconversion::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4SpecialCuts::PostStepDoIt(), G4StepLimiter::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNLoss::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4VAdjointReverseReaction::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4WeightWindowProcess::PostStepDoIt(), and G4XrayReflection::PostStepDoIt().

◆ operator=()

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

◆ ProposeCharge()

void G4ParticleChange::ProposeCharge ( G4double finalCharge)
inline

◆ ProposeEnergy()

◆ ProposeGlobalTime()

void G4ParticleChange::ProposeGlobalTime ( G4double t)
inline

◆ ProposeLocalTime()

void G4ParticleChange::ProposeLocalTime ( G4double t)
inline

◆ ProposeMagneticMoment()

void G4ParticleChange::ProposeMagneticMoment ( G4double finalMagneticMoment)
inline

◆ ProposeMass()

void G4ParticleChange::ProposeMass ( G4double finalMass)
inline

◆ ProposeMomentumDirection() [1/2]

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

◆ ProposeMomentumDirection() [2/2]

◆ ProposePolarization() [1/2]

void G4ParticleChange::ProposePolarization ( const G4ThreeVector & finalPoralization)
inline

◆ ProposePolarization() [2/2]

◆ ProposePosition() [1/2]

void G4ParticleChange::ProposePosition ( const G4ThreeVector & finalPosition)
inline

◆ ProposePosition() [2/2]

◆ ProposeProperTime()

void G4ParticleChange::ProposeProperTime ( G4double finalProperTime)
inline

◆ ProposeVelocity()

◆ UpdateStepForAlongStep()

G4Step * G4ParticleChange::UpdateStepForAlongStep ( G4Step * Step)
overridevirtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 137 of file G4ParticleChange.cc.

138{
139 // A physics process always calculates the final state of the
140 // particle relative to the initial state at the beginning
141 // of the Step, i.e., based on information of G4Track (or
142 // equivalently the PreStepPoint).
143 // So, the differences (delta) between these two states have to be
144 // calculated and be accumulated in PostStepPoint
145
146 // Take note that the return type of GetMomentumDirectionChange() is a
147 // pointer to G4ParticleMometum. Also it is a normalized
148 // momentum vector
149
150 const G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
151 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
152
153 // set Mass/Charge/MagneticMoment
154 pPostStepPoint->SetMass(theMassChange);
155 pPostStepPoint->SetCharge(theChargeChange);
157
158 // calculate new kinetic energy
159 G4double preEnergy = pPreStepPoint->GetKineticEnergy();
161 pPostStepPoint->GetKineticEnergy() + (theEnergyChange - preEnergy);
162
163 // update kinetic energy and momentum direction
164 if(energy > 0.0)
165 {
166 // calculate new momentum
167 G4ThreeVector pMomentum = pPostStepPoint->GetMomentum()
169 theMassChange) - pPreStepPoint->GetMomentum());
170 G4double tMomentum2 = pMomentum.mag2();
171 G4ThreeVector direction(1.0, 0.0, 0.0);
172 if(tMomentum2 > 0.)
173 {
174 direction = pMomentum / std::sqrt(tMomentum2);
175 }
176 pPostStepPoint->SetMomentumDirection(direction);
177 pPostStepPoint->SetKineticEnergy(energy);
178
179 // if velocity is not set it should be recomputed
180 //
182 {
183 if (theMassChange > 0.0)
184 {
185 theVelocityChange = CLHEP::c_light *
186 std::sqrt(energy*(energy + 2*theMassChange))/(energy + theMassChange);
187 }
188 else
189 {
190 // zero mass particle
191 theVelocityChange = CLHEP::c_light;
192 // optical photon case
194 {
195 G4Track* pTrack = pStep->GetTrack();
196 G4double e = pTrack->GetKineticEnergy();
197 pTrack->SetKineticEnergy(energy);
199 pTrack->SetKineticEnergy(e);
200 }
201 }
202 }
203 pPostStepPoint->SetVelocity(theVelocityChange);
204 }
205 else
206 {
207 // stop case
208 pPostStepPoint->SetKineticEnergy(0.0);
209 pPostStepPoint->SetVelocity(0.0);
210 }
211
212 // update polarization
213 pPostStepPoint->AddPolarization(thePolarizationChange -
214 pPreStepPoint->GetPolarization());
215
216 // update position and time
217 pPostStepPoint->AddPosition(thePositionChange - pPreStepPoint->GetPosition());
218 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
219 pPostStepPoint->AddLocalTime(theTimeChange - theLocalTime0);
220 pPostStepPoint->AddProperTime(theProperTimeChange -
221 pPreStepPoint->GetProperTime());
222
224 {
225 pPostStepPoint->SetWeight(theParentWeight);
226 }
227
228#ifdef G4VERBOSE
230#endif
231
232 // update the G4Step specific attributes
233 return UpdateStepInfo(pStep);
234}
double G4double
Definition G4Types.hh:83
double mag2() const
G4ThreeVector CalcMomentum(G4double energy, G4ThreeVector direction, G4double mass) const
G4Step * UpdateStepInfo(G4Step *Step)
void AddPolarization(const G4ThreeVector &aValue)
void SetMagneticMoment(G4double value)
void SetKineticEnergy(const G4double aValue)
void SetWeight(G4double aValue)
void SetMass(G4double value)
void SetCharge(G4double value)
G4double GetProperTime() const
void AddProperTime(const G4double aValue)
void SetVelocity(G4double v)
void AddPosition(const G4ThreeVector &aValue)
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetPosition() const
void AddGlobalTime(const G4double aValue)
const G4ThreeVector & GetPolarization() const
G4double GetKineticEnergy() const
void AddLocalTime(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
G4double CalculateVelocityForOpticalPhoton() const
Definition G4Track.cc:160
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
void SetKineticEnergy(const G4double aValue)
virtual G4bool CheckIt(const G4Track &)
G4double energy(const ThreeVector &p, const G4double m)

◆ UpdateStepForAtRest()

G4Step * G4ParticleChange::UpdateStepForAtRest ( G4Step * Step)
overridevirtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 296 of file G4ParticleChange.cc.

297{
298 // A physics process always calculates the final state of the particle
299
300 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
301
302 // set Mass/Charge
303 pPostStepPoint->SetMass(theMassChange);
304 pPostStepPoint->SetCharge(theChargeChange);
306
307 // update kinetic energy and momentum direction
309 pPostStepPoint->SetKineticEnergy(theEnergyChange);
312 pPostStepPoint->SetVelocity(theVelocityChange);
313
314 // update polarization
315 pPostStepPoint->SetPolarization(thePolarizationChange);
316
317 // update position and time
318 pPostStepPoint->SetPosition(thePositionChange);
319 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
320 pPostStepPoint->SetLocalTime(theTimeChange);
321 pPostStepPoint->SetProperTime(theProperTimeChange);
322
324 {
325 pPostStepPoint->SetWeight(theParentWeight);
326 }
327
328#ifdef G4VERBOSE
330#endif
331
332 // update the G4Step specific attributes
333 return UpdateStepInfo(pStep);
334}
void SetLocalTime(const G4double aValue)
void SetProperTime(const G4double aValue)
void SetPosition(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4double CalculateVelocity() const

◆ UpdateStepForPostStep()

G4Step * G4ParticleChange::UpdateStepForPostStep ( G4Step * Step)
overridevirtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 237 of file G4ParticleChange.cc.

238{
239 // A physics process always calculates the final state of the particle
240
241 // Take note that the return type of GetMomentumChange() is a
242 // pointer to G4ParticleMometum. Also it is a normalized
243 // momentum vector
244
245 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
246 G4Track* pTrack = pStep->GetTrack();
247
248 // set Mass/Charge
249 pPostStepPoint->SetMass(theMassChange);
250 pPostStepPoint->SetCharge(theChargeChange);
252
253 // update kinetic energy and momentum direction
255
256 // calculate velocity
257 if(theEnergyChange > 0.0)
258 {
259 pPostStepPoint->SetKineticEnergy(theEnergyChange);
262 {
264 }
265 pPostStepPoint->SetVelocity(theVelocityChange);
266 }
267 else
268 {
269 pPostStepPoint->SetKineticEnergy(0.0);
270 pPostStepPoint->SetVelocity(0.0);
271 }
272
273 // update polarization
274 pPostStepPoint->SetPolarization(thePolarizationChange);
275
276 // update position and time
277 pPostStepPoint->SetPosition(thePositionChange);
278 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
279 pPostStepPoint->SetLocalTime(theTimeChange);
280 pPostStepPoint->SetProperTime(theProperTimeChange);
281
283 {
284 pPostStepPoint->SetWeight(theParentWeight);
285 }
286
287#ifdef G4VERBOSE
289#endif
290
291 // update the G4Step specific attributes
292 return UpdateStepInfo(pStep);
293}

◆ UpdateStepInfo()

Member Data Documentation

◆ isVelocityChanged

G4bool G4ParticleChange::isVelocityChanged = false
protected

◆ theChargeChange

G4double G4ParticleChange::theChargeChange = 0.0
protected

◆ theEnergyChange

◆ theGlobalTime0

G4double G4ParticleChange::theGlobalTime0 = 0.0
protected

Definition at line 205 of file G4ParticleChange.hh.

Referenced by Initialize().

◆ theLocalTime0

G4double G4ParticleChange::theLocalTime0 = 0.0
protected

◆ theMagneticMomentChange

G4double G4ParticleChange::theMagneticMomentChange = 0.0
protected

◆ theMassChange

G4double G4ParticleChange::theMassChange = 0.0
protected

◆ theMomentumDirectionChange

◆ thePolarizationChange

◆ thePositionChange

◆ theProperTimeChange

G4double G4ParticleChange::theProperTimeChange = 0.0
protected

◆ theTimeChange

◆ theVelocityChange


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