Geant4 9.6.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
 
G4double theVelocityChange
 
G4bool isVelocityChanged
 
G4ThreeVector thePositionChange
 
G4double theGlobalTime0
 
G4double theLocalTime0
 
G4double theTimeChange
 
G4double theProperTimeChange
 
G4double theMassChange
 
G4double theChargeChange
 
G4double theMagneticMomentChange
 
const G4TracktheCurrentTrack
 
- Protected Attributes inherited from G4VParticleChange
G4TrackFastVectortheListOfSecondaries
 
G4int theNumberOfSecondaries
 
G4int theSizeOftheListOfSecondaries
 
G4TrackStatus theStatusChange
 
G4SteppingControl theSteppingControlFlag
 
G4double theLocalEnergyDeposit
 
G4double theNonIonizingEnergyDeposit
 
G4double theTrueStepLength
 
G4bool theFirstStepInVolume
 
G4bool theLastStepInVolume
 
G4double theParentWeight
 
G4bool isParentWeightProposed
 
G4bool fSetSecondaryWeightByProcess
 
G4double theParentGlobalTime
 
G4int verboseLevel
 
G4bool debugFlag
 

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

Constructor & Destructor Documentation

◆ G4ParticleChange() [1/2]

G4ParticleChange::G4ParticleChange ( )

Definition at line 54 of file G4ParticleChange.cc.

58 theEnergyChange(0.),
65{
66}
G4ThreeVector thePositionChange
G4ThreeVector theMomentumDirectionChange
G4double theProperTimeChange
G4ThreeVector thePolarizationChange
G4double theMagneticMomentChange
const G4Track * theCurrentTrack

◆ ~G4ParticleChange()

G4ParticleChange::~G4ParticleChange ( )
virtual

Definition at line 68 of file G4ParticleChange.cc.

69{
70#ifdef G4VERBOSE
71 if (verboseLevel>2) {
72 G4cout << "G4ParticleChange::~G4ParticleChange() " << G4endl;
73 }
74#endif
75}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ G4ParticleChange() [2/2]

G4ParticleChange::G4ParticleChange ( const G4ParticleChange right)
protected

Member Function Documentation

◆ AddSecondary() [1/4]

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

Definition at line 168 of file G4ParticleChange.cc.

170{
171 // create track
172 G4Track* aTrack = new G4Track(aParticle, GetGlobalTime(), thePositionChange);
173
174 // set IsGoodGorTrackingFlag
175 if (IsGoodForTracking) aTrack->SetGoodForTrackingFlag();
176
177 // Touchable handle is copied to keep the pointer
179
180 // add a secondary
182}
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 201 of file G4ParticleChange.cc.

204{
205 // create track
206 G4Track* aTrack = new G4Track(aParticle, newTime, thePositionChange);
207
208 // set IsGoodGorTrackingFlag
209 if (IsGoodForTracking) aTrack->SetGoodForTrackingFlag();
210
211 // Touchable handle is copied to keep the pointer
213
214 // add a secondary
216}

◆ AddSecondary() [3/4]

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

Definition at line 184 of file G4ParticleChange.cc.

187{
188 // create track
189 G4Track* aTrack = new G4Track(aParticle, GetGlobalTime(), newPosition);
190
191 // set IsGoodGorTrackingFlag
192 if (IsGoodForTracking) aTrack->SetGoodForTrackingFlag();
193
194 // Touchable is a temporary object, so you cannot keep the pointer
195 aTrack->SetTouchableHandle((G4VTouchable*)0);
196
197 // add a secondary
199}

◆ AddSecondary() [4/4]

void G4ParticleChange::AddSecondary ( G4Track aSecondary)

Definition at line 218 of file G4ParticleChange.cc.

219{
220 // add a secondary
222}

