G4FastStep Class Reference

#include <G4FastStep.hh>

Inheritance diagram for G4FastStep:

G4VParticleChange

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 &)

Protected Member Functions

 G4FastStep (const G4FastStep &right)
G4FastStepoperator= (const G4FastStep &right)

Detailed Description

Definition at line 91 of file G4FastStep.hh.


Constructor & Destructor Documentation

G4FastStep::G4FastStep (  ) 

Definition at line 297 of file G4FastStep.cc.

References G4cerr, G4endl, and G4VParticleChange::verboseLevel.

00298   : G4VParticleChange()
00299 {
00300   if (verboseLevel>2)
00301   {
00302     G4cerr << "G4FastStep::G4FastStep() " << G4endl;
00303   }
00304 }

G4FastStep::~G4FastStep (  )  [virtual]

Definition at line 306 of file G4FastStep.cc.

References G4cerr, G4endl, and G4VParticleChange::verboseLevel.

00307 {
00308   if (verboseLevel>2)
00309   {
00310     G4cerr << "G4FastStep::~G4FastStep() " << G4endl;
00311   }
00312 }

G4FastStep::G4FastStep ( const G4FastStep right  )  [protected]

Definition at line 315 of file G4FastStep.cc.

00316   : G4VParticleChange()
00317 {
00318    *this = right;
00319 }


Member Function Documentation

G4bool G4FastStep::CheckIt ( const G4Track  )  [virtual]

Reimplemented from G4VParticleChange.

Definition at line 479 of file G4FastStep.cc.

References G4VParticleChange::CheckIt(), DumpInfo(), FatalException, G4cout, G4endl, G4Exception(), G4VParticleChange::GetAccuracyForException(), G4VParticleChange::GetAccuracyForWarning(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4Track::GetProperTime(), JustWarning, and ns.

Referenced by UpdateStepForAtRest(), and UpdateStepForPostStep().

00480 {
00481   //
00482   //      In the G4FastStep::CheckIt
00483   //      We only check a bit
00484   //      
00485   //      If the user violates the energy,
00486   //      We don't care, we agree.
00487   //
00488   //      But for theMomentumDirectionChange,
00489   //      We do pay attention.
00490   //      And if too large is its range,
00491   //      We issue an Exception.
00492   //
00493   //
00494   // It means, the G4FastStep::CheckIt issues an exception only for the
00495   // theMomentumDirectionChange which should be an unit vector
00496   // and it corrects it because it could cause problems for the ulterior
00497   // tracking.For the rest, only warning are issued.
00498 
00499   G4bool    itsOK = true;
00500   G4bool    exitWithError = false;
00501   G4double  accuracy;
00502   
00503   // Energy should not be larger than the initial value
00504   accuracy = ( theEnergyChange - aTrack.GetKineticEnergy())/MeV;
00505   if (accuracy > GetAccuracyForWarning())
00506     {
00507       G4ExceptionDescription ed;
00508       ed << "The energy becomes larger than the initial value, difference = " <<  accuracy << " MeV" << G4endl;
00509       G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
00510                   "FastSim006",
00511                   JustWarning, ed);
00512       itsOK = false;
00513       if (accuracy > GetAccuracyForException())  {exitWithError = true;}
00514     }
00515   
00516   G4bool itsOKforMomentum = true;
00517   if ( theEnergyChange >0.)
00518     {
00519       accuracy = std::abs(theMomentumChange.mag2()-1.0);
00520       if (accuracy > GetAccuracyForWarning())
00521         {
00522           G4ExceptionDescription ed;
00523           ed << "The Momentum Change is not a unit vector, difference = " <<  accuracy << G4endl;
00524           G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
00525                       "FastSim007",
00526                       JustWarning, ed);
00527           itsOK = itsOKforMomentum = false;
00528           if (accuracy > GetAccuracyForException())  {exitWithError = true;}
00529         }
00530     }
00531   
00532   accuracy = (aTrack.GetGlobalTime()- theTimeChange)/ns;  
00533   if (accuracy > GetAccuracyForWarning())
00534     {
00535       G4ExceptionDescription ed;
00536       ed << "The global time is getting backward, difference = " <<  accuracy << " ns" << G4endl;
00537       G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
00538                   "FastSim008",
00539                   JustWarning, ed);
00540       itsOK = false;
00541     }
00542   
00543   accuracy = (aTrack.GetProperTime() - theProperTimeChange )/ns;
00544   if (accuracy >  GetAccuracyForWarning())
00545     {
00546       G4ExceptionDescription ed;
00547       ed << "The proper time is getting backward, difference = " <<  accuracy << " ns" << G4endl;
00548       G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
00549                   "FastSim009",
00550                   JustWarning, ed);
00551       itsOK = false;
00552     }
00553   
00554   if (!itsOK)
00555     { 
00556       G4cout << "ERROR - G4FastStep::CheckIt() " << G4endl;
00557       G4cout << "        Pointer : " << this << G4endl ;
00558       DumpInfo();
00559     }
00560   
00561   // Exit with error
00562   if (exitWithError)
00563     {
00564       G4ExceptionDescription ed;
00565       ed << "An inaccuracy in G4FastStep is beyond tolerance." << G4endl;
00566       G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)",
00567                   "FastSim010",
00568                   FatalException, ed);
00569     }
00570   
00571   //correction for Momentum only.
00572   if (!itsOKforMomentum) {
00573     G4double vmag = theMomentumChange.mag();
00574     theMomentumChange = (1./vmag)*theMomentumChange;
00575   }
00576   
00577   itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack); 
00578   return itsOK;
00579 }

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

