Geant4 10.7.0
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 ()
 
virtual ~G4ParticleChange ()
 
G4bool operator== (const G4ParticleChange &right) const
 
G4bool operator!= (const G4ParticleChange &right) const
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
virtual G4StepUpdateStepForAtRest (G4Step *Step)
 
virtual G4StepUpdateStepForPostStep (G4Step *Step)
 
virtual void Initialize (const G4Track &)
 
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)
 
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

 G4ParticleChange (const G4ParticleChange &right)
 
G4ParticleChangeoperator= (const G4ParticleChange &right)
 
G4StepUpdateStepInfo (G4Step *Step)
 
- 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
 

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
 
const G4TracktheCurrentTrack = nullptr
 
- 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
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 

Detailed Description

Definition at line 56 of file G4ParticleChange.hh.

Constructor & Destructor Documentation

◆ G4ParticleChange() [1/2]

G4ParticleChange::G4ParticleChange ( )

Definition at line 41 of file G4ParticleChange.cc.

◆ ~G4ParticleChange()

G4ParticleChange::~G4ParticleChange ( )
virtual

Definition at line 47 of file G4ParticleChange.cc.

48{
49}

◆ G4ParticleChange() [2/2]

G4ParticleChange::G4ParticleChange ( const G4ParticleChange right)
protected

Definition at line 52 of file G4ParticleChange.cc.

53 : G4VParticleChange(right)
54{
56
66 isVelocityChanged = true;
70}
G4ThreeVector thePositionChange
G4ThreeVector theMomentumDirectionChange
G4double theProperTimeChange
G4ThreeVector thePolarizationChange
G4double theMagneticMomentChange
const G4Track * theCurrentTrack

Member Function Documentation

◆ AddSecondary() [1/4]

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

Definition at line 132 of file G4ParticleChange.cc.

134{
135 // create track
136 G4Track* aTrack = new G4Track(aParticle, GetGlobalTime(), thePositionChange);
137
138 // set IsGoodGorTrackingFlag
139 if(IsGoodForTracking)
140 aTrack->SetGoodForTrackingFlag();
141
142 // touchable handle is copied to keep the pointer
144
145 // add a secondary
147}
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)

◆ AddSecondary() [2/4]

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

Definition at line 169 of file G4ParticleChange.cc.

171{
172 // create track
173 G4Track* aTrack = new G4Track(aParticle, newTime, thePositionChange);
174
175 // set IsGoodGorTrackingFlag
176 if(IsGoodForTracking)
177 aTrack->SetGoodForTrackingFlag();
178
179 // touchable handle is copied to keep the pointer
181
182 // add a secondary
184}

◆ AddSecondary() [3/4]

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

Definition at line 150 of file G4ParticleChange.cc.

153{
154 // create track
155 G4Track* aTrack = new G4Track(aParticle, GetGlobalTime(), newPosition);
156
157 // set IsGoodGorTrackingFlag
158 if(IsGoodForTracking)
159 aTrack->SetGoodForTrackingFlag();
160
161 // touchable is a temporary object, so you cannot keep the pointer
162 aTrack->SetTouchableHandle(nullptr);
163
164 // add a secondary
166}

◆ AddSecondary() [4/4]

void G4ParticleChange::AddSecondary ( G4Track aSecondary)

◆ CalcMomentum()

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

◆ CheckIt()

G4bool G4ParticleChange::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 467 of file G4ParticleChange.cc.