Referenced by G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4MuonMinusCaptureAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4HadronicProcess::FillResult(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4QSynchRad::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

◆ CheckIt()

G4bool G4ParticleChange::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 486 of file G4ParticleChange.cc.

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

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

◆ DumpInfo()

void G4ParticleChange::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 425 of file G4ParticleChange.cc.

426{
427// use base-class DumpInfo
429
430 G4int oldprc = G4cout.precision(3);
431
432 G4cout << " Mass (GeV) : "
433 << std::setw(20) << theMassChange/GeV
434 << G4endl;
435 G4cout << " Charge (eplus) : "
436 << std::setw(20) << theChargeChange/eplus
437 << G4endl;
438 G4cout << " MagneticMoment : "
439 << std::setw(20) << theMagneticMomentChange << G4endl;
440 G4cout << " : = " << std::setw(20)
441 << theMagneticMomentChange*2.*theMassChange/c_squared/eplus/hbar_Planck
442 << "*[e hbar]/[2 m]"
443 << G4endl;
444 G4cout << " Position - x (mm) : "
445 << std::setw(20) << thePositionChange.x()/mm
446 << G4endl;
447 G4cout << " Position - y (mm) : "
448 << std::setw(20) << thePositionChange.y()/mm
449 << G4endl;
450 G4cout << " Position - z (mm) : "
451 << std::setw(20) << thePositionChange.z()/mm
452 << G4endl;
453 G4cout << " Time (ns) : "
454 << std::setw(20) << theTimeChange/ns
455 << G4endl;
456 G4cout << " Proper Time (ns) : "
457 << std::setw(20) << theProperTimeChange/ns
458 << G4endl;
459 G4cout << " Momentum Direct - x : "
460 << std::setw(20) << theMomentumDirectionChange.x()
461 << G4endl;
462 G4cout << " Momentum Direct - y : "
463 << std::setw(20) << theMomentumDirectionChange.y()
464 << G4endl;
465 G4cout << " Momentum Direct - z : "
466 << std::setw(20) << theMomentumDirectionChange.z()
467 << G4endl;
468 G4cout << " Kinetic Energy (MeV): "
469 << std::setw(20) << theEnergyChange/MeV
470 << G4endl;
471 G4cout << " Velocity (/c): "
472 << std::setw(20) << theVelocityChange/c_light
473 << G4endl;
474 G4cout << " Polarization - x : "
475 << std::setw(20) << thePolarizationChange.x()
476 << G4endl;
477 G4cout << " Polarization - y : "
478 << std::setw(20) << thePolarizationChange.y()
479 << G4endl;
480 G4cout << " Polarization - z : "
481 << std::setw(20) << thePolarizationChange.z()
482 << G4endl;
483 G4cout.precision(oldprc);
484}
virtual void DumpInfo() const

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

◆ GetCharge()

G4double G4ParticleChange::GetCharge ( ) const

◆ GetEnergy()

G4double G4ParticleChange::GetEnergy ( ) const

◆ GetGlobalPosition()

G4ThreeVector G4ParticleChange::GetGlobalPosition ( const G4ThreeVector displacement) const

◆ GetGlobalTime()

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

◆ 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 228 of file G4ParticleChange.cc.

229{
230 // use base class's method at first
232 theCurrentTrack= &track;
233
234 // set Energy/Momentum etc. equal to those of the parent particle
235 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
236 theEnergyChange = pParticle->GetKineticEnergy();
238 isVelocityChanged = false;
241 theProperTimeChange = pParticle->GetProperTime();
242
243 // Set mass/charge/MagneticMoment of DynamicParticle
244 theMassChange = pParticle->GetMass();
245 theChargeChange = pParticle->GetCharge();
247
248 // set Position equal to those of the parent track
250
251 // set TimeChange equal to local time of the parent track
252 theTimeChange = track.GetLocalTime();
253
254 // set initial Local/Global time of the parent track
255 theLocalTime0 = track.GetLocalTime();
257
258}
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(), G4hImpactIonisation::AlongStepDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4MuonMinusCaptureAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4DNAMolecularDecay::DecayIt(), SpecialCuts::PostStepDoIt(), G4ImportanceProcess::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4WeightWindowProcess::PostStepDoIt(), G4VAdjointReverseReaction::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4StepLimiter::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpRayleigh::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), and G4QSynchRad::PostStepDoIt().