Definition at line 224 of file G4FastStep.cc.

References G4VParticleChange::AddSecondary(), G4FastTrack::GetInverseAffineTransformation(), G4DynamicParticle::GetPolarization(), G4DynamicParticle::SetMomentumDirection(), and G4DynamicParticle::SetPolarization().

00228 {
00229   // ----------------------------------------
00230   // Quantities in global coordinates system.
00231   //  
00232   // The allocated globalDynamics is deleted
00233   // by the destructor of the G4Track.
00234   // ----------------------------------------
00235   G4DynamicParticle* globalDynamics =
00236     new G4DynamicParticle(dynamics);
00237   G4ThreeVector globalPosition(position);
00238   
00239   // -----------------------------------
00240   // Convert to global system if needed:
00241   // -----------------------------------
00242   if (localCoordinates)
00243     {
00244       // -- Momentum Direction:
00245       globalDynamics->SetMomentumDirection(fFastTrack->
00246                                            GetInverseAffineTransformation()->
00247                                            TransformAxis(globalDynamics->
00248                                                          GetMomentumDirection()));
00249       // -- Polarization:
00250       G4ThreeVector globalPolarization;
00251       globalPolarization = fFastTrack->GetInverseAffineTransformation()->
00252         TransformAxis(globalDynamics->GetPolarization());
00253       globalDynamics->SetPolarization(
00254                                       globalPolarization.x(),
00255                                       globalPolarization.y(),
00256                                       globalPolarization.z()
00257                                       );
00258       
00259       // -- Position:
00260       globalPosition = fFastTrack->GetInverseAffineTransformation()->
00261         TransformPoint(globalPosition);
00262     }
00263   
00264   //-------------------------------------
00265   // Create the G4Track of the secondary:
00266   //-------------------------------------
00267   G4Track* secondary = new G4Track(
00268                                    globalDynamics,
00269                                    time,
00270                                    globalPosition
00271                                    );
00272   
00273   //-------------------------------
00274   // and feed the changes:
00275   //-------------------------------
00276   AddSecondary(secondary);
00277   
00278   //--------------------------------------
00279   // returns the pointer on the secondary:
00280   //--------------------------------------
00281   return secondary;
00282 }

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

