Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes
G4ParticleChangeForDecay Class Reference

#include <G4ParticleChangeForDecay.hh>

Inheritance diagram for G4ParticleChangeForDecay:
G4VParticleChange G4ParticleChangeForRadDecay

Public Member Functions

void AddSecondary (G4Track *aSecondary)
 
virtual G4bool CheckIt (const G4Track &)
 
void Clear ()
 
void ClearDebugFlag ()
 
virtual void DumpInfo () const
 
 G4ParticleChangeForDecay ()
 
G4bool GetDebugFlag () const
 
G4bool GetFirstStepInVolume () const
 
G4double GetGlobalTime (G4double timeDelay=0.0) const
 
G4bool GetLastStepInVolume () const
 
G4double GetLocalEnergyDeposit () const
 
G4double GetLocalTime (G4double timeDelay=0.0) const
 
G4double GetNonIonizingEnergyDeposit () const
 
G4int GetNumberOfSecondaries () const
 
G4double GetParentWeight () const
 
const G4ThreeVectorGetPolarization () const
 
G4TrackGetSecondary (G4int anIndex) const
 
G4SteppingControl GetSteppingControl () const
 
G4TrackStatus GetTrackStatus () const
 
G4double GetTrueStepLength () const
 
G4int GetVerboseLevel () const
 
G4double GetWeight () const
 
virtual void Initialize (const G4Track &)
 
G4bool IsParentWeightSetByProcess () const
 
G4bool IsSecondaryWeightSetByProcess () const
 
G4bool operator!= (const G4ParticleChangeForDecay &right) const
 
G4bool operator!= (const G4VParticleChange &right) const
 
G4bool operator== (const G4ParticleChangeForDecay &right) const
 
G4bool operator== (const G4VParticleChange &right) const
 
void ProposeFirstStepInVolume (G4bool flag)
 
void ProposeGlobalTime (G4double t)
 
void ProposeLastStepInVolume (G4bool flag)
 
void ProposeLocalEnergyDeposit (G4double anEnergyPart)
 
void ProposeLocalTime (G4double t)
 
void ProposeNonIonizingEnergyDeposit (G4double anEnergyPart)
 
void ProposeParentWeight (G4double finalWeight)
 
void ProposePolarization (const G4ThreeVector &finalPoralization)
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
void ProposeSteppingControl (G4SteppingControl StepControlFlag)
 
void ProposeTrackStatus (G4TrackStatus status)
 
void ProposeTrueStepLength (G4double truePathLength)
 
void ProposeWeight (G4double finalWeight)
 
void SetDebugFlag ()
 
void SetNumberOfSecondaries (G4int totSecondaries)
 
void SetParentWeightByProcess (G4bool)
 
void SetSecondaryWeightByProcess (G4bool)
 
void SetVerboseLevel (G4int vLevel)
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
virtual G4StepUpdateStepForAtRest (G4Step *Step)
 
virtual G4StepUpdateStepForPostStep (G4Step *Step)
 
virtual ~G4ParticleChangeForDecay ()
 

Protected Member Functions

G4bool CheckSecondary (G4Track &)
 
 G4ParticleChangeForDecay (const G4ParticleChangeForDecay &right)
 
G4double GetAccuracyForException () const
 
G4double GetAccuracyForWarning () const
 
void InitializeLocalEnergyDeposit (const G4Track &)
 
void InitializeParentGlobalTime (const G4Track &)
 
void InitializeParentWeight (const G4Track &)
 
