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

#include <G4ParticleChangeForDecay.hh>

+ Inheritance diagram for G4ParticleChangeForDecay:

Public Member Functions

 G4ParticleChangeForDecay ()
 
 ~G4ParticleChangeForDecay () override=default
 
 G4ParticleChangeForDecay (const G4ParticleChangeForDecay &right)=delete
 
G4ParticleChangeForDecayoperator= (const G4ParticleChangeForDecay &right)=delete
 
G4StepUpdateStepForAtRest (G4Step *Step) final
 
G4StepUpdateStepForPostStep (G4Step *Step) final
 
void Initialize (const G4Track &) final
 
void ProposeGlobalTime (G4double t)
 
void ProposeLocalTime (G4double t)
 
G4double GetGlobalTime (G4double timeDelay=0.0) const
 
G4double GetLocalTime (G4double timeDelay=0.0) const
 
const G4ThreeVectorGetPolarization () const
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
void ProposePolarization (const G4ThreeVector &finalPoralization)
 
void DumpInfo () const final
 
G4bool CheckIt (const G4Track &) final
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()=default
 
 G4VParticleChange (const G4VParticleChange &right)=delete
 
G4VParticleChangeoperator= (const G4VParticleChange &right)=delete
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
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
 
void ClearDebugFlag ()
 
void SetDebugFlag ()
 
G4bool GetDebugFlag () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VParticleChange
G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeLocalEnergyDeposit ()
 
void InitializeSteppingControl ()
 
