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

#include <G4VParticleChange.hh>

Inheritance diagram for G4VParticleChange:
Collaboration diagram for G4VParticleChange:

Public Member Functions

 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

 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

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
 

Static Protected Attributes

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

Detailed Description

Definition at line 94 of file G4VParticleChange.hh.

Constructor & Destructor Documentation

G4VParticleChange::G4VParticleChange ( )

Definition at line 48 of file G4VParticleChange.cc.

56  theTrueStepLength(0.0),
57  theFirstStepInVolume(false),
58  theLastStepInVolume(false),
59  theParentWeight(1.0),
63  verboseLevel(1),
64  debugFlag(false)
65 {
66 #ifdef G4VERBOSE
67  // activate CHeckIt if in VERBOSE mode
68  debugFlag = true;
69 #endif
71 }
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
const G4int G4TrackFastVectorSize
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit
G4VParticleChange::~G4VParticleChange ( )
virtual

Definition at line 73 of file G4VParticleChange.cc.

73  {
74  // check if tracks still exist in theListOfSecondaries
75  if (theNumberOfSecondaries>0) {
76 #ifdef G4VERBOSE
77  if (verboseLevel>0) {
78  G4cout << "G4VParticleChange::~G4VParticleChange() Warning ";
79  G4cout << "theListOfSecondaries is not empty ";
80  }
81 #endif
82  for (G4int index= 0; index<theNumberOfSecondaries; index++){
83  delete (*theListOfSecondaries)[index] ;
84  }
85  }
86  delete theListOfSecondaries;
87 }
G4TrackFastVector * theListOfSecondaries
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4VParticleChange::G4VParticleChange ( const G4VParticleChange right)
protected

Definition at line 89 of file G4VParticleChange.cc.

101  isParentWeightProposed(false),
103  theParentGlobalTime(0.0),
104  verboseLevel(right.verboseLevel),
105  debugFlag(right.debugFlag)
106 {
109  for (G4int index = 0; index<theNumberOfSecondaries; index++){
110  G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
111  theListOfSecondaries->SetElement(index, newTrack);
112  }
113 }
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4TrackFastVector * theListOfSecondaries
int G4int
Definition: G4Types.hh:78
G4SteppingControl theSteppingControlFlag
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
const G4int G4TrackFastVectorSize
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:

Member Function Documentation

void G4VParticleChange::AddSecondary ( G4Track aSecondary)

Definition at line 171 of file G4VParticleChange.cc.

172 {
173  if (debugFlag) CheckSecondary(*aTrack);
174 
175  // add a secondary after size check
177  // Set weight of secondary tracks
178  if (!fSetSecondaryWeightByProcess) aTrack->SetWeight(theParentWeight);
181  } else {
182  delete aTrack;
183 
184 #ifdef G4VERBOSE
185  if (verboseLevel>0) {
186  G4cout << "G4VParticleChange::AddSecondary() Warning ";
187  G4cout << "theListOfSecondaries is full !! " << G4endl;
188  G4cout << " The track is deleted " << G4endl;
189  }
190 #endif
191  G4Exception("G4VParticleChange::AddSecondary",
192  "TRACK101", JustWarning,
193  "Secondary Bug is full. The track is deleted");
194  }
195 }
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4TrackFastVector * theListOfSecondaries
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G4bool CheckSecondary(G4Track &)

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4VParticleChange::CheckIt ( const G4Track )
virtual

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForLoss, G4ParticleChangeForMSC, and G4ParticleChangeForGamma.

Definition at line 316 of file G4VParticleChange.cc.