void InitializeSecondaries (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeStepInVolumeFlags (const G4Track &)
 
void InitializeSteppingControl (const G4Track &)
 
void InitializeTrueStepLength (const G4Track &)
 
G4ParticleChangeForDecayoperator= (const G4ParticleChangeForDecay &right)
 
G4StepUpdateStepInfo (G4Step *Step)
 

Protected Attributes

G4bool debugFlag = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool isParentWeightProposed = false
 
G4bool theFirstStepInVolume = false
 
G4double theGlobalTime0 = 0.0
 
G4bool theLastStepInVolume = false
 
G4TrackFastVectortheListOfSecondaries = nullptr
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theLocalTime0 = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4double theParentGlobalTime = 0.0
 
G4double theParentWeight = 1.0
 
G4ThreeVector thePolarizationChange
 
G4int theSizeOftheListOfSecondaries = 0
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theTimeChange = 0.0
 
G4double theTrueStepLength = 0.0
 
G4int verboseLevel = 1
 

Static Protected Attributes

static const G4double accuracyForException = 0.001
 
static const G4double accuracyForWarning = 1.0e-9
 

Detailed Description

Definition at line 47 of file G4ParticleChangeForDecay.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForDecay() [1/2]

G4ParticleChangeForDecay::G4ParticleChangeForDecay ( )

Definition at line 40 of file G4ParticleChangeForDecay.cc.

◆ ~G4ParticleChangeForDecay()

G4ParticleChangeForDecay::~G4ParticleChangeForDecay ( )
virtual

Definition at line 46 of file G4ParticleChangeForDecay.cc.

47{
48}

◆ G4ParticleChangeForDecay() [2/2]

G4ParticleChangeForDecay::G4ParticleChangeForDecay ( const G4ParticleChangeForDecay right)
protected

Member Function Documentation

◆ AddSecondary()

void G4VParticleChange::AddSecondary ( G4Track aSecondary)
inherited

Definition at line 146 of file G4VParticleChange.cc.

147{
148 if(debugFlag)
149 CheckSecondary(*aTrack);
150
151 // add a secondary after size check
153 {
154 // set weight of secondary tracks
156 aTrack->SetWeight(theParentWeight);
159 }
160 else
161 {
162 delete aTrack;
163
164 G4Exception("G4VParticleChange::AddSecondary()", "TRACK101", JustWarning,
165 "Secondary buffer is full. The track is deleted!");
166 }
167}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
G4TrackFastVector * theListOfSecondaries
G4bool CheckSecondary(G4Track &)

References G4VParticleChange::CheckSecondary(), G4VParticleChange::debugFlag, G4VParticleChange::fSetSecondaryWeightByProcess, G4Exception(), JustWarning, G4FastVector< Type, N >::SetElement(), G4Track::SetWeight(), G4VParticleChange::theListOfSecondaries, G4VParticleChange::theNumberOfSecondaries, G4VParticleChange::theParentWeight, and G4VParticleChange::theSizeOftheListOfSecondaries.

Referenced by G4ParticleChangeForGamma::AddSecondary(), G4ParticleChange::AddSecondary(), G4eplusAnnihilation::AtRestDoIt(), G4FastStep::CreateSecondaryTrack(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4VEnergyLossProcess::FillSecondariesAlongStep(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), and G4ParticleChangeForOccurenceBiasing::StealSecondaries().

◆ CheckIt()

G4bool G4ParticleChangeForDecay::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 190 of file G4ParticleChangeForDecay.cc.

191{
192 G4bool exitWithError = false;
193
194 G4double accuracy;
195
196 // local time should not go back
197 G4bool itsOK = true;
198 accuracy = -1.0 * (theTimeChange - theLocalTime0) / ns;
199 if(accuracy > accuracyForWarning)
200 {
201 itsOK = false;
202 exitWithError = (accuracy > accuracyForException);
203#ifdef G4VERBOSE
204 G4cout << " G4ParticleChangeForDecay::CheckIt : ";
205 G4cout << "the local time goes back !!"
206 << " Difference: " << accuracy << "[ns] " << G4endl;
207 G4cout << "initial local time " << theLocalTime0 / ns << "[ns] "
208 << "initial global time " << theGlobalTime0 / ns << "[ns] "
209 << G4endl;
211 << " E=" << aTrack.GetKineticEnergy() / MeV
212 << " pos=" << aTrack.GetPosition().x() / m << ", "
213 << aTrack.GetPosition().y() / m << ", "
214 << aTrack.GetPosition().z() / m << 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 {
226 G4Exception("G4ParticleChangeForDecay::CheckIt()", "TRACK005",
227 EventMustBeAborted, "time was illegal");
228 }
229
230 // correction
231 if(!itsOK)
232 {
233 theTimeChange = aTrack.GetLocalTime();
234 }
235
236 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
237 return itsOK;
238}
@ EventMustBeAborted
static constexpr double m
Definition: G4SIunits.hh:109
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
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
#define ns
Definition: xmlparse.cc:614

References G4VParticleChange::accuracyForException, G4VParticleChange::accuracyForWarning, G4VParticleChange::CheckIt(), DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), G4Track::GetDefinition(), G4Track::GetKineticEnergy(), G4Track::GetLocalTime(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), m, MeV, ns, theGlobalTime0, theLocalTime0, theTimeChange, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by UpdateStepForAtRest().

◆ CheckSecondary()

G4bool G4VParticleChange::CheckSecondary ( G4Track aTrack)
protectedinherited

Definition at line 385 of file G4VParticleChange.cc.

386{
387 G4bool exitWithError = false;
388 G4double accuracy;
389 static G4ThreadLocal G4int nError = 0;
390#ifdef G4VERBOSE
391 const G4int maxError = 30;
392#endif
393
394 // MomentumDirection should be unit vector
395 G4bool itsOKforMomentum = true;
396 if(aTrack.GetKineticEnergy() > 0.)
397 {
398 accuracy = std::fabs((aTrack.GetMomentumDirection()).mag2() - 1.0);
399 if(accuracy > accuracyForWarning)
400 {
401 itsOKforMomentum = false;
402 nError += 1;
403 exitWithError = exitWithError || (accuracy > accuracyForException);
404#ifdef G4VERBOSE
405 if(nError < maxError)
406 {
407 G4cout << " G4VParticleChange::CheckSecondary : ";
408 G4cout << "the Momentum direction is not unit vector !! "
409 << " Difference: " << accuracy << G4endl;
411 << " E=" << aTrack.GetKineticEnergy() / MeV
412 << " pos=" << aTrack.GetPosition().x() / m << ", "
413 << aTrack.GetPosition().y() / m << ", "
414 << aTrack.GetPosition().z() / m << G4endl;
415 }
416#endif
417 }
418 }
419
420 // Kinetic Energy should not be negative
421 G4bool itsOKforEnergy = true;
422 accuracy = -1.0 * (aTrack.GetKineticEnergy()) / MeV;
423 if(accuracy > accuracyForWarning)
424 {
425 itsOKforEnergy = false;
426 nError += 1;
427 exitWithError = exitWithError || (accuracy > accuracyForException);
428#ifdef G4VERBOSE
429 if(nError < maxError)
430 {
431 G4cout << " G4VParticleChange::CheckSecondary : ";
432 G4cout << "the kinetic energy is negative !!"
433 << " Difference: " << accuracy << "[MeV] " << G4endl;
434 G4cout << " G4VParticleChange::CheckSecondary : ";
435 G4cout << "the global time of secondary is earlier than the parent !!"
436 << " Difference: " << accuracy << "[ns] " << G4endl;
438 << " E=" << aTrack.GetKineticEnergy() / MeV
439 << " pos=" << aTrack.GetPosition().x() / m << ", "
440 << aTrack.GetPosition().y() / m << ", "
441 << aTrack.GetPosition().z() / m << G4endl;
442 }
443#endif
444 }
445 // Check timing of secondaries
446 G4bool itsOKforTiming = true;
447
448 accuracy = (theParentGlobalTime - aTrack.GetGlobalTime()) / ns;
449 if(accuracy > accuracyForWarning)
450 {
451 itsOKforTiming = false;
452 nError += 1;
453 exitWithError = (accuracy > accuracyForException);
454#ifdef G4VERBOSE
455 if(nError < maxError)
456 {
457 G4cout << " G4VParticleChange::CheckSecondary : ";
458 G4cout
459 << "the global time of secondary goes back comapared to the parent !!"
460 << " Difference: " << accuracy << "[ns] " << G4endl;
462 << " E=" << aTrack.GetKineticEnergy() / MeV
463 << " pos=" << aTrack.GetPosition().x() / m << ", "
464 << aTrack.GetPosition().y() / m << ", "
465 << aTrack.GetPosition().z() / m
466 << " time=" << aTrack.GetGlobalTime() / ns
467 << " parent time=" << theParentGlobalTime / ns << G4endl;
468 }
469#endif
470 }
471
472 // Exit with error
473 if(exitWithError)
474 {
475 G4Exception("G4VParticleChange::CheckSecondary()", "TRACK001",
476 EventMustBeAborted, "Secondary with illegal energy/momentum ");
477 }
478
479 G4bool itsOK = itsOKforMomentum && itsOKforEnergy && itsOKforTiming;
480
481 // correction
482 if(!itsOKforMomentum)
483 {
484 G4double vmag = (aTrack.GetMomentumDirection()).mag();
485 aTrack.SetMomentumDirection((1. / vmag) * aTrack.GetMomentumDirection());
486 }
487 if(!itsOKforEnergy)
488 {
489 aTrack.SetKineticEnergy(0.0);
490 }
491
492 if(!itsOK)
493 {
494 this->DumpInfo();
495 }
496
497 return itsOK;
498}
int G4int
Definition: G4Types.hh:85
G4double GetGlobalTime() const
const G4ThreeVector & GetMomentumDirection() const
void SetKineticEnergy(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
virtual void DumpInfo() const
#define G4ThreadLocal
Definition: tls.hh:77

References G4VParticleChange::accuracyForException, G4VParticleChange::accuracyForWarning, G4VParticleChange::DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), G4ThreadLocal, G4Track::GetDefinition(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4Track::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), m, MeV, ns, G4Track::SetKineticEnergy(), G4Track::SetMomentumDirection(), G4VParticleChange::theParentGlobalTime, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4VParticleChange::AddSecondary().

◆ Clear()

void G4VParticleChange::Clear ( )
inherited

◆ ClearDebugFlag()

void G4VParticleChange::ClearDebugFlag ( )
inherited

◆ DumpInfo()

void G4ParticleChangeForDecay::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 174 of file G4ParticleChangeForDecay.cc.

175{
176 // Show header
178
179 G4int oldprc = G4cout.precision(3);
180 G4cout << " proposed local Time (ns) : " << std::setw(20)
181 << theTimeChange / ns << G4endl;
182 G4cout << " initial local Time (ns) : " << std::setw(20)
183 << theLocalTime0 / ns << G4endl;
184 G4cout << " initial global Time (ns) : " << std::setw(20)
185 << theGlobalTime0 / ns << G4endl;
186 G4cout.precision(oldprc);
187}

References G4VParticleChange::DumpInfo(), G4cout, G4endl, ns, theGlobalTime0, theLocalTime0, and theTimeChange.

Referenced by CheckIt().

◆ GetAccuracyForException()

G4double G4VParticleChange::GetAccuracyForException ( ) const
protectedinherited

Definition at line 507 of file G4VParticleChange.cc.

508{
510}

References G4VParticleChange::accuracyForException.

Referenced by G4FastStep::CheckIt().

◆ GetAccuracyForWarning()

G4double G4VParticleChange::GetAccuracyForWarning ( ) const
protectedinherited

Definition at line 501 of file G4VParticleChange.cc.

502{
503 return accuracyForWarning;
504}

References G4VParticleChange::accuracyForWarning.

Referenced by G4FastStep::CheckIt().

◆ GetDebugFlag()

G4bool G4VParticleChange::GetDebugFlag ( ) const
inherited

◆ GetFirstStepInVolume()

G4bool G4VParticleChange::GetFirstStepInVolume ( ) const
inherited

◆ GetGlobalTime()

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

Definition at line 123 of file G4ParticleChangeForDecay.hh.

124{
125 // Convert the time delay to the global time.
126 return theGlobalTime0 + (theTimeChange - theLocalTime0) + timeDelay;
127}

References theGlobalTime0, theLocalTime0, and theTimeChange.

Referenced by UpdateStepForAtRest().

◆ GetLastStepInVolume()

G4bool G4VParticleChange::GetLastStepInVolume ( ) const
inherited

◆ GetLocalEnergyDeposit()

G4double G4VParticleChange::GetLocalEnergyDeposit ( ) const
inherited

◆ GetLocalTime()

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

Definition at line 136 of file G4ParticleChangeForDecay.hh.

137{
138 // Convert the time delay to the local time.
139 return theTimeChange + timeDelay;
140}

References theTimeChange.

◆ GetNonIonizingEnergyDeposit()

G4double G4VParticleChange::GetNonIonizingEnergyDeposit ( ) const
inherited

◆ GetNumberOfSecondaries()

G4int G4VParticleChange::GetNumberOfSecondaries ( ) const
inherited

◆ GetParentWeight()

G4double G4VParticleChange::GetParentWeight ( ) const
inherited

◆ GetPolarization()

const G4ThreeVector * G4ParticleChangeForDecay::GetPolarization ( ) const
inline

Definition at line 143 of file G4ParticleChangeForDecay.hh.

144{
145 return &thePolarizationChange;
146}

References thePolarizationChange.

◆ GetSecondary()

G4Track * G4VParticleChange::GetSecondary ( G4int  anIndex) const
inherited

◆ GetSteppingControl()

G4SteppingControl G4VParticleChange::GetSteppingControl ( ) const
inherited

◆ GetTrackStatus()

G4TrackStatus G4VParticleChange::GetTrackStatus ( ) const
inherited

◆ GetTrueStepLength()

G4double G4VParticleChange::GetTrueStepLength ( ) const
inherited

◆ GetVerboseLevel()

G4int G4VParticleChange::GetVerboseLevel ( ) const
inherited

◆ GetWeight()

G4double G4VParticleChange::GetWeight ( ) const
inherited

◆ Initialize()

void G4ParticleChangeForDecay::Initialize ( const G4Track track)
virtual

Reimplemented from G4VParticleChange.

Definition at line 112 of file G4ParticleChangeForDecay.cc.

113{
114 // use base class's method at first
116
117 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
118
119 // set TimeChange equal to local time of the parent track
120 theTimeChange = track.GetLocalTime();
121
122 // set initial Local/Global time of the parent track
123 theLocalTime0 = track.GetLocalTime();
125
126 // set the Polarization equal to those of the parent track
128}
const G4ThreeVector & GetPolarization() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void Initialize(const G4Track &)

References G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4Track::GetLocalTime(), G4DynamicParticle::GetPolarization(), G4VParticleChange::Initialize(), theGlobalTime0, theLocalTime0, thePolarizationChange, and theTimeChange.

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

◆ InitializeLocalEnergyDeposit()

void G4VParticleChange::InitializeLocalEnergyDeposit ( const G4Track )
protectedinherited

◆ InitializeParentGlobalTime()

void G4VParticleChange::InitializeParentGlobalTime ( const G4Track )
protectedinherited

◆ InitializeParentWeight()

void G4VParticleChange::InitializeParentWeight ( const G4Track )
protectedinherited

◆ InitializeSecondaries()

void G4VParticleChange::InitializeSecondaries ( const G4Track )
protectedinherited

◆ InitializeStatusChange()

void G4VParticleChange::InitializeStatusChange ( const G4Track )
protectedinherited

◆ InitializeStepInVolumeFlags()

void G4VParticleChange::InitializeStepInVolumeFlags ( const G4Track )
protectedinherited

◆ InitializeSteppingControl()

void G4VParticleChange::InitializeSteppingControl ( const G4Track )
protectedinherited

◆ InitializeTrueStepLength()

void G4VParticleChange::InitializeTrueStepLength ( const G4Track )
protectedinherited

◆ IsParentWeightSetByProcess()

G4bool G4VParticleChange::IsParentWeightSetByProcess ( ) const
inherited

Definition at line 516 of file G4VParticleChange.cc.

516{ return true; }

◆ IsSecondaryWeightSetByProcess()

G4bool G4VParticleChange::IsSecondaryWeightSetByProcess ( ) const
inherited

◆ operator!=() [1/2]

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

Definition at line 105 of file G4ParticleChangeForDecay.cc.

107{
108 return ((G4VParticleChange*) this != (G4VParticleChange*) &right);
109}

◆ operator!=() [2/2]

G4bool G4VParticleChange::operator!= ( const G4VParticleChange right) const
inherited

Definition at line 140 of file G4VParticleChange.cc.

141{
142 return (this != (G4VParticleChange*) &right);
143}

◆ operator=()

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

Definition at line 63 of file G4ParticleChangeForDecay.cc.

64{
65 if(this != &right)
66 {
68 {
69 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
70 {
71 if((*theListOfSecondaries)[index])
72 delete(*theListOfSecondaries)[index];
73 }
74 }
78 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
79 {
80 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index]));
81 theListOfSecondaries->SetElement(index, newTrack);
82 }
83
88
93 }
94 return *this;
95}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
G4TrackStatus theStatusChange
G4SteppingControl theSteppingControlFlag

