Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FastStep Class Reference

#include <G4FastStep.hh>

Inheritance diagram for G4FastStep:
Collaboration diagram for G4FastStep:

Public Member Functions

void KillPrimaryTrack ()
 
void ProposePrimaryTrackFinalPosition (const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalPosition (const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackFinalTime (G4double)
 
void SetPrimaryTrackFinalTime (G4double)
 
void ProposePrimaryTrackFinalProperTime (G4double)
 
void SetPrimaryTrackFinalProperTime (G4double)
 
void ProposePrimaryTrackFinalMomentumDirection (const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalMomentum (const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackFinalKineticEnergy (G4double)
 
void SetPrimaryTrackFinalKineticEnergy (G4double)
 
void ProposePrimaryTrackFinalKineticEnergyAndDirection (G4double, const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalKineticEnergyAndDirection (G4double, const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackFinalPolarization (const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalPolarization (const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackPathLength (G4double)
 
void SetPrimaryTrackPathLength (G4double)
 
void ProposePrimaryTrackFinalEventBiasingWeight (G4double)
 
void SetPrimaryTrackFinalEventBiasingWeight (G4double)
 
void SetNumberOfSecondaryTracks (G4int)
 
G4int GetNumberOfSecondaryTracks ()
 
G4TrackCreateSecondaryTrack (const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)
 
G4TrackCreateSecondaryTrack (const G4DynamicParticle &, G4ThreeVector, G4double, G4bool localCoordinates=true)
 
G4TrackGetSecondaryTrack (G4int)
 
void ProposeTotalEnergyDeposited (G4double anEnergyPart)
 
void SetTotalEnergyDeposited (G4double anEnergyPart)
 
G4double GetTotalEnergyDeposited () const
 
void ForceSteppingHitInvocation ()
 
 G4FastStep ()
 
virtual ~G4FastStep ()
 
G4bool operator== (const G4FastStep &right) const
 
G4bool operator!= (const G4FastStep &right) const
 
G4StepUpdateStepForAtRest (G4Step *Step)
 
G4StepUpdateStepForPostStep (G4Step *Step)
 
void Initialize (const G4FastTrack &)
 
void DumpInfo () const
 
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

 G4FastStep (const G4FastStep &right)
 
G4FastStepoperator= (const G4FastStep &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
 

Additional Inherited Members

- 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
 
- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 

Detailed Description

Definition at line 91 of file G4FastStep.hh.

Constructor & Destructor Documentation

G4FastStep::G4FastStep ( )

Definition at line 298 of file G4FastStep.cc.

299  : G4VParticleChange(),
300  theEnergyChange ( 0.0 ),
301  theTimeChange ( 0.0 ),
302  theProperTimeChange( 0.0 ),
303  fFastTrack ( nullptr ),
304  theWeightChange ( 0.0 )
305 {
306  if (verboseLevel>2)
307  {
308  G4cerr << "G4FastStep::G4FastStep()" << G4endl;
309  }
310 }
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr
G4FastStep::~G4FastStep ( )
virtual

Definition at line 312 of file G4FastStep.cc.

313 {
314  if (verboseLevel>2)
315  {
316  G4cerr << "G4FastStep::~G4FastStep()" << G4endl;
317  }
318 }
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr
G4FastStep::G4FastStep ( const G4FastStep right)
protected

Definition at line 321 of file G4FastStep.cc.

323 {
324  *this = right;
325 }

Member Function Documentation

G4bool G4FastStep::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 466 of file G4FastStep.cc.

467 {
468  //
469  // In the G4FastStep::CheckIt
470  // We only check a bit
471  //
472  // If the user violates the energy,
473  // We don't care, we agree.
474  //
475  // But for theMomentumDirectionChange,
476  // We do pay attention.
477  // And if too large is its range,
478  // We issue an Exception.
479  //
480  //
481  // It means, the G4FastStep::CheckIt issues an exception only for the
482  // theMomentumDirectionChange which should be an unit vector
483  // and it corrects it because it could cause problems for the ulterior
484  // tracking.For the rest, only warning are issued.
485 
486  G4bool itsOK = true;
487  G4bool exitWithError = false;
488  G4double accuracy;
489 
490  // Energy should not be larger than the initial value
491  accuracy = ( theEnergyChange - aTrack.GetKineticEnergy())/MeV;
492  if (accuracy > GetAccuracyForWarning())
493  {
495  ed << "The energy becomes larger than the initial value, difference = " << accuracy << " MeV" << G4endl;
496  G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
497  "FastSim006",
498  JustWarning, ed);
499  itsOK = false;
500  if (accuracy > GetAccuracyForException()) {exitWithError = true;}
501  }
502 
503  G4bool itsOKforMomentum = true;
504  if ( theEnergyChange >0.)
505  {
506  accuracy = std::abs(theMomentumChange.mag2()-1.0);
507  if (accuracy > GetAccuracyForWarning())
508  {
510  ed << "The Momentum Change is not a unit vector, difference = " << accuracy << G4endl;
511  G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
512  "FastSim007",
513  JustWarning, ed);
514  itsOK = itsOKforMomentum = false;
515  if (accuracy > GetAccuracyForException()) {exitWithError = true;}
516  }
517  }
518 
519  accuracy = (aTrack.GetGlobalTime()- theTimeChange)/ns;
520  if (accuracy > GetAccuracyForWarning())
521  {
523  ed << "The global time is getting backward, difference = " << accuracy << " ns" << G4endl;
524  G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
525  "FastSim008",
526  JustWarning, ed);
527  itsOK = false;
528  }
529 
530  accuracy = (aTrack.GetProperTime() - theProperTimeChange )/ns;
531  if (accuracy > GetAccuracyForWarning())
532  {
534  ed << "The proper time is getting backward, difference = " << accuracy << " ns" << G4endl;
535  G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
536  "FastSim009",
537  JustWarning, ed);
538  itsOK = false;
539  }
540 
541  if (!itsOK)
542  {
543  G4cout << "ERROR - G4FastStep::CheckIt() " << G4endl;
544  G4cout << " Pointer : " << this << G4endl ;
545  DumpInfo();
546  }
547 
548  // Exit with error
549  if (exitWithError)
550  {
552  ed << "An inaccuracy in G4FastStep is beyond tolerance." << G4endl;
553  G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
554  "FastSim010",
555  FatalException, ed);
556  }
557 
558  //correction for Momentum only.
559  if (!itsOKforMomentum) {
560  G4double vmag = theMomentumChange.mag();
561  theMomentumChange = (1./vmag)*theMomentumChange;
562  }
563 
564  itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
565  return itsOK;
566 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetProperTime() const
void DumpInfo() const
Definition: G4FastStep.cc:443
virtual G4bool CheckIt(const G4Track &)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
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
double mag2() const
G4double GetAccuracyForException() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
G4double GetAccuracyForWarning() const
double mag() const
#define ns
Definition: xmlparse.cc:614

Here is the call graph for this function:

Here is the caller graph for this function:

G4Track * G4FastStep::CreateSecondaryTrack ( const G4DynamicParticle dynamics,
G4ThreeVector  polarization,
G4ThreeVector  position,
G4double  time,
G4bool  localCoordinates = true 
)

Definition at line 203 of file G4FastStep.cc.

208 {
209  G4DynamicParticle dummyDynamics(dynamics);
210 
211  // ------------------------------------------
212  // Add the polarization to the dummyDynamics:
213  // ------------------------------------------
214  dummyDynamics.SetPolarization(polarization.x(),
215  polarization.y(),
216  polarization.z());
217 
218  return CreateSecondaryTrack(dummyDynamics, position, time, localCoordinates);
219 }
double x() const
G4Track * CreateSecondaryTrack(const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)
Definition: G4FastStep.cc:203
double z() const
double y() const

Here is the call graph for this function:

G4Track * G4FastStep::CreateSecondaryTrack ( const G4DynamicParticle dynamics,
G4ThreeVector  position,
G4double  time,
G4bool  localCoordinates = true 
)

Definition at line 225 of file G4FastStep.cc.

229 {
230  // ----------------------------------------
231  // Quantities in global coordinates system.
232  //
233  // The allocated globalDynamics is deleted
234  // by the destructor of the G4Track.
235  // ----------------------------------------
236  G4DynamicParticle* globalDynamics =
237  new G4DynamicParticle(dynamics);
238  G4ThreeVector globalPosition(position);
239 
240  // -----------------------------------
241  // Convert to global system if needed:
242  // -----------------------------------
243  if (localCoordinates)
244  {
245  // -- Momentum Direction:
246  globalDynamics->SetMomentumDirection(fFastTrack->
247  GetInverseAffineTransformation()->
248  TransformAxis(globalDynamics->
249  GetMomentumDirection()));
250  // -- Polarization:
251  G4ThreeVector globalPolarization;
252  globalPolarization = fFastTrack->GetInverseAffineTransformation()->
253  TransformAxis(globalDynamics->GetPolarization());
254  globalDynamics->SetPolarization(
255  globalPolarization.x(),
256  globalPolarization.y(),
257  globalPolarization.z()
258  );
259 
260  // -- Position:
261  globalPosition = fFastTrack->GetInverseAffineTransformation()->
262  TransformPoint(globalPosition);
263  }
264 
265  //-------------------------------------
266  // Create the G4Track of the secondary:
267  //-------------------------------------
268  G4Track* secondary = new G4Track(
269  globalDynamics,
270  time,
271  globalPosition
272  );
273 
274  //-------------------------------
275  // and feed the changes:
276  //-------------------------------
277  AddSecondary(secondary);
278 
279  //--------------------------------------
280  // returns the pointer on the secondary:
281  //--------------------------------------
282  return secondary;
283 }
double x() const
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:238
void SetMomentumDirection(const G4ThreeVector &aDirection)
double z() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void AddSecondary(G4Track *aSecondary)
double y() const
const G4ThreeVector & GetPolarization() const

Here is the call graph for this function:

void G4FastStep::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 443 of file G4FastStep.cc.

444 {
445 // use base-class DumpInfo
447 
448  G4cout << " Position - x (mm) : " << G4BestUnit( thePositionChange.x(), "Length" ) << G4endl;
449  G4cout << " Position - y (mm) : " << G4BestUnit( thePositionChange.y(), "Length" ) << G4endl;
450  G4cout << " Position - z (mm) : " << G4BestUnit( thePositionChange.z(), "Length" ) << G4endl;
451  G4cout << " Time (ns) : " << G4BestUnit( theTimeChange, "Time" ) << G4endl;
452  G4cout << " Proper Time (ns) : " << G4BestUnit( theProperTimeChange, "Time" ) << G4endl;
453  G4int olprc = G4cout.precision(3);
454  G4cout << " Momentum Direct - x : " << std::setw(20) << theMomentumChange.x() << G4endl;
455  G4cout << " Momentum Direct - y : " << std::setw(20) << theMomentumChange.y() << G4endl;
456  G4cout << " Momentum Direct - z : " << std::setw(20) << theMomentumChange.z() << G4endl;
457  G4cout.precision(olprc);
458  G4cout << " Kinetic Energy (MeV): " << G4BestUnit( theEnergyChange, "Energy" ) << G4endl;
459  G4cout.precision(3);
460  G4cout << " Polarization - x : " << std::setw(20) << thePolarizationChange.x() << G4endl;
461  G4cout << " Polarization - y : " << std::setw(20) << thePolarizationChange.y() << G4endl;
462  G4cout << " Polarization - z : " << std::setw(20) << thePolarizationChange.z() << G4endl;
463  G4cout.precision(olprc);
464 }
double x() const
virtual void DumpInfo() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
double z() const
G4GLOB_DLL std::ostream G4cout
double y() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FastStep::ForceSteppingHitInvocation ( )
G4int G4FastStep::GetNumberOfSecondaryTracks ( )
G4Track* G4FastStep::GetSecondaryTrack ( G4int  )
G4double G4FastStep::GetTotalEnergyDeposited ( ) const
void G4FastStep::Initialize ( const G4FastTrack fastTrack)

Definition at line 54 of file G4FastStep.cc.

55 {
56  // keeps the fastTrack reference
57  fFastTrack=&fastTrack;
58 
59  // currentTrack will be used to Initialize the other data members
60  const G4Track& currentTrack = *(fFastTrack->GetPrimaryTrack());
61 
62  // use base class's method at first
63  G4VParticleChange::Initialize(currentTrack);
64 
65  // set Energy/Momentum etc. equal to those of the parent particle
66  const G4DynamicParticle* pParticle = currentTrack.GetDynamicParticle();
67  theEnergyChange = pParticle->GetKineticEnergy();
68  theMomentumChange = pParticle->GetMomentumDirection();
69  thePolarizationChange = pParticle->GetPolarization();
70  theProperTimeChange = pParticle->GetProperTime();
71 
72  // set Position/Time etc. equal to those of the parent track
73  thePositionChange = currentTrack.GetPosition();
74  theTimeChange = currentTrack.GetGlobalTime();
75 
76  // switch off stepping hit invokation by default:
78 
79  // event biasing weigth:
80  theWeightChange = currentTrack.GetWeight();
81 }
virtual void Initialize(const G4Track &)
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetProperTime() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4SteppingControl theSteppingControlFlag
G4double GetGlobalTime() const
const G4ThreeVector & GetPolarization() const
G4double GetWeight() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FastStep::KillPrimaryTrack ( )

Definition at line 88 of file G4FastStep.cc.

89 {
92 }
void SetPrimaryTrackFinalKineticEnergy(G4double)
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

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

Definition at line 361 of file G4FastStep.cc.

362 {
363  return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
364 }
G4FastStep & G4FastStep::operator= ( const G4FastStep right)
protected

Definition at line 328 of file G4FastStep.cc.

329 {
330  if (this != &right)
331  {
337  theMomentumChange = right.theMomentumChange;
338  thePolarizationChange = right.thePolarizationChange;
339  thePositionChange = right.thePositionChange;
340  theProperTimeChange = right.theProperTimeChange;
341  theTimeChange = right.theTimeChange;
342  theEnergyChange = right.theEnergyChange;
346  theWeightChange = right.theWeightChange;
347  fFastTrack = right.fFastTrack;
348  }
349  return *this;
350 }
G4VParticleChange & operator=(const G4VParticleChange &right)
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag
G4TrackStatus theStatusChange

Here is the call graph for this function:

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

Definition at line 356 of file G4FastStep.cc.

357 {
358  return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
359 }
void G4FastStep::ProposePrimaryTrackFinalEventBiasingWeight ( G4double  )
void G4FastStep::ProposePrimaryTrackFinalKineticEnergy ( G4double  )
void G4FastStep::ProposePrimaryTrackFinalKineticEnergyAndDirection ( G4double  kineticEnergy,
const G4ThreeVector direction,
G4bool  localCoordinates = true 
)

Definition at line 151 of file G4FastStep.cc.

154 {
155  // Compute global direction if needed...
156  G4ThreeVector globalDirection = direction;
157  if (localCoordinates)
158  globalDirection =fFastTrack->GetInverseAffineTransformation()->
159  TransformAxis(direction);
160  // ...and feed the globalMomentum (ensuring unitarity)
161  SetMomentumChange(globalDirection.unit());
162  SetPrimaryTrackFinalKineticEnergy(kineticEnergy);
163 }
void SetPrimaryTrackFinalKineticEnergy(G4double)
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:238
Hep3Vector unit() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FastStep::ProposePrimaryTrackFinalMomentumDirection ( const G4ThreeVector momentum,
G4bool  localCoordinates = true 
)

Definition at line 125 of file G4FastStep.cc.

127 {
128  // Compute the momentum in global reference
129  // system if needed ...
130  G4ThreeVector globalMomentum = momentum;
131  if (localCoordinates)
132  globalMomentum = fFastTrack->GetInverseAffineTransformation()->
133  TransformAxis(momentum);
134  // ...and feed the globalMomentum (ensuring unitarity)
135  SetMomentumChange(globalMomentum.unit());
136 }
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:238
Hep3Vector unit() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FastStep::ProposePrimaryTrackFinalPolarization ( const G4ThreeVector polarization,
G4bool  localCoordinates = true 
)

Definition at line 179 of file G4FastStep.cc.

181 {
182  // Compute polarization in global system if needed:
183  G4ThreeVector globalPolarization(polarization);
184  if (localCoordinates)
185  globalPolarization = fFastTrack->GetInverseAffineTransformation()->
186  TransformAxis(globalPolarization);
187  // Feed the particle globalPolarization:
188  thePolarizationChange = globalPolarization;
189 }
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:238

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FastStep::ProposePrimaryTrackFinalPosition ( const G4ThreeVector position,
G4bool  localCoordinates = true 
)

Definition at line 99 of file G4FastStep.cc.

101 {
102  // Compute the position coordinate in global
103  // reference system if needed ...
104  G4ThreeVector globalPosition = position;
105  if (localCoordinates)
106  globalPosition = fFastTrack->GetInverseAffineTransformation()->
107  TransformPoint(position);
108  // ...and feed the globalPosition:
109  thePositionChange = globalPosition;
110 }
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:238
#define position
Definition: xmlparse.cc:622

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FastStep::ProposePrimaryTrackFinalProperTime ( G4double  )
void G4FastStep::ProposePrimaryTrackFinalTime ( G4double  )
void G4FastStep::ProposePrimaryTrackPathLength ( G4double  )
void G4FastStep::ProposeTotalEnergyDeposited ( G4double  anEnergyPart)
void G4FastStep::SetNumberOfSecondaryTracks ( G4int  )
void G4FastStep::SetPrimaryTrackFinalEventBiasingWeight ( G4double  )
void G4FastStep::SetPrimaryTrackFinalKineticEnergy ( G4double  )

Here is the caller graph for this function:

void G4FastStep::SetPrimaryTrackFinalKineticEnergyAndDirection ( G4double  kineticEnergy,
const G4ThreeVector direction,
G4bool  localCoordinates = true 
)

Definition at line 167 of file G4FastStep.cc.

170 {
171  ProposePrimaryTrackFinalKineticEnergyAndDirection(kineticEnergy, direction, localCoordinates);
172 }
void ProposePrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:151

Here is the call graph for this function:

void G4FastStep::SetPrimaryTrackFinalMomentum ( const G4ThreeVector momentum,
G4bool  localCoordinates = true 
)

Definition at line 140 of file G4FastStep.cc.

142 {
143  ProposePrimaryTrackFinalMomentumDirection(momentum, localCoordinates);
144 }
void ProposePrimaryTrackFinalMomentumDirection(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:125

Here is the call graph for this function:

void G4FastStep::SetPrimaryTrackFinalPolarization ( const G4ThreeVector polarization,
G4bool  localCoordinates = true 
)

Definition at line 193 of file G4FastStep.cc.

195 {
196  ProposePrimaryTrackFinalPolarization(polarization, localCoordinates);
197 }
void ProposePrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:179

Here is the call graph for this function:

void G4FastStep::SetPrimaryTrackFinalPosition ( const G4ThreeVector position,
G4bool  localCoordinates = true 
)

Definition at line 114 of file G4FastStep.cc.

116 {
117  ProposePrimaryTrackFinalPosition(position, localCoordinates);
118 }
void ProposePrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:99

Here is the call graph for this function:

void G4FastStep::SetPrimaryTrackFinalProperTime ( G4double  )
void G4FastStep::SetPrimaryTrackFinalTime ( G4double  )
void G4FastStep::SetPrimaryTrackPathLength ( G4double  )
void G4FastStep::SetTotalEnergyDeposited ( G4double  anEnergyPart)
G4Step * G4FastStep::UpdateStepForAtRest ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 407 of file G4FastStep.cc.

408 {
409  // A physics process always calculates the final state of the particle
410 
411  // G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
412  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
413  G4Track* aTrack = pStep->GetTrack();
414  // G4double mass = aTrack->GetDynamicParticle()->GetMass();
415 
416  // update kinetic energy and momentum direction
417  pPostStepPoint->SetMomentumDirection(theMomentumChange);
418  pPostStepPoint->SetKineticEnergy( theEnergyChange );
419 
420  // update polarization
421  pPostStepPoint->SetPolarization( thePolarizationChange );
422 
423  // update position and time
424  pPostStepPoint->SetPosition( thePositionChange );
425  pPostStepPoint->SetGlobalTime( theTimeChange );
426  pPostStepPoint->AddLocalTime( theTimeChange
427  - aTrack->GetGlobalTime());
428  pPostStepPoint->SetProperTime( theProperTimeChange );
429 
430  // update weight
431  pPostStepPoint->SetWeight( theWeightChange );
432 
433  if (debugFlag) CheckIt(*aTrack);
434 
435  // Update the G4Step specific attributes
436  return UpdateStepInfo(pStep);
437 }
void SetPosition(const G4ThreeVector &aValue)
void SetWeight(G4double aValue)
G4bool CheckIt(const G4Track &)
Definition: G4FastStep.cc:466
void AddLocalTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4double GetGlobalTime() const
void SetProperTime(const G4double aValue)
G4Step * UpdateStepInfo(G4Step *Step)
void SetKineticEnergy(const G4double aValue)

Here is the call graph for this function:

G4Step * G4FastStep::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 370 of file G4FastStep.cc.

371 {
372  // A physics process always calculates the final state of the particle
373 
374  // Take note that the return type of GetMomentumChange is a
375  // pointer to G4ParticleMometum. Also it is a normalized
376  // momentum vector.
377 
378  // G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
379  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
380  G4Track* aTrack = pStep->GetTrack();
381  // G4double mass = aTrack->GetDynamicParticle()->GetMass();
382 
383  // update kinetic energy and momentum direction
384  pPostStepPoint->SetMomentumDirection(theMomentumChange);
385  pPostStepPoint->SetKineticEnergy( theEnergyChange );
386 
387  // update polarization
388  pPostStepPoint->SetPolarization( thePolarizationChange );
389 
390  // update position and time
391  pPostStepPoint->SetPosition( thePositionChange );
392  pPostStepPoint->SetGlobalTime( theTimeChange );
393  pPostStepPoint->AddLocalTime( theTimeChange
394  - aTrack->GetGlobalTime());
395  pPostStepPoint->SetProperTime( theProperTimeChange );
396 
397  // update weight
398  pPostStepPoint->SetWeight( theWeightChange );
399 
400  if (debugFlag) CheckIt(*aTrack);
401 
402 
403  // Update the G4Step specific attributes
404  return UpdateStepInfo(pStep);
405 }
void SetPosition(const G4ThreeVector &aValue)
void SetWeight(G4double aValue)
G4bool CheckIt(const G4Track &)
Definition: G4FastStep.cc:466
void AddLocalTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4double GetGlobalTime() const
void SetProperTime(const G4double aValue)
G4Step * UpdateStepInfo(G4Step *Step)
void SetKineticEnergy(const G4double aValue)

Here is the call graph for this function:


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