468{
469 G4bool exitWithError = false;
470 G4double accuracy;
471 static G4ThreadLocal G4int nError = 0;
472#ifdef G4VERBOSE
473 const G4int maxError = 30;
474#endif
475
476 // no check in case of "fStopAndKill"
478 return G4VParticleChange::CheckIt(aTrack);
479
480 // MomentumDirection should be unit vector
481 G4bool itsOKforMomentum = true;
482 if(theEnergyChange > 0.)
483 {
484 accuracy = std::fabs(theMomentumDirectionChange.mag2() - 1.0);
485 if(accuracy > accuracyForWarning)
486 {
487 itsOKforMomentum = false;
488 nError += 1;
489 exitWithError = exitWithError || (accuracy > accuracyForException);
490#ifdef G4VERBOSE
491 if(nError < maxError)
492 {
493 G4cout << " G4ParticleChange::CheckIt : ";
494 G4cout << "the Momentum Change is not unit vector !!"
495 << " Difference: " << accuracy << G4endl;
497 << " E=" << aTrack.GetKineticEnergy() / MeV
498 << " pos=" << aTrack.GetPosition().x() / m << ", "
499 << aTrack.GetPosition().y() / m << ", "
500 << aTrack.GetPosition().z() / m << G4endl;
501 }
502#endif
503 }
504 }
505
506 // both global and proper time should not go back
507 G4bool itsOKforGlobalTime = true;
508 accuracy = (aTrack.GetLocalTime() - theTimeChange) / ns;
509 if(accuracy > accuracyForWarning)
510 {
511 itsOKforGlobalTime = false;
512 nError += 1;
513 exitWithError = exitWithError || (accuracy > accuracyForException);
514#ifdef G4VERBOSE
515 if(nError < maxError)
516 {
517 G4cout << " G4ParticleChange::CheckIt : ";
518 G4cout << "the local time goes back !!"
519 << " Difference: " << accuracy << "[ns] " << G4endl;
521 << " E=" << aTrack.GetKineticEnergy() / MeV
522 << " pos=" << aTrack.GetPosition().x() / m << ", "
523 << aTrack.GetPosition().y() / m << ", "
524 << aTrack.GetPosition().z() / m
525 << " global time=" << aTrack.GetGlobalTime() / ns
526 << " local time=" << aTrack.GetLocalTime() / ns
527 << " proper time=" << aTrack.GetProperTime() / ns << G4endl;
528 }
529#endif
530 }
531
532 G4bool itsOKforProperTime = true;
533 accuracy = (aTrack.GetProperTime() - theProperTimeChange) / ns;
534 if(accuracy > accuracyForWarning)
535 {
536 itsOKforProperTime = false;
537 nError += 1;
538 exitWithError = exitWithError || (accuracy > accuracyForException);
539#ifdef G4VERBOSE
540 if(nError < maxError)
541 {
542 G4cout << " G4ParticleChange::CheckIt : ";
543 G4cout << "the proper time goes back !!"
544 << " Difference: " << accuracy << "[ns] " << G4endl;
546 << " E=" << aTrack.GetKineticEnergy() / MeV
547 << " pos=" << aTrack.GetPosition().x() / m << ", "
548 << aTrack.GetPosition().y() / m << ", "
549 << aTrack.GetPosition().z() / m
550 << " global time=" << aTrack.GetGlobalTime() / ns
551 << " local time=" << aTrack.GetLocalTime() / ns
552 << " proper time=" << aTrack.GetProperTime() / ns << G4endl;
553 }
554#endif
555 }
556
557 // kinetic Energy should not be negative
558 G4bool itsOKforEnergy = true;
559 accuracy = -1.0 * theEnergyChange / MeV;
560 if(accuracy > accuracyForWarning)
561 {
562 itsOKforEnergy = false;
563 nError += 1;
564 exitWithError = exitWithError || (accuracy > accuracyForException);
565#ifdef G4VERBOSE
566 if(nError < maxError)
567 {
568 G4cout << " G4ParticleChange::CheckIt : ";
569 G4cout << "the kinetic energy is negative !!"
570 << " Difference: " << accuracy << "[MeV] " << G4endl;
572 << " E=" << aTrack.GetKineticEnergy() / MeV
573 << " pos=" << aTrack.GetPosition().x() / m << ", "
574 << aTrack.GetPosition().y() / m << ", "
575 << aTrack.GetPosition().z() / m << G4endl;
576 }
577#endif
578 }
579
580 // velocity should not be less than c_light
581 G4bool itsOKforVelocity = true;
582 if(theVelocityChange < 0.)
583 {
584 itsOKforVelocity = false;
585 nError += 1;
586 exitWithError = true;
587#ifdef G4VERBOSE
588 if(nError < maxError)
589 {
590 G4cout << " G4ParticleChange::CheckIt : ";
591 G4cout << "the velocity is negative !!"
592 << " Velocity: " << theVelocityChange / c_light << G4endl;
594 << " E=" << aTrack.GetKineticEnergy() / MeV
595 << " pos=" << aTrack.GetPosition().x() / m << ", "
596 << aTrack.GetPosition().y() / m << ", "
597 << aTrack.GetPosition().z() / m << G4endl;
598 }
599#endif
600 }
601
602 accuracy = theVelocityChange / c_light - 1.0;
603 if(accuracy > accuracyForWarning)
604 {
605 itsOKforVelocity = false;
606 nError += 1;
607 exitWithError = exitWithError || (accuracy > accuracyForException);
608#ifdef G4VERBOSE
609 if(nError < maxError)
610 {
611 G4cout << " G4ParticleChange::CheckIt : ";
612 G4cout << "the velocity is greater than c_light !!" << G4endl;
613 G4cout << " Velocity: " << theVelocityChange / c_light << G4endl;
615 << " E=" << aTrack.GetKineticEnergy() / MeV
616 << " pos=" << aTrack.GetPosition().x() / m << ", "
617 << aTrack.GetPosition().y() / m << ", "
618 << aTrack.GetPosition().z() / m << G4endl;
619 }
620#endif
621 }
622
623 G4bool itsOK = itsOKforMomentum && itsOKforEnergy && itsOKforVelocity &&
624 itsOKforProperTime && itsOKforGlobalTime;
625 // dump out information of this particle change
626#ifdef G4VERBOSE
627 if(!itsOK)
628 {
629 DumpInfo();
630 }
631#endif
632
633 // Exit with error
634 if(exitWithError)
635 {
636 G4Exception("G4ParticleChange::CheckIt()", "TRACK003", EventMustBeAborted,
637 "momentum, energy, and/or time was illegal");
638 }
639 // correction
640 if(!itsOKforMomentum)
641 {
644 }
645 if(!itsOKforGlobalTime)
646 {
647 theTimeChange = aTrack.GetLocalTime();
648 }
649 if(!itsOKforProperTime)
650 {
652 }
653 if(!itsOKforEnergy)
654 {
655 theEnergyChange = 0.0;
656 }
657 if(!itsOKforVelocity)
658 {
659 theVelocityChange = c_light;
660 }
661
662 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
663 return itsOK;
664}
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double mag2() const
double y() const
double mag() const
virtual void DumpInfo() const
const G4String & GetParticleName() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4double GetLocalTime() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
virtual G4bool CheckIt(const G4Track &)
static const G4double accuracyForException
static const G4double accuracyForWarning
G4TrackStatus GetTrackStatus() const
#define G4ThreadLocal
Definition: tls.hh:77
#define ns
Definition: xmlparse.cc:614