Definition at line 202 of file G4FastStep.cc.

References G4DynamicParticle::SetPolarization().

00207 {
00208   G4DynamicParticle dummyDynamics(dynamics);
00209   
00210   // ------------------------------------------
00211   // Add the polarization to the dummyDynamics:
00212   // ------------------------------------------
00213   dummyDynamics.SetPolarization(polarization.x(),
00214                                 polarization.y(),
00215                                 polarization.z());
00216   
00217   return CreateSecondaryTrack(dummyDynamics, position, time, localCoordinates);
00218 }

void G4FastStep::DumpInfo (  )  const [virtual]

Reimplemented from G4VParticleChange.

Definition at line 435 of file G4FastStep.cc.

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

Referenced by CheckIt().

00436 {
00437 // use base-class DumpInfo
00438   G4VParticleChange::DumpInfo();
00439 
00440   G4cout.precision(3);
00441   G4cout << "        Position - x (mm)   : " 
00442        << std::setw(20) << thePositionChange.x()/mm
00443        << G4endl; 
00444   G4cout << "        Position - y (mm)   : " 
00445        << std::setw(20) << thePositionChange.y()/mm
00446        << G4endl; 
00447   G4cout << "        Position - z (mm)   : " 
00448        << std::setw(20) << thePositionChange.z()/mm
00449        << G4endl;
00450   G4cout << "        Time (ns)           : " 
00451        << std::setw(20) << theTimeChange/ns
00452        << G4endl;
00453   G4cout << "        Proper Time (ns)    : " 
00454        << std::setw(20) << theProperTimeChange/ns
00455        << G4endl;
00456   G4cout << "        Momentum Direct - x : " 
00457        << std::setw(20) << theMomentumChange.x()
00458        << G4endl;
00459   G4cout << "        Momentum Direct - y : " 
00460        << std::setw(20) << theMomentumChange.y()
00461        << G4endl;
00462   G4cout << "        Momentum Direct - z : " 
00463        << std::setw(20) << theMomentumChange.z()
00464        << G4endl;
00465   G4cout << "        Kinetic Energy (MeV): " 
00466        << std::setw(20) << theEnergyChange/MeV
00467        << G4endl;
00468   G4cout << "        Polarization - x    : " 
00469        << std::setw(20) << thePolarizationChange.x()
00470        << G4endl;
00471   G4cout << "        Polarization - y    : " 
00472        << std::setw(20) << thePolarizationChange.y()
00473        << G4endl;
00474   G4cout << "        Polarization - z    : " 
00475        << std::setw(20) <<  thePolarizationChange.z()
00476        << G4endl;
00477 }

void G4FastStep::ForceSteppingHitInvocation (  )  [inline]

Definition at line 123 of file G4FastStep.icc.

References NormalCondition, and G4VParticleChange::ProposeSteppingControl().

00124 {
00125   ProposeSteppingControl(NormalCondition);
00126 }

G4int G4FastStep::GetNumberOfSecondaryTracks (  )  [inline]

Definition at line 93 of file G4FastStep.icc.

References G4VParticleChange::GetNumberOfSecondaries().

00094 {
00095   return GetNumberOfSecondaries();
00096 }

G4Track * G4FastStep::GetSecondaryTrack ( G4int   )  [inline]

Definition at line 98 of file G4FastStep.icc.

References G4VParticleChange::GetSecondary().

00099 {
00100   return GetSecondary(i);
00101 }

G4double G4FastStep::GetTotalEnergyDeposited (  )  const [inline]

Definition at line 117 of file G4FastStep.icc.

References G4VParticleChange::GetLocalEnergyDeposit().

00118 {
00119   return GetLocalEnergyDeposit();
00120 }

void G4FastStep::Initialize ( const G4FastTrack  ) 

Definition at line 53 of file G4FastStep.cc.