References G4FastVector< Type, N >::SetElement(), theGlobalTime0, G4VParticleChange::theListOfSecondaries, G4VParticleChange::theLocalEnergyDeposit, theLocalTime0, G4VParticleChange::theNumberOfSecondaries, thePolarizationChange, G4VParticleChange::theStatusChange, G4VParticleChange::theSteppingControlFlag, theTimeChange, and G4VParticleChange::theTrueStepLength.

◆ operator==() [1/2]

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

Definition at line 98 of file G4ParticleChangeForDecay.cc.

100{
101 return ((G4VParticleChange*) this == (G4VParticleChange*) &right);
102}

◆ operator==() [2/2]

G4bool G4VParticleChange::operator== ( const G4VParticleChange right) const
inherited

Definition at line 134 of file G4VParticleChange.cc.

135{
136 return (this == (G4VParticleChange*) &right);
137}

◆ ProposeFirstStepInVolume()

void G4VParticleChange::ProposeFirstStepInVolume ( G4bool  flag)
inherited

◆ ProposeGlobalTime()

void G4ParticleChangeForDecay::ProposeGlobalTime ( G4double  t)
inline

Definition at line 117 of file G4ParticleChangeForDecay.hh.

118{
120}

References theGlobalTime0, theLocalTime0, and theTimeChange.

