Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParticleChangeForDecay Class Reference

#include <G4ParticleChangeForDecay.hh>

Inheritance diagram for G4ParticleChangeForDecay:
Collaboration 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 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)
 
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
 

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::G4ParticleChangeForDecay ( )

Definition at line 48 of file G4ParticleChangeForDecay.cc.

49  : G4VParticleChange(),
51 {
52 #ifdef G4VERBOSE
53  if (verboseLevel>2) {
54  G4cout << "G4ParticleChangeForDecay::G4ParticleChangeForDecay() " << G4endl;
55  }
56 #endif
57 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
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 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForDecay::G4ParticleChangeForDecay ( const G4ParticleChangeForDecay right)
protected

Member Function Documentation

G4bool G4ParticleChangeForDecay::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 201 of file G4ParticleChangeForDecay.cc.

202 {
203  G4bool exitWithError = false;
204 
205  G4double accuracy;
206 
207  // local time should not go back
208  G4bool itsOK =true;
209  accuracy = -1.0*(theTimeChange - theLocalTime0)/ns;
210  if (accuracy > accuracyForWarning) {
211  itsOK = false;
212  exitWithError = (accuracy > accuracyForException);
213 #ifdef G4VERBOSE
214  G4cout << " G4ParticleChangeForDecay::CheckIt : ";
215  G4cout << "the local time goes back !!"
216  << " Difference: " << accuracy << "[ns] " <<G4endl;
217  G4cout << "initial local time "<< theLocalTime0/ns << "[ns] "
218  << "initial global time "<< theGlobalTime0/ns << "[ns] " <<G4endl;
219  G4cout << aTrack.GetDefinition()->GetParticleName()
220  << " E=" << aTrack.GetKineticEnergy()/MeV
221  << " pos=" << aTrack.GetPosition().x()/m
222  << ", " << aTrack.GetPosition().y()/m
223  << ", " << aTrack.GetPosition().z()/m
224  <<G4endl;
225 #endif
226  }
227 
228  // dump out information of this particle change
229 #ifdef G4VERBOSE
230  if (!itsOK) DumpInfo();
231 #endif
232 
233  // Exit with error
234  if (exitWithError) {
235  G4Exception("G4ParticleChangeForDecay::CheckIt",
236  "TRACK005",EventMustBeAborted,
237  "time was illegal");
238  }
239 
240  // correction
241  if (!itsOK) {
242  theTimeChange = aTrack.GetLocalTime();
243  }
244 
245  itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
246  return itsOK;
247 }
G4ParticleDefinition * GetDefinition() const
G4double GetLocalTime() const
double x() const
const G4ThreeVector & GetPosition() const
const G4String & GetParticleName() const
double z() const
static const G4double accuracyForException
virtual G4bool CheckIt(const G4Track &)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double y() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
static const G4double accuracyForWarning
#define ns
Definition: xmlparse.cc:614

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParticleChangeForDecay::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 186 of file G4ParticleChangeForDecay.cc.

187 {
188 // Show header
190 
191  G4int oldprc = G4cout.precision(3);
192  G4cout << " proposed local Time (ns) : "
193  << std::setw(20) << theTimeChange/ns << G4endl;
194  G4cout << " initial local Time (ns) : "
195  << std::setw(20) << theLocalTime0/ns << G4endl;
196  G4cout << " initial global Time (ns) : "
197  << std::setw(20) << theGlobalTime0/ns << G4endl;
198  G4cout.precision(oldprc);
199 }
virtual void DumpInfo() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
#define ns
Definition: xmlparse.cc:614

Here is the call graph for this function:

Here is the caller graph for this function:

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 }

Here is the caller graph for this function:

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 }
const G4ThreeVector * G4ParticleChangeForDecay::GetPolarization ( void  ) const
inline

Definition at line 152 of file G4ParticleChangeForDecay.hh.

153 {
154  return &thePolarizationChange;
155 }
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();
136  theGlobalTime0 = track.GetGlobalTime();
137 
138  // set the Polarization equal to those of the parent track
139  thePolarizationChange = pParticle->GetPolarization();
140 }
virtual void Initialize(const G4Track &)
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetGlobalTime() const
const G4ThreeVector & GetPolarization() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 116 of file G4ParticleChangeForDecay.cc.

117 {
118  return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
119 }
G4ParticleChangeForDecay & G4ParticleChangeForDecay::operator= ( const G4ParticleChangeForDecay right)
protected

Definition at line 77 of file G4ParticleChangeForDecay.cc.

78 {
79  if (this != &right){
80  if (theNumberOfSecondaries>0) {
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  }
91  delete theListOfSecondaries;
93  theNumberOfSecondaries = right.theNumberOfSecondaries;
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 }
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4TrackFastVector * theListOfSecondaries
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4SteppingControl theSteppingControlFlag
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
G4TrackStatus theStatusChange

Here is the call graph for this function:

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

Definition at line 111 of file G4ParticleChangeForDecay.cc.

112 {
113  return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
114 }
void G4ParticleChangeForDecay::ProposeGlobalTime ( G4double  t)
inline

Definition at line 126 of file G4ParticleChangeForDecay.hh.

Here is the caller graph for this function:

void G4ParticleChangeForDecay::ProposeLocalTime ( G4double  t)
inline

Definition at line 139 of file G4ParticleChangeForDecay.hh.

140 {
141  theTimeChange = t;
142 }

Here is the caller graph for this function:

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)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParticleChangeForDecay::ProposePolarization ( const G4ThreeVector finalPoralization)
inline

Definition at line 158 of file G4ParticleChangeForDecay.hh.

159 {
160  thePolarizationChange = finalPoralization;
161 }
G4Step * G4ParticleChangeForDecay::UpdateStepForAtRest ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 161 of file G4ParticleChangeForDecay.cc.

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

Here is the call graph for this function:

G4Step * G4ParticleChangeForDecay::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 146 of file G4ParticleChangeForDecay.cc.

147 {
148  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
149 
151  pPostStepPoint->SetWeight( theParentWeight );
152  }
153  // update polarization
154  pPostStepPoint->SetPolarization( thePolarizationChange );
155 
156  // Update the G4Step specific attributes
157  return UpdateStepInfo(pStep);
158 }
void SetWeight(G4double aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4Step * UpdateStepInfo(G4Step *Step)

Here is the call graph for this function:

Member Data Documentation

G4double G4ParticleChangeForDecay::theGlobalTime0
protected

Definition at line 109 of file G4ParticleChangeForDecay.hh.

G4double G4ParticleChangeForDecay::theLocalTime0
protected

Definition at line 111 of file G4ParticleChangeForDecay.hh.

G4ThreeVector G4ParticleChangeForDecay::thePolarizationChange
protected

Definition at line 117 of file G4ParticleChangeForDecay.hh.

G4double G4ParticleChangeForDecay::theTimeChange
protected

Definition at line 114 of file G4ParticleChangeForDecay.hh.


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