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

#include <G4ParticleChangeForDecay.hh>

+ Inheritance diagram for G4ParticleChangeForDecay:

Public Member Functions

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

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

Protected Attributes

G4double theGlobalTime0
 
G4double theLocalTime0
 
G4double theTimeChange
 
G4ThreeVector thePolarizationChange
 
- 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 54 of file G4ParticleChangeForDecay.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForDecay() [1/2]

G4ParticleChangeForDecay::G4ParticleChangeForDecay ( )

Definition at line 48 of file G4ParticleChangeForDecay.cc.

51{
52#ifdef G4VERBOSE
53 if (verboseLevel>2) {
54 G4cout << "G4ParticleChangeForDecay::G4ParticleChangeForDecay() " << G4endl;
55 }
56#endif
57}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ ~G4ParticleChangeForDecay()

G4ParticleChangeForDecay::~G4ParticleChangeForDecay ( )
virtual

Definition at line 59 of file G4ParticleChangeForDecay.cc.

60{
61#ifdef G4VERBOSE
62 if (verboseLevel>2) {
63 G4cout << "G4ParticleChangeForDecay::~G4ParticleChangeForDecay() " << G4endl;
64 }
65#endif
66}

◆ G4ParticleChangeForDecay() [2/2]

G4ParticleChangeForDecay::G4ParticleChangeForDecay ( const G4ParticleChangeForDecay right)
protected

Member Function Documentation

◆ CheckIt()

G4bool G4ParticleChangeForDecay::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 193 of file G4ParticleChangeForDecay.cc.

194{
195 G4bool exitWithError = false;
196
197 G4double accuracy;
198
199 // local time should not go back
200 G4bool itsOK =true;
201 accuracy = -1.0*(theTimeChange - theLocalTime0)/ns;
202 if (accuracy > accuracyForWarning) {
203 itsOK = false;
204 exitWithError = (accuracy > accuracyForException);
205#ifdef G4VERBOSE
206 G4cout << " G4ParticleChangeForDecay::CheckIt : ";
207 G4cout << "the local time goes back !!"
208 << " Difference: " << accuracy << "[ns] " <<G4endl;
210 << " E=" << aTrack.GetKineticEnergy()/MeV
211 << " pos=" << aTrack.GetPosition().x()/m
212 << ", " << aTrack.GetPosition().y()/m
213 << ", " << aTrack.GetPosition().z()/m
214 <<G4endl;
215#endif
216 }
217
218 // dump out information of this particle change
219#ifdef G4VERBOSE
220 if (!itsOK) DumpInfo();
221#endif
222
223 // Exit with error
224 if (exitWithError) {
225 G4Exception("G4ParticleChangeForDecay::CheckIt",
226 "TRACK005",EventMustBeAborted,
227 "time was illegal");
228 }
229
230 // correction
231 if (!itsOK) {
232 theTimeChange = aTrack.GetLocalTime();
233 }
234
235 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
236 return itsOK;
237}
@ EventMustBeAborted
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
double z() const
double x() const
double y() const
const G4String & GetParticleName() const
const G4ThreeVector & GetPosition() const
G4double GetLocalTime() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
virtual G4bool CheckIt(const G4Track &)
static const G4double accuracyForException
static const G4double accuracyForWarning
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 UpdateStepForAtRest().

◆ DumpInfo()

void G4ParticleChangeForDecay::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 182 of file G4ParticleChangeForDecay.cc.

183{
184// Show header
186
187 G4int oldprc = G4cout.precision(3);
188 G4cout << " Time (ns) : "
189 << std::setw(20) << theTimeChange/ns << G4endl;
190 G4cout.precision(oldprc);
191}
int G4int
Definition: G4Types.hh:66
virtual void DumpInfo() const

Referenced by CheckIt().

◆ GetGlobalTime()

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

Definition at line 132 of file G4ParticleChangeForDecay.hh.

133{
134 // Convert the time delay to the global time.
135 return theGlobalTime0 + (theTimeChange-theLocalTime0) + timeDelay;
136}

Referenced by UpdateStepForAtRest().

◆ GetLocalTime()

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

Definition at line 145 of file G4ParticleChangeForDecay.hh.

146{
147 // Convert the time delay to the local time.
148 return theTimeChange + timeDelay;
149}

◆ GetPolarization()

const G4ThreeVector * G4ParticleChangeForDecay::GetPolarization ( ) const
inline

Definition at line 152 of file G4ParticleChangeForDecay.hh.

153{
154 return &thePolarizationChange;
155}