Referenced by UpdateStepForAlongStep(), G4ParticleChangeForTransport::UpdateStepForAlongStep(), UpdateStepForAtRest(), and UpdateStepForPostStep().

◆ DumpInfo()

void G4ParticleChange::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 420 of file G4ParticleChange.cc.

421{
422 // use base-class DumpInfo
424
425 G4int oldprc = G4cout.precision(3);
426
427 G4cout << " Mass (GeV) : " << std::setw(20) << theMassChange / GeV
428 << G4endl;
429 G4cout << " Charge (eplus) : " << std::setw(20)
430 << theChargeChange / eplus << G4endl;
431 G4cout << " MagneticMoment : " << std::setw(20)
433 G4cout << " : = " << std::setw(20)
435 * 2. * theMassChange / c_squared / eplus / hbar_Planck
436 << "*[e hbar]/[2 m]" << G4endl;
437 G4cout << " Position - x (mm) : " << std::setw(20)
438 << thePositionChange.x() / mm << G4endl;
439 G4cout << " Position - y (mm) : " << std::setw(20)
440 << thePositionChange.y() / mm << G4endl;
441 G4cout << " Position - z (mm) : " << std::setw(20)
442 << thePositionChange.z() / mm << G4endl;
443 G4cout << " Time (ns) : " << std::setw(20)
444 << theTimeChange / ns << G4endl;
445 G4cout << " Proper Time (ns) : " << std::setw(20)
447 G4cout << " Momentum Direct - x : " << std::setw(20)
449 G4cout << " Momentum Direct - y : " << std::setw(20)
451 G4cout << " Momentum Direct - z : " << std::setw(20)
453 G4cout << " Kinetic Energy (MeV): " << std::setw(20)
454 << theEnergyChange / MeV << G4endl;
455 G4cout << " Velocity (/c): " << std::setw(20)
456 << theVelocityChange / c_light << G4endl;
457 G4cout << " Polarization - x : " << std::setw(20)
459 G4cout << " Polarization - y : " << std::setw(20)
461 G4cout << " Polarization - z : " << std::setw(20)
463 G4cout.precision(oldprc);
464}
virtual void DumpInfo() const