References AvoidHitInvocation, G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4DynamicParticle::GetPolarization(), G4Track::GetPosition(), G4FastTrack::GetPrimaryTrack(), G4DynamicParticle::GetProperTime(), G4Track::GetWeight(), G4VParticleChange::Initialize(), and G4VParticleChange::theSteppingControlFlag.

Referenced by G4FastSimulationManager::AtRestGetFastSimulationManagerTrigger(), and G4FastSimulationManager::PostStepGetFastSimulationManagerTrigger().

00054 {
00055   // keeps the fastTrack reference
00056   fFastTrack=&fastTrack;
00057 
00058   // currentTrack will be used to Initialize the other data members
00059   const G4Track& currentTrack = *(fFastTrack->GetPrimaryTrack());
00060 
00061   // use base class's method at first
00062   G4VParticleChange::Initialize(currentTrack);
00063 
00064   // set Energy/Momentum etc. equal to those of the parent particle
00065   const G4DynamicParticle*  pParticle = currentTrack.GetDynamicParticle();
00066   theEnergyChange          = pParticle->GetKineticEnergy();
00067   theMomentumChange        = pParticle->GetMomentumDirection();
00068   thePolarizationChange    = pParticle->GetPolarization();
00069   theProperTimeChange      = pParticle->GetProperTime();
00070 
00071   // set Position/Time etc. equal to those of the parent track
00072   thePositionChange      = currentTrack.GetPosition();
00073   theTimeChange          = currentTrack.GetGlobalTime();
00074 
00075   // switch off stepping hit invokation by default:
00076   theSteppingControlFlag = AvoidHitInvocation;
00077 
00078   // event biasing weigth:
00079   theWeightChange        = currentTrack.GetWeight();
00080 }  

void G4FastStep::KillPrimaryTrack (  ) 

Definition at line 87 of file G4FastStep.cc.

References fStopAndKill, G4VParticleChange::ProposeTrackStatus(), and SetPrimaryTrackFinalKineticEnergy().

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

Definition at line 353 of file G4FastStep.cc.

00354 {
00355    return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
00356 }

G4FastStep & G4FastStep::operator= ( const G4FastStep right  )  [protected]

Definition at line 322 of file G4FastStep.cc.

References G4VParticleChange::operator=(), G4VParticleChange::theListOfSecondaries, G4VParticleChange::theLocalEnergyDeposit, G4VParticleChange::theNumberOfSecondaries, G4VParticleChange::theSizeOftheListOfSecondaries, G4VParticleChange::theStatusChange, G4VParticleChange::theSteppingControlFlag, and G4VParticleChange::theTrueStepLength.

00323 {
00324    if (this != &right)
00325    {
00326      G4VParticleChange::operator=(right);
00327      theListOfSecondaries          = right.theListOfSecondaries;
00328      theSizeOftheListOfSecondaries = right.theSizeOftheListOfSecondaries;
00329      theNumberOfSecondaries        = right.theNumberOfSecondaries;
00330      theStatusChange               = right.theStatusChange;
00331      theMomentumChange             = right.theMomentumChange;
00332      thePolarizationChange         = right.thePolarizationChange;
00333      thePositionChange             = right.thePositionChange;
00334      theTimeChange                 = right.theTimeChange;
00335      theEnergyChange               = right.theEnergyChange;
00336      theTrueStepLength             = right.theTrueStepLength;
00337      theLocalEnergyDeposit         = right.theLocalEnergyDeposit;
00338      theSteppingControlFlag        = right.theSteppingControlFlag;
00339      theWeightChange               = right.theWeightChange;
00340    }
00341    return *this;
00342 }

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

Definition at line 348 of file G4FastStep.cc.

00349 {
00350    return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
00351 }

void G4FastStep::ProposePrimaryTrackFinalEventBiasingWeight ( G4double   )  [inline]

Definition at line 148 of file G4FastStep.icc.

Referenced by SetPrimaryTrackFinalEventBiasingWeight().