◆ operator!=()

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

Definition at line 158 of file G4ParticleChange.cc.

159{
160 return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
161}

◆ operator=()

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

Definition at line 102 of file G4ParticleChange.cc.

103{
104#ifdef G4VERBOSE
105 if (verboseLevel>1) {
106 G4cout << "G4ParticleChange:: assignment operator is called " << G4endl;
107 }
108#endif
109 if (this != &right){
111#ifdef G4VERBOSE
112 if (verboseLevel>0) {
113 G4cout << "G4ParticleChange: assignment operator Warning ";
114 G4cout << "theListOfSecondaries is not empty ";
115 }
116#endif
117 for (G4int index= 0; index<theNumberOfSecondaries; index++){
118 if ( (*theListOfSecondaries)[index] ) delete (*theListOfSecondaries)[index] ;
119 }
120 }
121 delete theListOfSecondaries;
122
125 for (G4int index = 0; index<theNumberOfSecondaries; index++){
126 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
127 theListOfSecondaries->SetElement(index, newTrack); }
128
131
141 isVelocityChanged = true;
145
149 }
150 return *this;
151}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag

◆ operator==()

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

Definition at line 153 of file G4ParticleChange.cc.

154{
155 return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
156}

◆ 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()

void G4ParticleChange::ProposeVelocity ( G4double  finalVelocity)

◆ UpdateStepForAlongStep()

G4Step * G4ParticleChange::UpdateStepForAlongStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 264 of file G4ParticleChange.cc.

265{
266 // A physics process always calculates the final state of the
267 // particle relative to the initial state at the beginning
268 // of the Step, i.e., based on information of G4Track (or
269 // equivalently the PreStepPoint).
270 // So, the differences (delta) between these two states have to be
271 // calculated and be accumulated in PostStepPoint.
272
273 // Take note that the return type of GetMomentumDirectionChange is a
274 // pointer to G4ParticleMometum. Also it is a normalized
275 // momentum vector.
276
277 G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
278 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
279 G4double mass = theMassChange;
280
281 // Set Mass/Charge/MagneticMoment
282 pPostStepPoint->SetMass(theMassChange);
283 pPostStepPoint->SetCharge(theChargeChange);
285
286 // calculate new kinetic energy
287 G4double energy = pPostStepPoint->GetKineticEnergy()
288 + (theEnergyChange - pPreStepPoint->GetKineticEnergy());
289
290 // update kinetic energy and momentum direction
291 if (energy > 0.0) {
292 // calculate new momentum
293 G4ThreeVector pMomentum = pPostStepPoint->GetMomentum()
295 - pPreStepPoint->GetMomentum());
296 G4double tMomentum = pMomentum.mag();
297 G4ThreeVector direction(1.0,0.0,0.0);
298 if( tMomentum > 0. ){
299 G4double inv_Momentum= 1.0 / tMomentum;
300 direction= pMomentum * inv_Momentum;
301 }
302 pPostStepPoint->SetMomentumDirection(direction);
303 pPostStepPoint->SetKineticEnergy( energy );
304 } else {
305 // stop case
306 //pPostStepPoint->SetMomentumDirection(G4ThreeVector(1., 0., 0.));
307 pPostStepPoint->SetKineticEnergy(0.0);
308 }
309
310 if (!isVelocityChanged) theVelocityChange = pStep->GetTrack()->CalculateVelocity();
311 pPostStepPoint->SetVelocity(theVelocityChange);
312
313 // update polarization
315 - pPreStepPoint->GetPolarization());
316
317 // update position and time
318 pPostStepPoint->AddPosition( thePositionChange
319 - pPreStepPoint->GetPosition() );
320 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
321 pPostStepPoint->AddLocalTime( theTimeChange - theLocalTime0 );
322 pPostStepPoint->AddProperTime( theProperTimeChange
323 - pPreStepPoint->GetProperTime());
324
326 pPostStepPoint->SetWeight( theParentWeight );
327 }
328
329#ifdef G4VERBOSE
330 G4Track* aTrack = pStep->GetTrack();
331 if (debugFlag) CheckIt(*aTrack);
332#endif
333
334 // Update the G4Step specific attributes
335 return UpdateStepInfo(pStep);
336}
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)