Referenced by CheckIt(), and G4ParticleChangeForTransport::DumpInfo().

◆ GetCharge()

G4double G4ParticleChange::GetCharge ( ) const

◆ GetEnergy()

◆ GetGlobalPosition()

G4ThreeVector G4ParticleChange::GetGlobalPosition ( const G4ThreeVector displacement) const

◆ GetGlobalTime()

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

Referenced by AddSecondary().

◆ GetLocalTime()

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

◆ GetMagneticMoment()

G4double G4ParticleChange::GetMagneticMoment ( ) const

◆ GetMass()

G4double G4ParticleChange::GetMass ( ) const

◆ GetMomentumDirection()

◆ GetPolarization()

const G4ThreeVector * G4ParticleChange::GetPolarization ( ) const

◆ GetPosition()

const G4ThreeVector * G4ParticleChange::GetPosition ( ) const

◆ GetProperTime()

G4double G4ParticleChange::GetProperTime ( ) const

◆ GetVelocity()

G4double G4ParticleChange::GetVelocity ( ) const

◆ Initialize()

void G4ParticleChange::Initialize ( const G4Track track)
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 194 of file G4ParticleChange.cc.

195{
196 // use base class's method at first
198 theCurrentTrack = &track;
199
200 // set Energy/Momentum etc. equal to those of the parent particle
201 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
202 theEnergyChange = pParticle->GetKineticEnergy();
204 isVelocityChanged = false;
207 theProperTimeChange = pParticle->GetProperTime();
208
209 // set mass/charge/MagneticMoment of DynamicParticle
210 theMassChange = pParticle->GetMass();
211 theChargeChange = pParticle->GetCharge();
213
214 // set Position equal to those of the parent track
216
217 // set TimeChange equal to local time of the parent track
218 theTimeChange = track.GetLocalTime();
219
220 // set initial Local/Global time of the parent track
221 theLocalTime0 = track.GetLocalTime();
223}
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 G4DynamicParticle * GetDynamicParticle() const
virtual void Initialize(const G4Track &)

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

◆ operator!=()

G4bool G4ParticleChange::operator!= ( const G4ParticleChange right) const

Definition at line 126 of file G4ParticleChange.cc.

127{
128 return ((G4VParticleChange*) this != (G4VParticleChange*) &right);
129}

◆ operator=()

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

Definition at line 73 of file G4ParticleChange.cc.

74{
75 if(this != &right)
76 {
78 {
79 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
80 {
81 if((*theListOfSecondaries)[index])
82 delete(*theListOfSecondaries)[index];
83 }
84 }
86
89 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
90 {
91 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index]));
92 theListOfSecondaries->SetElement(index, newTrack);
93 }
94
97
107 isVelocityChanged = true;
111
115 }
116 return *this;
117}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag

◆ operator==()

G4bool G4ParticleChange::operator== ( const G4ParticleChange right) const

Definition at line 120 of file G4ParticleChange.cc.

121{
122 return ((G4VParticleChange*) this == (G4VParticleChange*) &right);
123}

◆ ProposeCharge()

void G4ParticleChange::ProposeCharge ( G4double  finalCharge)

◆ ProposeEnergy()

◆ ProposeGlobalTime()

◆ ProposeLocalTime()

◆ ProposeMagneticMoment()

void G4ParticleChange::ProposeMagneticMoment ( G4double  finalMagneticMoment)

◆ ProposeMass()