00149 {
00150   theWeightChange = w;
00151 }

void G4FastStep::ProposePrimaryTrackFinalKineticEnergy ( G4double   )  [inline]

Definition at line 57 of file G4FastStep.icc.

Referenced by SetPrimaryTrackFinalKineticEnergy().

00058 {
00059   theEnergyChange = kineticEnergy;
00060 }

void G4FastStep::ProposePrimaryTrackFinalKineticEnergyAndDirection ( G4double  ,
const G4ThreeVector ,
G4bool  localCoordinates = true 
)

Definition at line 150 of file G4FastStep.cc.

References G4FastTrack::GetInverseAffineTransformation(), and SetPrimaryTrackFinalKineticEnergy().

Referenced by SetPrimaryTrackFinalKineticEnergyAndDirection().

00153 {
00154   // Compute global direction if needed...
00155   G4ThreeVector globalDirection = direction;
00156   if (localCoordinates)
00157     globalDirection =fFastTrack->GetInverseAffineTransformation()-> 
00158       TransformAxis(direction);
00159   // ...and feed the globalMomentum (ensuring unitarity)
00160   SetMomentumChange(globalDirection.unit());
00161   SetPrimaryTrackFinalKineticEnergy(kineticEnergy);
00162 }

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

Definition at line 124 of file G4FastStep.cc.

References G4FastTrack::GetInverseAffineTransformation().

Referenced by SetPrimaryTrackFinalMomentum().

00126 {
00127   // Compute the momentum in global reference
00128   // system if needed ...
00129   G4ThreeVector globalMomentum = momentum;
00130   if (localCoordinates)
00131     globalMomentum = fFastTrack->GetInverseAffineTransformation()->
00132       TransformAxis(momentum);
00133   // ...and feed the globalMomentum (ensuring unitarity)
00134   SetMomentumChange(globalMomentum.unit());
00135 }

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

Definition at line 178 of file G4FastStep.cc.

References G4FastTrack::GetInverseAffineTransformation().

Referenced by SetPrimaryTrackFinalPolarization().

00180 {
00181   // Compute polarization in global system if needed:
00182   G4ThreeVector globalPolarization(polarization);
00183   if (localCoordinates)
00184     globalPolarization = fFastTrack->GetInverseAffineTransformation()->
00185       TransformAxis(globalPolarization);  
00186   // Feed the particle globalPolarization:
00187   thePolarizationChange = globalPolarization;
00188 }

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

Definition at line 98 of file G4FastStep.cc.

References G4FastTrack::GetInverseAffineTransformation().

Referenced by SetPrimaryTrackFinalPosition().

00100 {
00101   // Compute the position coordinate in global
00102   // reference system if needed ...
00103   G4ThreeVector globalPosition = position;
00104   if (localCoordinates) 
00105     globalPosition = fFastTrack->GetInverseAffineTransformation()->
00106       TransformPoint(position);
00107   // ...and feed the globalPosition:
00108   thePositionChange = globalPosition;
00109 }

void G4FastStep::ProposePrimaryTrackFinalProperTime ( G4double   )  [inline]

Definition at line 44 of file G4FastStep.icc.

Referenced by SetPrimaryTrackFinalProperTime().

00045 {
00046   theProperTimeChange = properTime;
00047 }

void G4FastStep::ProposePrimaryTrackFinalTime ( G4double   )  [inline]

Definition at line 32 of file G4FastStep.icc.

Referenced by SetPrimaryTrackFinalTime().

00033 {
00034   theTimeChange = time;
00035 }

void G4FastStep::ProposePrimaryTrackPathLength ( G4double   )  [inline]

Definition at line 70 of file G4FastStep.icc.

References G4VParticleChange::ProposeTrueStepLength().

Referenced by SetPrimaryTrackPathLength().

00071 {
00072   ProposeTrueStepLength(pathLength);
00073 }