321 {
322 
323  G4bool exitWithError = false;
324  G4double accuracy;
325  static G4ThreadLocal G4int nError = 0;
326 #ifdef G4VERBOSE
327  const G4int maxError = 30;
328 #endif
329 
330  // Energy deposit should not be negative
331  G4bool itsOKforEnergy = true;
332  accuracy = -1.0*theLocalEnergyDeposit/MeV;
333  if (accuracy > accuracyForWarning) {
334  itsOKforEnergy = false;
335  nError += 1;
336  exitWithError = (accuracy > accuracyForException);
337 #ifdef G4VERBOSE
338  if (nError < maxError) {
339  G4cout << " G4VParticleChange::CheckIt : ";
340  G4cout << "the energy deposit is negative !!"
341  << " Difference: " << accuracy << "[MeV] " <<G4endl;
342  G4cout << aTrack.GetDefinition()->GetParticleName()
343  << " E=" << aTrack.GetKineticEnergy()/MeV
344  << " pos=" << aTrack.GetPosition().x()/m
345  << ", " << aTrack.GetPosition().y()/m
346  << ", " << aTrack.GetPosition().z()/m
347  <<G4endl;
348  }
349 #endif
350  }
351 
352  // true path length should not be negative
353  G4bool itsOKforStepLength = true;
354  accuracy = -1.0*theTrueStepLength/mm;
355  if (accuracy > accuracyForWarning) {
356  itsOKforStepLength = false;
357  nError += 1;
358  exitWithError = (accuracy > accuracyForException);
359 #ifdef G4VERBOSE
360  if (nError < maxError) {
361  G4cout << " G4VParticleChange::CheckIt : ";
362  G4cout << "the true step length is negative !!"
363  << " Difference: " << accuracy << "[MeV] " <<G4endl;
364  G4cout << aTrack.GetDefinition()->GetParticleName()
365  << " E=" << aTrack.GetKineticEnergy()/MeV
366  << " pos=" << aTrack.GetPosition().x()/m
367  << ", " << aTrack.GetPosition().y()/m
368  << ", " << aTrack.GetPosition().z()/m
369  <<G4endl;
370  }
371 #endif
372 
373  }
374 #ifdef G4VERBOSE
375  if (!itsOKforStepLength || !itsOKforEnergy ){
376  // dump out information of this particle change
377  DumpInfo();
378  }
379 #endif
380 
381  // Exit with error
382  if (exitWithError) {
383  G4Exception("G4VParticleChange::CheckIt",
384  "TRACK001", EventMustBeAborted,
385  "Step length and/or energy deposit was illegal");
386  }
387 
388  // correction
389  if ( !itsOKforStepLength ) {
390  theTrueStepLength = (1.e-12)*mm;
391  }
392  if ( !itsOKforEnergy ) {
393  theLocalEnergyDeposit = 0.0;
394  }
395  return (itsOKforStepLength && itsOKforEnergy );
396 }
static constexpr double mm
Definition: G4SIunits.hh:115
virtual void DumpInfo() const
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
static const G4double accuracyForException
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
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
static const G4double accuracyForWarning

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4VParticleChange::CheckSecondary ( G4Track aTrack)
protected

Definition at line 398 of file G4VParticleChange.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VParticleChange::Clear ( )

Here is the caller graph for this function:

void G4VParticleChange::ClearDebugFlag ( )

Here is the caller graph for this function:

void G4VParticleChange::DumpInfo ( ) const
virtual

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForLoss, G4ParticleChangeForTransport, G4ParticleChangeForDecay, G4ParticleChangeForMSC, and G4ParticleChangeForGamma.

Definition at line 252 of file G4VParticleChange.cc.