void G4ParticleChange::ProposeMass ( G4double  finalMass)

◆ ProposeMomentumDirection() [1/2]

void G4ParticleChange::ProposeMomentumDirection ( const G4ThreeVector Pfinal)

◆ ProposeMomentumDirection() [2/2]

◆ ProposePolarization() [1/2]

void G4ParticleChange::ProposePolarization ( const G4ThreeVector finalPoralization)

◆ ProposePolarization() [2/2]

◆ ProposePosition() [1/2]

void G4ParticleChange::ProposePosition ( const G4ThreeVector finalPosition)

◆ ProposePosition() [2/2]

◆ ProposeProperTime()

void G4ParticleChange::ProposeProperTime ( G4double  finalProperTime)

◆ ProposeVelocity()

◆ UpdateStepForAlongStep()

G4Step * G4ParticleChange::UpdateStepForAlongStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 226 of file G4ParticleChange.cc.

227{
228 // A physics process always calculates the final state of the
229 // particle relative to the initial state at the beginning
230 // of the Step, i.e., based on information of G4Track (or
231 // equivalently the PreStepPoint).
232 // So, the differences (delta) between these two states have to be
233 // calculated and be accumulated in PostStepPoint
234
235 // Take note that the return type of GetMomentumDirectionChange() is a
236 // pointer to G4ParticleMometum. Also it is a normalized
237 // momentum vector
238
239 G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
240 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
241 G4Track* pTrack = pStep->GetTrack();
242 G4double mass = theMassChange;
243
244 // set Mass/Charge/MagneticMoment
245 pPostStepPoint->SetMass(theMassChange);
246 pPostStepPoint->SetCharge(theChargeChange);
248
249 // calculate new kinetic energy
250 G4double preEnergy = pPreStepPoint->GetKineticEnergy();
252 pPostStepPoint->GetKineticEnergy() + (theEnergyChange - preEnergy);
253
254 // update kinetic energy and momentum direction
255 if(energy > 0.0)
256 {
257 // calculate new momentum
258 G4ThreeVector pMomentum = pPostStepPoint->GetMomentum()
260 - pPreStepPoint->GetMomentum());
261 G4double tMomentum = pMomentum.mag();
262 G4ThreeVector direction(1.0, 0.0, 0.0);
263 if(tMomentum > 0.)
264 {
265 G4double inv_Momentum = 1.0 / tMomentum;
266 direction = pMomentum * inv_Momentum;
267 }
268 pPostStepPoint->SetMomentumDirection(direction);
269 pPostStepPoint->SetKineticEnergy(energy);
270 }
271 else
272 {
273 // stop case
274 // pPostStepPoint->SetMomentumDirection(G4ThreeVector(1., 0., 0.));
275 pPostStepPoint->SetKineticEnergy(0.0);
276 }
277 // calculate velocity
279 {
280 if(energy > 0.0)
281 {
282 pTrack->SetKineticEnergy(energy);
284 pTrack->SetKineticEnergy(preEnergy);
285 }
286 else if(theMassChange > 0.0)
287 {
288 theVelocityChange = 0.0;
289 }
290 }
291 pPostStepPoint->SetVelocity(theVelocityChange);
292
293 // update polarization
294 pPostStepPoint->AddPolarization(thePolarizationChange -
295 pPreStepPoint->GetPolarization());
296
297 // update position and time
298 pPostStepPoint->AddPosition(thePositionChange - pPreStepPoint->GetPosition());
299 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
300 pPostStepPoint->AddLocalTime(theTimeChange - theLocalTime0);
301 pPostStepPoint->AddProperTime(theProperTimeChange -
302 pPreStepPoint->GetProperTime());
303
305 {
306 pPostStepPoint->SetWeight(theParentWeight);
307 }
308
309#ifdef G4VERBOSE
310 G4Track* aTrack = pStep->GetTrack();
311 if(debugFlag) { CheckIt(*aTrack); }
312#endif
313
314 // update the G4Step specific attributes
315 return UpdateStepInfo(pStep);
316}
G4ThreeVector CalcMomentum(G4double energy, G4ThreeVector direction, G4double mass) const
G4Step * UpdateStepInfo(G4Step *Step)
virtual G4bool CheckIt(const G4Track &)
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 CalculateVelocity() const
void SetKineticEnergy(const G4double aValue)
G4double energy(const ThreeVector &p, const G4double m)