void G4FastStep::ProposeTotalEnergyDeposited ( G4double  anEnergyPart  )  [inline]

Definition at line 107 of file G4FastStep.icc.

References G4VParticleChange::ProposeLocalEnergyDeposit().

Referenced by SetTotalEnergyDeposited().

00108 {
00109   ProposeLocalEnergyDeposit(anEnergyPart);
00110 }

void G4FastStep::SetNumberOfSecondaryTracks ( G4int   )  [inline]

Definition at line 87 of file G4FastStep.icc.

References G4VParticleChange::SetNumberOfSecondaries().

00088 {
00089   SetNumberOfSecondaries(nSecondaries);
00090 }

void G4FastStep::SetPrimaryTrackFinalEventBiasingWeight ( G4double   )  [inline]

Definition at line 153 of file G4FastStep.icc.

References ProposePrimaryTrackFinalEventBiasingWeight().

void G4FastStep::SetPrimaryTrackFinalKineticEnergy ( G4double   )  [inline]

Definition at line 63 of file G4FastStep.icc.

References ProposePrimaryTrackFinalKineticEnergy().

Referenced by KillPrimaryTrack(), and ProposePrimaryTrackFinalKineticEnergyAndDirection().

00064 {
00065   ProposePrimaryTrackFinalKineticEnergy(kineticEnergy);
00066 }

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

Definition at line 166 of file G4FastStep.cc.

References ProposePrimaryTrackFinalKineticEnergyAndDirection().

00169 {
00170   ProposePrimaryTrackFinalKineticEnergyAndDirection(kineticEnergy, direction, localCoordinates);
00171 }

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

Definition at line 139 of file G4FastStep.cc.

References ProposePrimaryTrackFinalMomentumDirection().

00141 {
00142   ProposePrimaryTrackFinalMomentumDirection(momentum, localCoordinates);
00143 }

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

Definition at line 192 of file G4FastStep.cc.

References ProposePrimaryTrackFinalPolarization().

00194 {
00195   ProposePrimaryTrackFinalPolarization(polarization, localCoordinates);
00196 }

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

Definition at line 113 of file G4FastStep.cc.

References ProposePrimaryTrackFinalPosition().

00115 {
00116   ProposePrimaryTrackFinalPosition(position, localCoordinates);
00117 }

void G4FastStep::SetPrimaryTrackFinalProperTime ( G4double   )  [inline]

Definition at line 49 of file G4FastStep.icc.

References ProposePrimaryTrackFinalProperTime().

00050 {
00051   ProposePrimaryTrackFinalProperTime(properTime);
00052 }

void G4FastStep::SetPrimaryTrackFinalTime ( G4double   )  [inline]

Definition at line 37 of file G4FastStep.icc.

References ProposePrimaryTrackFinalTime().

00038 {
00039   ProposePrimaryTrackFinalTime(time);
00040 }

void G4FastStep::SetPrimaryTrackPathLength ( G4double   )  [inline]

Definition at line 75 of file G4FastStep.icc.

References ProposePrimaryTrackPathLength().

00076 {
00077   ProposePrimaryTrackPathLength(pathLength);
00078 }

void G4FastStep::SetTotalEnergyDeposited ( G4double  anEnergyPart  )  [inline]

Definition at line 111 of file G4FastStep.icc.

References ProposeTotalEnergyDeposited().

00112 {
00113   ProposeTotalEnergyDeposited(anEnergyPart);
00114 }

G4Step * G4FastStep::UpdateStepForAtRest ( G4Step Step  )  [virtual]

Reimplemented from G4VParticleChange.

Definition at line 399 of file G4FastStep.cc.

References G4StepPoint::AddLocalTime(), CheckIt(), G4VParticleChange::debugFlag, G4Track::GetGlobalTime(), G4Step::GetPostStepPoint(), G4Step::GetTrack(), G4StepPoint::SetGlobalTime(), G4StepPoint::SetKineticEnergy(), G4StepPoint::SetMomentumDirection(), G4StepPoint::SetPolarization(), G4StepPoint::SetPosition(), G4StepPoint::SetProperTime(), G4StepPoint::SetWeight(), and G4VParticleChange::UpdateStepInfo().