253 {
254 
255 // Show header
256  G4int olprc = G4cout.precision(3);
257  G4cout << " -----------------------------------------------"
258  << G4endl;
259  G4cout << " G4ParticleChange Information " << std::setw(20) << G4endl;
260  G4cout << " -----------------------------------------------"
261  << G4endl;
262 
263  G4cout << " # of 2ndaries : "
264  << std::setw(20) << theNumberOfSecondaries
265  << G4endl;
266 
267  if (theNumberOfSecondaries >0) {
268  G4cout << " Pointer to 2ndaries : "
269  << std::setw(20) << GetSecondary(0)
270  << G4endl;
271  G4cout << " (Showed only 1st one)"
272  << G4endl;
273  }
274  G4cout << " -----------------------------------------------"
275  << G4endl;
276 
277  G4cout << " Energy Deposit (MeV): "
278  << std::setw(20) << theLocalEnergyDeposit/MeV
279  << G4endl;
280 
281  G4cout << " Non-ionizing Energy Deposit (MeV): "
282  << std::setw(20) << theNonIonizingEnergyDeposit/MeV
283  << G4endl;
284 
285  G4cout << " Track Status : "
286  << std::setw(20);
287  if( theStatusChange == fAlive ){
288  G4cout << " Alive";
289  } else if( theStatusChange == fStopButAlive ){
290  G4cout << " StopButAlive";
291  } else if( theStatusChange == fStopAndKill ){
292  G4cout << " StopAndKill";
294  G4cout << " KillTrackAndSecondaries";
295  } else if( theStatusChange == fSuspend ){
296  G4cout << " Suspend";
297  } else if( theStatusChange == fPostponeToNextEvent ){
298  G4cout << " PostponeToNextEvent";
299  }
300  G4cout << G4endl;
301  G4cout << " True Path Length (mm) : "
302  << std::setw(20) << theTrueStepLength/mm
303  << G4endl;
304  G4cout << " Stepping Control : "
305  << std::setw(20) << theSteppingControlFlag
306  << G4endl;
307  if (theFirstStepInVolume) {
308  G4cout << " First Step In the voulme : " << G4endl;
309  }
310  if (theLastStepInVolume) {
311  G4cout << " Last Step In the voulme : " << G4endl;
312  }
313  G4cout.precision(olprc);
314 }
G4Track * GetSecondary(G4int anIndex) const
static constexpr double mm
Definition: G4SIunits.hh:115
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4SteppingControl theSteppingControlFlag
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VParticleChange::GetAccuracyForException ( ) const
protected

Definition at line 513 of file G4VParticleChange.cc.

514 {
515  return accuracyForException;
516 }
static const G4double accuracyForException

Here is the caller graph for this function:

G4double G4VParticleChange::GetAccuracyForWarning ( ) const
protected

Definition at line 508 of file G4VParticleChange.cc.

509 {
510  return accuracyForWarning;
511 }
static const G4double accuracyForWarning

Here is the caller graph for this function:

G4bool G4VParticleChange::GetDebugFlag ( ) const
G4bool G4VParticleChange::GetFirstStepInVolume ( ) const

Here is the caller graph for this function:

G4bool G4VParticleChange::GetLastStepInVolume ( ) const

Here is the caller graph for this function:

G4double G4VParticleChange::GetLocalEnergyDeposit ( ) const

Here is the caller graph for this function:

G4double G4VParticleChange::GetNonIonizingEnergyDeposit ( ) const
G4int G4VParticleChange::GetNumberOfSecondaries ( ) const

Here is the caller graph for this function:

G4double G4VParticleChange::GetParentWeight ( ) const

Here is the caller graph for this function:

G4Track* G4VParticleChange::GetSecondary ( G4int  anIndex) const

Here is the caller graph for this function:

G4SteppingControl G4VParticleChange::GetSteppingControl ( ) const
G4TrackStatus G4VParticleChange::GetTrackStatus ( ) const

Here is the caller graph for this function:

G4double G4VParticleChange::GetTrueStepLength ( ) const
G4int G4VParticleChange::GetVerboseLevel ( ) const
G4double G4VParticleChange::GetWeight ( ) const
virtual void G4VParticleChange::Initialize ( const G4Track )
virtual
void G4VParticleChange::InitializeLocalEnergyDeposit ( const G4Track )
protected
void G4VParticleChange::InitializeParentGlobalTime ( const G4Track )
protected
void G4VParticleChange::InitializeParentWeight ( const G4Track )
protected
void G4VParticleChange::InitializeSecondaries ( const G4Track )
protected

Here is the caller graph for this function:

void G4VParticleChange::InitializeStatusChange ( const G4Track )
protected
void G4VParticleChange::InitializeStepInVolumeFlags ( const G4Track )
protected
void G4VParticleChange::InitializeSteppingControl ( const G4Track )
protected
void G4VParticleChange::InitializeTrueStepLength ( const G4Track )
protected
G4bool G4VParticleChange::IsParentWeightSetByProcess ( ) const

Definition at line 527 of file G4VParticleChange.cc.

528 {
529  return true;
530 }
G4bool G4VParticleChange::IsSecondaryWeightSetByProcess ( ) const
G4bool G4VParticleChange::operator!= ( const G4VParticleChange right) const