◆ UpdateStepForAtRest()

G4Step * G4ParticleChange::UpdateStepForAtRest ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 378 of file G4ParticleChange.cc.

379{
380 // A physics process always calculates the final state of the particle
381
382 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
383
384 // set Mass/Charge
385 pPostStepPoint->SetMass(theMassChange);
386 pPostStepPoint->SetCharge(theChargeChange);
388
389 // update kinetic energy and momentum direction
391 pPostStepPoint->SetKineticEnergy(theEnergyChange);
393 theVelocityChange = pStep->GetTrack()->CalculateVelocity();
394 pPostStepPoint->SetVelocity(theVelocityChange);
395
396 // update polarization
397 pPostStepPoint->SetPolarization(thePolarizationChange);
398
399 // update position and time
400 pPostStepPoint->SetPosition(thePositionChange);
401 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
402 pPostStepPoint->SetLocalTime(theTimeChange);
403 pPostStepPoint->SetProperTime(theProperTimeChange);
404
406 {
407 pPostStepPoint->SetWeight(theParentWeight);
408 }
409
410#ifdef G4VERBOSE
411 G4Track* aTrack = pStep->GetTrack();
412 if(debugFlag) { CheckIt(*aTrack); }
413#endif
414
415 // update the G4Step specific attributes
416 return UpdateStepInfo(pStep);
417}
void SetLocalTime(const G4double aValue)
void SetProperTime(const G4double aValue)
void SetPosition(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)

◆ UpdateStepForPostStep()

G4Step * G4ParticleChange::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 319 of file G4ParticleChange.cc.

320{
321 // A physics process always calculates the final state of the particle
322
323 // Take note that the return type of GetMomentumChange() is a
324 // pointer to G4ParticleMometum. Also it is a normalized
325 // momentum vector
326
327 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
328 G4Track* pTrack = pStep->GetTrack();
329
330 // set Mass/Charge
331 pPostStepPoint->SetMass(theMassChange);
332 pPostStepPoint->SetCharge(theChargeChange);
334
335 // update kinetic energy and momentum direction
337 pPostStepPoint->SetKineticEnergy(theEnergyChange);
338
339 // calculate velocity
342 {
343 if(theEnergyChange > 0.0)
344 {
346 }
347 else if(theMassChange > 0.0)
348 {
349 theVelocityChange = 0.0;
350 }
351 }
352 pPostStepPoint->SetVelocity(theVelocityChange);
353
354 // update polarization
355 pPostStepPoint->SetPolarization(thePolarizationChange);
356
357 // update position and time
358 pPostStepPoint->SetPosition(thePositionChange);
359 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
360 pPostStepPoint->SetLocalTime(theTimeChange);
361 pPostStepPoint->SetProperTime(theProperTimeChange);
362
364 {
365 pPostStepPoint->SetWeight(theParentWeight);
366 }
367
368#ifdef G4VERBOSE
369 G4Track* aTrack = pStep->GetTrack();
370 if(debugFlag) { CheckIt(*aTrack); }
371#endif
372
373 // update the G4Step specific attributes
374 return UpdateStepInfo(pStep);
375}

◆ UpdateStepInfo()

Member Data Documentation

◆ isVelocityChanged

◆ theChargeChange

G4double G4ParticleChange::theChargeChange = 0.0
protected

◆ theCurrentTrack

const G4Track* G4ParticleChange::theCurrentTrack = nullptr
protected

Definition at line 232 of file G4ParticleChange.hh.

Referenced by AddSecondary(), G4ParticleChange(), Initialize(), and operator=().

◆ theEnergyChange

◆ theGlobalTime0

G4double G4ParticleChange::theGlobalTime0 = 0.0
protected

Definition at line 212 of file G4ParticleChange.hh.

Referenced by G4ParticleChange(), Initialize(), and operator=().

◆ 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

◆ theTimeChange

◆ theVelocityChange


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