Referenced by G4UnknownDecay::DecayIt().

◆ ProposeLastStepInVolume()

void G4VParticleChange::ProposeLastStepInVolume ( G4bool  flag)
inherited

◆ ProposeLocalEnergyDeposit()

void G4VParticleChange::ProposeLocalEnergyDeposit ( G4double  anEnergyPart)
inherited

Referenced by G4VEnergyLossProcess::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4OpBoundaryProcess::DoAbsorption(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4SpecialCuts::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4DNABornExcitationModel1::SampleSecondaries(), G4DNABornExcitationModel2::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNADiracRMatrixExcitationModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNAQuinnPlasmonExcitationModel::SampleSecondaries(), G4DNARelativisticIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4PAIPhotModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreIonisationModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4PenelopeBremsstrahlungModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopeIonisationModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4LEPTSAttachmentModel::SampleSecondaries(), G4LEPTSDissociationModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSExcitationModel::SampleSecondaries(), G4LEPTSIonisationModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4LEPTSRotExcitationModel::SampleSecondaries(), G4LEPTSVibExcitationModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ ProposeLocalTime()

void G4ParticleChangeForDecay::ProposeLocalTime ( G4double  t)
inline

◆ ProposeNonIonizingEnergyDeposit()

void G4VParticleChange::ProposeNonIonizingEnergyDeposit ( G4double  anEnergyPart)
inherited

◆ ProposeParentWeight()

void G4VParticleChange::ProposeParentWeight ( G4double  finalWeight)
inherited

◆ ProposePolarization() [1/2]

void G4ParticleChangeForDecay::ProposePolarization ( const G4ThreeVector finalPoralization)
inline

Definition at line 149 of file G4ParticleChangeForDecay.hh.

151{
152 thePolarizationChange = finalPoralization;
153}

References thePolarizationChange.

◆ ProposePolarization() [2/2]

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

◆ ProposeSteppingControl()

void G4VParticleChange::ProposeSteppingControl ( G4SteppingControl  StepControlFlag)
inherited

◆ ProposeTrackStatus()

void G4VParticleChange::ProposeTrackStatus ( G4TrackStatus  status)
inherited

Referenced by G4BiasingProcessInterface::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4OpBoundaryProcess::DoAbsorption(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4FastStep::KillPrimaryTrack(), G4ImportanceProcess::KillTrack(), G4WeightWindowProcess::KillTrack(), G4DNAElectronHoleRecombination::MakeReaction(), G4SpecialCuts::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4FastSimulationManagerProcess::PostStepDoIt(), G4PhononDownconversion::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNLoss::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4BiasingProcessInterface::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4PenelopeAnnihilationModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4eeToTwoGammaModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4eplusTo2GammaOKVIModel::SampleSecondaries(), G4eplusTo3GammaOKVIModel::SampleSecondaries(), G4eeToHadronsMultiModel::SampleSecondaries(), G4LEPTSAttachmentModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), and G4BetheHeitler5DModel::SampleSecondaries().

◆ ProposeTrueStepLength()

void G4VParticleChange::ProposeTrueStepLength ( G4double  truePathLength)
inherited

◆ ProposeWeight()

void G4VParticleChange::ProposeWeight ( G4double  finalWeight)
inherited

◆ SetDebugFlag()

void G4VParticleChange::SetDebugFlag ( )
inherited

◆ SetNumberOfSecondaries()

void G4VParticleChange::SetNumberOfSecondaries ( G4int  totSecondaries)
inherited

◆ SetParentWeightByProcess()

void G4VParticleChange::SetParentWeightByProcess ( G4bool  )
inherited

◆ SetSecondaryWeightByProcess()

void G4VParticleChange::SetSecondaryWeightByProcess ( G4bool  )
inherited

◆ SetVerboseLevel()

void G4VParticleChange::SetVerboseLevel ( G4int  vLevel)
inherited

◆ UpdateStepForAlongStep()

G4Step * G4VParticleChange::UpdateStepForAlongStep ( G4Step Step)
virtualinherited

◆ UpdateStepForAtRest()

G4Step * G4ParticleChangeForDecay::UpdateStepForAtRest ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 147 of file G4ParticleChangeForDecay.cc.

148{
149 // A physics process always calculates the final state of the particle
150
151 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
152
153 // update polarization
154 pPostStepPoint->SetPolarization(thePolarizationChange);
155
156 // update time
157 pPostStepPoint->SetGlobalTime(GetGlobalTime());
158 pPostStepPoint->SetLocalTime(theTimeChange);
159 pPostStepPoint->AddProperTime(theTimeChange - theLocalTime0);
160
161#ifdef G4VERBOSE
162 G4Track* aTrack = pStep->GetTrack();
163 if(debugFlag) { CheckIt(*aTrack); }
164#endif
165
167 pPostStepPoint->SetWeight(theParentWeight);
168
169 // Update the G4Step specific attributes
170 return UpdateStepInfo(pStep);
171}
G4double GetGlobalTime(G4double timeDelay=0.0) const
virtual G4bool CheckIt(const G4Track &)
void SetLocalTime(const G4double aValue)
void AddProperTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetPolarization(const G4ThreeVector &aValue)

References G4StepPoint::AddProperTime(), CheckIt(), G4VParticleChange::debugFlag, GetGlobalTime(), G4Step::GetPostStepPoint(), G4Step::GetTrack(), G4VParticleChange::isParentWeightProposed, G4StepPoint::SetGlobalTime(), G4StepPoint::SetLocalTime(), G4StepPoint::SetPolarization(), G4StepPoint::SetWeight(), theLocalTime0, G4VParticleChange::theParentWeight, thePolarizationChange, theTimeChange, and G4VParticleChange::UpdateStepInfo().

◆ UpdateStepForPostStep()

G4Step * G4ParticleChangeForDecay::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 131 of file G4ParticleChangeForDecay.cc.

132{
133 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
134
136 {
137 pPostStepPoint->SetWeight(theParentWeight);
138 }
139 // update polarization
140 pPostStepPoint->SetPolarization(thePolarizationChange);
141
142 // update the G4Step specific attributes
143 return UpdateStepInfo(pStep);
144}

References G4Step::GetPostStepPoint(), G4VParticleChange::isParentWeightProposed, G4StepPoint::SetPolarization(), G4StepPoint::SetWeight(), G4VParticleChange::theParentWeight, thePolarizationChange, and G4VParticleChange::UpdateStepInfo().

◆ UpdateStepInfo()

G4Step * G4VParticleChange::UpdateStepInfo ( G4Step Step)
protectedinherited

Definition at line 170 of file G4VParticleChange.cc.

171{
172 // Update the G4Step specific attributes
173 pStep->SetStepLength(theTrueStepLength);
174 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
175 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
176 pStep->SetControlFlag(theSteppingControlFlag);
177
179 {
180 pStep->SetFirstStepFlag();
181 }
182 else
183 {
184 pStep->ClearFirstStepFlag();
185 }
187 {
188 pStep->SetLastStepFlag();
189 }
190 else
191 {
192 pStep->ClearLastStepFlag();
193 }
194
195 return pStep;
196}
G4double theNonIonizingEnergyDeposit

References G4Step::AddNonIonizingEnergyDeposit(), G4Step::AddTotalEnergyDeposit(), G4Step::ClearFirstStepFlag(), G4Step::ClearLastStepFlag(), G4Step::SetControlFlag(), G4Step::SetFirstStepFlag(), G4Step::SetLastStepFlag(), G4Step::SetStepLength(), G4VParticleChange::theFirstStepInVolume, G4VParticleChange::theLastStepInVolume, G4VParticleChange::theLocalEnergyDeposit, G4VParticleChange::theNonIonizingEnergyDeposit, G4VParticleChange::theSteppingControlFlag, and G4VParticleChange::theTrueStepLength.

Referenced by G4VParticleChange::UpdateStepForAlongStep(), G4FastStep::UpdateStepForAtRest(), UpdateStepForAtRest(), G4VParticleChange::UpdateStepForAtRest(), G4FastStep::UpdateStepForPostStep(), UpdateStepForPostStep(), and G4VParticleChange::UpdateStepForPostStep().

Field Documentation

◆ accuracyForException

const G4double G4VParticleChange::accuracyForException = 0.001
staticprotectedinherited

◆ accuracyForWarning

const G4double G4VParticleChange::accuracyForWarning = 1.0e-9
staticprotectedinherited

◆ debugFlag

G4bool G4VParticleChange::debugFlag = false
protectedinherited

◆ fSetSecondaryWeightByProcess

G4bool G4VParticleChange::fSetSecondaryWeightByProcess = false
protectedinherited

◆ isParentWeightProposed

G4bool G4VParticleChange::isParentWeightProposed = false
protectedinherited

◆ theFirstStepInVolume

G4bool G4VParticleChange::theFirstStepInVolume = false
protectedinherited

◆ theGlobalTime0

G4double G4ParticleChangeForDecay::theGlobalTime0 = 0.0
protected

◆ theLastStepInVolume

G4bool G4VParticleChange::theLastStepInVolume = false
protectedinherited

◆ theListOfSecondaries

G4TrackFastVector* G4VParticleChange::theListOfSecondaries = nullptr
protectedinherited

◆ theLocalEnergyDeposit

G4double G4VParticleChange::theLocalEnergyDeposit = 0.0
protectedinherited

◆ theLocalTime0

G4double G4ParticleChangeForDecay::theLocalTime0 = 0.0
protected

◆ theNonIonizingEnergyDeposit

G4double G4VParticleChange::theNonIonizingEnergyDeposit = 0.0
protectedinherited

◆ theNumberOfSecondaries

G4int G4VParticleChange::theNumberOfSecondaries = 0
protectedinherited

◆ theParentGlobalTime

G4double G4VParticleChange::theParentGlobalTime = 0.0
protectedinherited

◆ theParentWeight

G4double G4VParticleChange::theParentWeight = 1.0
protectedinherited

◆ thePolarizationChange

G4ThreeVector G4ParticleChangeForDecay::thePolarizationChange
protected

◆ theSizeOftheListOfSecondaries

G4int G4VParticleChange::theSizeOftheListOfSecondaries = 0
protectedinherited

◆ theStatusChange

G4TrackStatus G4VParticleChange::theStatusChange = fAlive
protectedinherited

◆ theSteppingControlFlag

G4SteppingControl G4VParticleChange::theSteppingControlFlag = NormalCondition
protectedinherited

◆ theTimeChange

G4double G4ParticleChangeForDecay::theTimeChange = 0.0
protected

◆ theTrueStepLength

G4double G4VParticleChange::theTrueStepLength = 0.0
protectedinherited

◆ verboseLevel

G4int G4VParticleChange::verboseLevel = 1
protectedinherited

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