◆ Initialize()

void G4ParticleChangeForDecay::Initialize ( const G4Track track)
virtual

Reimplemented from G4VParticleChange.

Definition at line 124 of file G4ParticleChangeForDecay.cc.

125{
126 // use base class's method at first
128
129 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
130
131 // set TimeChange equal to local time of the parent track
132 theTimeChange = track.GetLocalTime();
133
134 // set initial Local/Global time of the parent track
135 theLocalTime0 = track.GetLocalTime();
137
138 // set the Polarization equal to those of the parent track
140}
const G4ThreeVector & GetPolarization() const
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void Initialize(const G4Track &)

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

◆ operator!=()

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

Definition at line 116 of file G4ParticleChangeForDecay.cc.

117{
118 return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
119}

◆ operator=()

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

Definition at line 77 of file G4ParticleChangeForDecay.cc.

78{
79 if (this != &right){
81#ifdef G4VERBOSE
82 if (verboseLevel>0) {
83 G4cout << "G4ParticleChangeForDecay: assignment operator Warning ";
84 G4cout << "theListOfSecondaries is not empty ";
85 }
86#endif
87 for (G4int index= 0; index<theNumberOfSecondaries; index++){
88 if ( (*theListOfSecondaries)[index] ) delete (*theListOfSecondaries)[index] ;
89 }
90 }
94 for (G4int index = 0; index<theNumberOfSecondaries; index++){
95 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
96 theListOfSecondaries->SetElement(index, newTrack); }
97
102
107 }
108 return *this;
109}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag

◆ operator==()

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

Definition at line 111 of file G4ParticleChangeForDecay.cc.

112{
113 return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
114}

◆ ProposeGlobalTime()

void G4ParticleChangeForDecay::ProposeGlobalTime ( G4double  t)
inline

Definition at line 126 of file G4ParticleChangeForDecay.hh.

127{
129}

Referenced by G4UnknownDecay::DecayIt().

◆ ProposeLocalTime()

void G4ParticleChangeForDecay::ProposeLocalTime ( G4double  t)
inline

Definition at line 139 of file G4ParticleChangeForDecay.hh.

140{
141 theTimeChange = t;
142}

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

◆ ProposePolarization() [1/2]

void G4ParticleChangeForDecay::ProposePolarization ( const G4ThreeVector finalPoralization)
inline

Definition at line 158 of file G4ParticleChangeForDecay.hh.

159{
160 thePolarizationChange = finalPoralization;
161}

◆ ProposePolarization() [2/2]

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

Definition at line 164 of file G4ParticleChangeForDecay.hh.

168{
172}
void setY(double)
void setZ(double)
void setX(double)

Referenced by G4DecayWithSpin::DecayIt().

◆ UpdateStepForAtRest()

G4Step * G4ParticleChangeForDecay::UpdateStepForAtRest ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 157 of file G4ParticleChangeForDecay.cc.

158{
159 // A physics process always calculates the final state of the particle
160
161 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
162
163 // update polarization
164 pPostStepPoint->SetPolarization( thePolarizationChange );
165
166 // update time
167 pPostStepPoint->SetGlobalTime( GetGlobalTime() );
168 pPostStepPoint->SetLocalTime( theTimeChange );
169 pPostStepPoint->AddProperTime (theTimeChange-theLocalTime0);
170
171#ifdef G4VERBOSE
172 G4Track* aTrack = pStep->GetTrack();
173 if (debugFlag) CheckIt(*aTrack);
174#endif
175
176 if (isParentWeightProposed )pPostStepPoint->SetWeight( theParentWeight );
177
178 // Update the G4Step specific attributes
179 return UpdateStepInfo(pStep);
180}
G4double GetGlobalTime(G4double timeDelay=0.0) const
virtual G4bool CheckIt(const G4Track &)
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)
virtual

Reimplemented from G4VParticleChange.

Definition at line 146 of file G4ParticleChangeForDecay.cc.

147{
149 pStep->GetPostStepPoint()->SetWeight( theParentWeight );
150 }
151
152 // Update the G4Step specific attributes
153 return UpdateStepInfo(pStep);
154}

Member Data Documentation

◆ theGlobalTime0

G4double G4ParticleChangeForDecay::theGlobalTime0
protected

◆ theLocalTime0

G4double G4ParticleChangeForDecay::theLocalTime0
protected

◆ thePolarizationChange

G4ThreeVector G4ParticleChangeForDecay::thePolarizationChange
protected

◆ theTimeChange


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