00400 { 
00401   // A physics process always calculates the final state of the particle
00402 
00403   // G4StepPoint* pPreStepPoint  = pStep->GetPreStepPoint(); 
00404   G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint(); 
00405   G4Track*     aTrack  = pStep->GetTrack();
00406   // G4double     mass = aTrack->GetDynamicParticle()->GetMass();
00407  
00408   // update kinetic energy and momentum direction
00409   pPostStepPoint->SetMomentumDirection(theMomentumChange);
00410   pPostStepPoint->SetKineticEnergy( theEnergyChange );
00411 
00412   // update polarization
00413   pPostStepPoint->SetPolarization( thePolarizationChange );
00414       
00415   // update position and time
00416   pPostStepPoint->SetPosition( thePositionChange  );
00417   pPostStepPoint->SetGlobalTime( theTimeChange  );
00418   pPostStepPoint->AddLocalTime( theTimeChange 
00419                                  - aTrack->GetGlobalTime());
00420   pPostStepPoint->SetProperTime( theProperTimeChange  );
00421 
00422   // update weight
00423   pPostStepPoint->SetWeight( theWeightChange );
00424 
00425   if (debugFlag) CheckIt(*aTrack);
00426 
00427   //  Update the G4Step specific attributes 
00428   return UpdateStepInfo(pStep);
00429 }

G4Step * G4FastStep::UpdateStepForPostStep ( G4Step Step  )  [virtual]

Reimplemented from G4VParticleChange.

Definition at line 362 of file G4FastStep.cc.

References G4StepPoint::AddLocalTime(), CheckIt(), G4VParticleChange::debugFlag, G4Track::GetGlobalTime(), G4Step::GetPostStepPoint(), G4Step::GetTrack(), G4StepPoint::SetGlobalTime(), G4StepPoint::SetKineticEnergy(), G4StepPoint::SetMomentumDirection(), G4StepPoint::SetPolarization(), G4StepPoint::SetPosition(), G4StepPoint::SetProperTime(), G4StepPoint::SetWeight(), and G4VParticleChange::UpdateStepInfo().

00363 { 
00364   // A physics process always calculates the final state of the particle
00365 
00366   // Take note that the return type of GetMomentumChange is a
00367   // pointer to G4ParticleMometum. Also it is a normalized 
00368   // momentum vector.
00369 
00370   //  G4StepPoint* pPreStepPoint  = pStep->GetPreStepPoint(); 
00371   G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint(); 
00372   G4Track*     aTrack  = pStep->GetTrack();
00373   //  G4double     mass = aTrack->GetDynamicParticle()->GetMass();
00374  
00375   // update kinetic energy and momentum direction
00376   pPostStepPoint->SetMomentumDirection(theMomentumChange);
00377   pPostStepPoint->SetKineticEnergy( theEnergyChange );
00378 
00379    // update polarization
00380   pPostStepPoint->SetPolarization( thePolarizationChange );
00381       
00382   // update position and time
00383   pPostStepPoint->SetPosition( thePositionChange  );
00384   pPostStepPoint->SetGlobalTime( theTimeChange  );
00385   pPostStepPoint->AddLocalTime( theTimeChange 
00386                                  - aTrack->GetGlobalTime());
00387   pPostStepPoint->SetProperTime( theProperTimeChange  );
00388 
00389   // update weight
00390   pPostStepPoint->SetWeight( theWeightChange );
00391 
00392   if (debugFlag) CheckIt(*aTrack);
00393 
00394   
00395   //  Update the G4Step specific attributes 
00396   return UpdateStepInfo(pStep);
00397 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:58 2013 for Geant4 by  doxygen 1.4.7