Definition at line 166 of file G4VParticleChange.cc.

167 {
168  return (this != (G4VParticleChange *) &right);
169 }
G4VParticleChange & G4VParticleChange::operator= ( const G4VParticleChange right)
protected

Definition at line 116 of file G4VParticleChange.cc.

117 {
118  if (this != &right){
119  if (theNumberOfSecondaries>0) {
120 #ifdef G4VERBOSE
121  if (verboseLevel>0) {
122  G4cout << "G4VParticleChange: assignment operator Warning ";
123  G4cout << "theListOfSecondaries is not empty ";
124  }
125 #endif
126  for (G4int index = 0; index<theNumberOfSecondaries; index++){
127  if ( (*theListOfSecondaries)[index] ) delete ((*theListOfSecondaries)[index]) ;
128  }
129  }
130  delete theListOfSecondaries;
131 
133  theNumberOfSecondaries = right.theNumberOfSecondaries;
134  for (G4int index = 0; index<theNumberOfSecondaries; index++){
135  G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
136  theListOfSecondaries->SetElement(index, newTrack);
137  }
143 
146 
150 
152 
153  verboseLevel = right.verboseLevel;
154  debugFlag = right.debugFlag;
155 
156  }
157  return *this;
158 }
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
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 160 of file G4VParticleChange.cc.

161 {
162  return (this == (G4VParticleChange *) &right);
163 }
void G4VParticleChange::ProposeFirstStepInVolume ( G4bool  flag)

Here is the caller graph for this function:

void G4VParticleChange::ProposeLastStepInVolume ( G4bool  flag)

Here is the caller graph for this function:

void G4VParticleChange::ProposeLocalEnergyDeposit ( G4double  anEnergyPart)
void G4VParticleChange::ProposeNonIonizingEnergyDeposit ( G4double  anEnergyPart)

Here is the caller graph for this function:

void G4VParticleChange::ProposeParentWeight ( G4double  finalWeight)

Here is the caller graph for this function:

void G4VParticleChange::ProposeSteppingControl ( G4SteppingControl  StepControlFlag)

Here is the caller graph for this function:

void G4VParticleChange::ProposeTrackStatus ( G4TrackStatus  status)
void G4VParticleChange::ProposeTrueStepLength ( G4double  truePathLength)

Here is the caller graph for this function:

void G4VParticleChange::ProposeWeight ( G4double  finalWeight)

Here is the caller graph for this function:

void G4VParticleChange::SetDebugFlag ( )
void G4VParticleChange::SetNumberOfSecondaries ( G4int  totSecondaries)

Here is the caller graph for this function:

void G4VParticleChange::SetParentWeightByProcess ( G4bool  )

Definition at line 522 of file G4VParticleChange.cc.

523 {
524 }

Here is the caller graph for this function:

void G4VParticleChange::SetSecondaryWeightByProcess ( G4bool  )

Here is the caller graph for this function:

void G4VParticleChange::SetVerboseLevel ( G4int  vLevel)
G4Step * G4VParticleChange::UpdateStepForAlongStep ( G4Step Step)
virtual

Reimplemented in G4ParticleChange, G4ParticleChangeForMSC, G4ParticleChangeForTransport, G4ParticleChangeForLoss, G4ParticleChangeForOccurenceBiasing, and G4ParticleChangeForNothing.

Definition at line 228 of file G4VParticleChange.cc.

229 {
230  if (isParentWeightProposed ){
231  G4double initialWeight = Step->GetPreStepPoint()->GetWeight();
232  G4double currentWeight = Step->GetPostStepPoint()->GetWeight();
233  G4double finalWeight = (theParentWeight/initialWeight)*currentWeight;
234  Step->GetPostStepPoint()->SetWeight( finalWeight );
235  }
236  return UpdateStepInfo(Step);
237 }
G4double GetWeight() const
void SetWeight(G4double aValue)
G4StepPoint * GetPreStepPoint() const
G4Step * UpdateStepInfo(G4Step *Step)
G4StepPoint * GetPostStepPoint() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4Step * G4VParticleChange::UpdateStepForAtRest ( G4Step Step)
virtual

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForTransport, G4ParticleChangeForGamma, G4ParticleChangeForOccurenceBiasing, and G4ParticleChangeForNothing.