void InitializeParentWeight (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries ()
 
void InitializeFromStep (const G4Step *)
 
G4double ComputeBeta (G4double kinEnergy)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 
- Protected Attributes inherited from G4VParticleChange
const G4TracktheCurrentTrack = nullptr
 
std::vector< G4Track * > theListOfSecondaries
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4double theTrueStepLength = 0.0
 
G4double theParentWeight = 1.0
 
G4double theParentGlobalTime = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4int theSizeOftheListOfSecondaries = 0
 
G4int verboseLevel = 1
 
G4int nError = 0
 
G4bool theFirstStepInVolume = false
 
G4bool theLastStepInVolume = false
 
G4bool isParentWeightProposed = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool debugFlag = false
 
- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 
static const G4int maxError = 10
 

Detailed Description

Definition at line 47 of file G4ParticleChangeForDecay.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForDecay() [1/2]

G4ParticleChangeForDecay::G4ParticleChangeForDecay ( )

Definition at line 39 of file G4ParticleChangeForDecay.cc.

40{}

◆ ~G4ParticleChangeForDecay()

G4ParticleChangeForDecay::~G4ParticleChangeForDecay ( )
overridedefault

◆ G4ParticleChangeForDecay() [2/2]

G4ParticleChangeForDecay::G4ParticleChangeForDecay ( const G4ParticleChangeForDecay & right)
delete

Member Function Documentation

◆ CheckIt()

G4bool G4ParticleChangeForDecay::CheckIt ( const G4Track & aTrack)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 121 of file G4ParticleChangeForDecay.cc.

122{
123 // local time should not go back
124 G4bool itsOK = true;
125 if(theTimeChange < theLocalTime0)
126 {
127 itsOK = false;
128 ++nError;
129#ifdef G4VERBOSE
130 if(nError < maxError)
131 {
132 G4cout << " G4ParticleChangeForDecay::CheckIt : ";
133 G4cout << "the local time goes back !!"
134 << " Difference: " << (theTimeChange - theLocalTime0) / ns
135 << "[ns] " << G4endl;
136 G4cout << "initial local time " << theLocalTime0 / ns << "[ns] "
137 << "initial global time " << theGlobalTime0 / ns << "[ns] "
138 << G4endl;
139 }
140#endif
141 theTimeChange = theLocalTime0;
142 }
143
144 if(!itsOK)
145 {
146 if(nError < maxError)
147 {
148#ifdef G4VERBOSE
149 // dump out information of this particle change
150 DumpInfo();
151#endif
152 G4Exception("G4ParticleChangeForDecay::CheckIt()", "TRACK005",
153 JustWarning, "time is illegal");
154 }
155 }
156
157 return (itsOK) && G4VParticleChange::CheckIt(aTrack);
158}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual G4bool CheckIt(const G4Track &)
static const G4int maxError

Referenced by UpdateStepForAtRest().

◆ DumpInfo()

void G4ParticleChangeForDecay::DumpInfo ( ) const
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 101 of file G4ParticleChangeForDecay.cc.

102{
103 // Show header
105
106 G4long oldprc = G4cout.precision(8);
107 G4cout << " -----------------------------------------------" << G4endl;
108 G4cout << " G4ParticleChangeForDecay proposes: " << G4endl;
109 G4cout << " Proposed local Time (ns): " << std::setw(20)
110 << theTimeChange / ns << G4endl;
111 G4cout << " Initial local Time (ns) : " << std::setw(20)
112 << theLocalTime0 / ns << G4endl;
113 G4cout << " Initial global Time (ns): " << std::setw(20)
114 << theGlobalTime0 / ns << G4endl;
115 G4cout << " Current global Time (ns): " << std::setw(20)
117 G4cout.precision(oldprc);
118}
long G4long
Definition G4Types.hh:87
G4double GetGlobalTime() const
virtual void DumpInfo() const
const G4Track * theCurrentTrack

Referenced by CheckIt().

◆ GetGlobalTime()

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

Definition at line 116 of file G4ParticleChangeForDecay.hh.

117{
118 // Convert the time delay to the global time.
119 return theGlobalTime0 + (theTimeChange - theLocalTime0) + timeDelay;
120}

Referenced by UpdateStepForAtRest().

◆ GetLocalTime()

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

Definition at line 129 of file G4ParticleChangeForDecay.hh.

130{
131 // Convert the time delay to the local time.
132 return theTimeChange + timeDelay;
133}

◆ GetPolarization()

const G4ThreeVector * G4ParticleChangeForDecay::GetPolarization ( ) const
inline

Definition at line 136 of file G4ParticleChangeForDecay.hh.

137{
138 return &thePolarizationChange;
139}

◆ Initialize()

void G4ParticleChangeForDecay::Initialize ( const G4Track & track)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 43 of file G4ParticleChangeForDecay.cc.

44{
45 // use base class's method at first
47
48 // set TimeChange equal to local time of the parent track
49 theLocalTime0 = theTimeChange = track.GetLocalTime();
50
51 // set initial Local/Global time of the parent track
52 theGlobalTime0 = track.GetGlobalTime();
53
54 // set the Polarization equal to those of the parent track
55 thePolarizationChange = track.GetDynamicParticle()->GetPolarization();
56}
const G4ThreeVector & GetPolarization() const
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void Initialize(const G4Track &)

Referenced by G4Decay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4UnknownDecay::DecayIt(), G4Decay::PostStepDoIt(), and G4DecayWithSpin::PostStepDoIt().

◆ operator=()

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

◆ ProposeGlobalTime()

void G4ParticleChangeForDecay::ProposeGlobalTime ( G4double t)
inline

Definition at line 110 of file G4ParticleChangeForDecay.hh.

111{
112 theTimeChange = (t - theGlobalTime0) + theLocalTime0;
113}

Referenced by G4UnknownDecay::DecayIt().

◆ ProposeLocalTime()

void G4ParticleChangeForDecay::ProposeLocalTime ( G4double t)
inline

Definition at line 123 of file G4ParticleChangeForDecay.hh.

124{
125 theTimeChange = t;
126}

Referenced by G4RadioactiveDecay::DecayAnalog(), G4Decay::DecayIt(), and G4Radioactivation::DecayIt().

◆ ProposePolarization() [1/2]

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

Definition at line 142 of file G4ParticleChangeForDecay.hh.

144{
145 thePolarizationChange = finalPoralization;
146}

◆ ProposePolarization() [2/2]

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

Definition at line 149 of file G4ParticleChangeForDecay.hh.

152{
153 thePolarizationChange.setX(Px);
154 thePolarizationChange.setY(Py);
155 thePolarizationChange.setZ(Pz);
156}
void setY(double)
void setZ(double)
void setX(double)

Referenced by G4DecayWithSpin::AtRestDoIt(), and G4DecayWithSpin::PostStepDoIt().

◆ UpdateStepForAtRest()

G4Step * G4ParticleChangeForDecay::UpdateStepForAtRest ( G4Step * Step)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 75 of file G4ParticleChangeForDecay.cc.

76{
77 // A physics process always calculates the final state of the particle
78
79 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
80
81 // update polarization
82 pPostStepPoint->SetPolarization(thePolarizationChange);
83
84 // update time
85 pPostStepPoint->SetGlobalTime(GetGlobalTime());
86 pPostStepPoint->SetLocalTime(theTimeChange);
87 pPostStepPoint->AddProperTime(theTimeChange - theLocalTime0);
88
89#ifdef G4VERBOSE
91#endif
92
94 pPostStepPoint->SetWeight(theParentWeight);
95
96 // Update the G4Step specific attributes
97 return UpdateStepInfo(pStep);
98}
G4double GetGlobalTime(G4double timeDelay=0.0) const
G4bool CheckIt(const G4Track &) final
void SetLocalTime(const G4double aValue)
void SetWeight(G4double aValue)
void AddProperTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4Step * UpdateStepInfo(G4Step *Step)

◆ UpdateStepForPostStep()

G4Step * G4ParticleChangeForDecay::UpdateStepForPostStep ( G4Step * Step)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 59 of file G4ParticleChangeForDecay.cc.

60{
61 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
62
64 {
65 pPostStepPoint->SetWeight(theParentWeight);
66 }
67 // update polarization
68 pPostStepPoint->SetPolarization(thePolarizationChange);
69
70 // update the G4Step specific attributes
71 return UpdateStepInfo(pStep);
72}

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