◆ UpdateStepForAtRest()

G4Step * G4ParticleChange::UpdateStepForAtRest ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Reimplemented in G4ParticleChangeForTransport.

Definition at line 382 of file G4ParticleChange.cc.

383{
384 // A physics process always calculates the final state of the particle
385
386 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
387
388 // Set Mass/Charge
389 pPostStepPoint->SetMass(theMassChange);
390 pPostStepPoint->SetCharge(theChargeChange);
392
393 // update kinetic energy and momentum direction
395 pPostStepPoint->SetKineticEnergy( theEnergyChange );
396 if (!isVelocityChanged) theVelocityChange = pStep->GetTrack()->CalculateVelocity();
397 pPostStepPoint->SetVelocity(theVelocityChange);
398
399 // update polarization
400 pPostStepPoint->SetPolarization( thePolarizationChange );
401
402 // update position and time
403 pPostStepPoint->SetPosition( thePositionChange );
404 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
405 pPostStepPoint->SetLocalTime( theTimeChange );
406 pPostStepPoint->SetProperTime( theProperTimeChange );
407
409 pPostStepPoint->SetWeight( theParentWeight );
410 }
411
412#ifdef G4VERBOSE
413 G4Track* aTrack = pStep->GetTrack();
414 if (debugFlag) CheckIt(*aTrack);
415#endif
416
417 // Update the G4Step specific attributes
418 return UpdateStepInfo(pStep);
419}
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 338 of file G4ParticleChange.cc.

339{
340 // A physics process always calculates the final state of the particle
341
342 // Take note that the return type of GetMomentumChange is a
343 // pointer to G4ParticleMometum. Also it is a normalized
344 // momentum vector.
345
346 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
347
348 // Set Mass/Charge
349 pPostStepPoint->SetMass(theMassChange);
350 pPostStepPoint->SetCharge(theChargeChange);
352
353 // update kinetic energy and momentum direction
355 pPostStepPoint->SetKineticEnergy( theEnergyChange );
356 if (!isVelocityChanged) theVelocityChange = pStep->GetTrack()->CalculateVelocity();
357 pPostStepPoint->SetVelocity(theVelocityChange);
358
359 // update polarization
360 pPostStepPoint->SetPolarization( thePolarizationChange );
361
362 // update position and time
363 pPostStepPoint->SetPosition( thePositionChange );
364 pPostStepPoint->AddGlobalTime(theTimeChange - theLocalTime0);
365 pPostStepPoint->SetLocalTime( theTimeChange );
366 pPostStepPoint->SetProperTime( theProperTimeChange );
367
369 pPostStepPoint->SetWeight( theParentWeight );
370 }
371
372#ifdef G4VERBOSE
373 G4Track* aTrack = pStep->GetTrack();
374 if (debugFlag) CheckIt(*aTrack);
375#endif
376
377 // Update the G4Step specific attributes
378 return UpdateStepInfo(pStep);
379}

◆ UpdateStepInfo()

Member Data Documentation

◆ isVelocityChanged

◆ theChargeChange

G4double G4ParticleChange::theChargeChange
protected

◆ theCurrentTrack

const G4Track* G4ParticleChange::theCurrentTrack
protected

Definition at line 261 of file G4ParticleChange.hh.

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

◆ theEnergyChange

◆ theGlobalTime0

G4double G4ParticleChange::theGlobalTime0
protected

Definition at line 241 of file G4ParticleChange.hh.

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

◆ theLocalTime0

G4double G4ParticleChange::theLocalTime0
protected

◆ theMagneticMomentChange

G4double G4ParticleChange::theMagneticMomentChange
protected

◆ theMassChange

G4double G4ParticleChange::theMassChange
protected

◆ theMomentumDirectionChange

◆ thePolarizationChange

◆ thePositionChange

◆ theProperTimeChange

◆ theTimeChange

◆ theVelocityChange


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