Definition at line 219 of file G4VParticleChange.cc.

220 {
223  }
224  return UpdateStepInfo(Step);
225 }
void SetWeight(G4double aValue)
G4Step * UpdateStepInfo(G4Step *Step)
G4StepPoint * GetPostStepPoint() const

Here is the call graph for this function:

G4Step * G4VParticleChange::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented in G4FastStep, G4ParticleChange, G4ParticleChangeForDecay, G4ParticleChangeForTransport, G4ParticleChangeForMSC, G4ParticleChangeForLoss, G4ParticleChangeForGamma, G4ParticleChangeForOccurenceBiasing, and G4ParticleChangeForNothing.

Definition at line 239 of file G4VParticleChange.cc.

240 {
241  if (isParentWeightProposed ){
243  }
244  return UpdateStepInfo(Step);
245 }
void SetWeight(G4double aValue)
G4Step * UpdateStepInfo(G4Step *Step)
G4StepPoint * GetPostStepPoint() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4Step * G4VParticleChange::UpdateStepInfo ( G4Step Step)
protected

Definition at line 202 of file G4VParticleChange.cc.

203 {
204  // Update the G4Step specific attributes
205  pStep->SetStepLength( theTrueStepLength );
206  pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
207  pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
208  pStep->SetControlFlag( theSteppingControlFlag );
209 
210  if (theFirstStepInVolume) {pStep->SetFirstStepFlag();}
211  else {pStep->ClearFirstStepFlag();}
212  if (theLastStepInVolume) {pStep->SetLastStepFlag();}
213  else {pStep->ClearLastStepFlag();}
214 
215  return pStep;
216 }
G4SteppingControl theSteppingControlFlag
G4double theNonIonizingEnergyDeposit

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

const G4double G4VParticleChange::accuracyForException = 0.001
staticprotected

Definition at line 318 of file G4VParticleChange.hh.

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

Definition at line 317 of file G4VParticleChange.hh.

G4bool G4VParticleChange::debugFlag
protected

Definition at line 314 of file G4VParticleChange.hh.

G4bool G4VParticleChange::fSetSecondaryWeightByProcess
protected

Definition at line 286 of file G4VParticleChange.hh.

G4bool G4VParticleChange::isParentWeightProposed
protected

Definition at line 284 of file G4VParticleChange.hh.

G4bool G4VParticleChange::theFirstStepInVolume
protected

Definition at line 278 of file G4VParticleChange.hh.

G4bool G4VParticleChange::theLastStepInVolume
protected

Definition at line 279 of file G4VParticleChange.hh.

G4TrackFastVector* G4VParticleChange::theListOfSecondaries
protected

Definition at line 245 of file G4VParticleChange.hh.

G4double G4VParticleChange::theLocalEnergyDeposit
protected

Definition at line 260 of file G4VParticleChange.hh.

G4double G4VParticleChange::theNonIonizingEnergyDeposit
protected

Definition at line 269 of file G4VParticleChange.hh.

G4int G4VParticleChange::theNumberOfSecondaries
protected

Definition at line 248 of file G4VParticleChange.hh.

G4double G4VParticleChange::theParentGlobalTime
protected

Definition at line 289 of file G4VParticleChange.hh.

G4double G4VParticleChange::theParentWeight
protected

Definition at line 282 of file G4VParticleChange.hh.

G4int G4VParticleChange::theSizeOftheListOfSecondaries
protected

Definition at line 251 of file G4VParticleChange.hh.

G4TrackStatus G4VParticleChange::theStatusChange
protected

Definition at line 254 of file G4VParticleChange.hh.

G4SteppingControl G4VParticleChange::theSteppingControlFlag
protected

Definition at line 257 of file G4VParticleChange.hh.

G4double G4VParticleChange::theTrueStepLength
protected

Definition at line 274 of file G4VParticleChange.hh.

G4int G4VParticleChange::verboseLevel
protected

Definition at line 293 of file G4VParticleChange.hh.


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