Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FastStep Class Reference

#include <G4FastStep.hh>

+ Inheritance diagram for G4FastStep:

Public Member Functions

 G4FastStep ()=default
 
 ~G4FastStep () override=default
 
 G4FastStep (const G4FastStep &right)=delete
 
G4FastStepoperator= (const G4FastStep &right)=delete
 
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 ()
 
G4StepUpdateStepForAtRest (G4Step *Step) override
 
G4StepUpdateStepForPostStep (G4Step *Step) override
 
void Initialize (const G4FastTrack &)
 
void DumpInfo () const override
 
G4bool CheckIt (const G4Track &) override
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()=default
 
 G4VParticleChange (const G4VParticleChange &right)=delete
 
G4VParticleChangeoperator= (const G4VParticleChange &right)=delete
 
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)
 
const G4TrackGetCurrentTrack () const
 
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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VParticleChange
G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeLocalEnergyDeposit ()
 
void InitializeSteppingControl ()
 
void InitializeParentWeight (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries ()
 
void InitializeFromStep (const G4Step *)
 
G4double ComputeBeta (G4double kinEnergy)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 
- Protected Attributes inherited from G4VParticleChange
const G4TracktheCurrentTrack = nullptr
 
std::vector< G4Track * > theListOfSecondaries
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4double theTrueStepLength = 0.0
 
G4double theParentWeight = 1.0
 
G4double theParentGlobalTime = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4int theSizeOftheListOfSecondaries = 0
 
G4int verboseLevel = 1
 
G4int nError = 0
 
G4bool theFirstStepInVolume = false
 
G4bool theLastStepInVolume = false
 
G4bool isParentWeightProposed = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool debugFlag = false
 
- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 
static const G4int maxError = 10
 

Detailed Description

Definition at line 89 of file G4FastStep.hh.

Constructor & Destructor Documentation

◆ G4FastStep() [1/2]

G4FastStep::G4FastStep ( )
default

◆ ~G4FastStep()

G4FastStep::~G4FastStep ( )
overridedefault

◆ G4FastStep() [2/2]

G4FastStep::G4FastStep ( const G4FastStep & right)
delete

Member Function Documentation

◆ CheckIt()

G4bool G4FastStep::CheckIt ( const G4Track & aTrack)
overridevirtual

Reimplemented from G4VParticleChange.

Definition at line 359 of file G4FastStep.cc.

360{
361 //
362 // In the G4FastStep::CheckIt
363 // We only check a bit
364 //
365 // If the user violates the energy,
366 // We don't care, we agree.
367 //
368 // But for theMomentumDirectionChange,
369 // We do pay attention.
370 // And if too large is its range,
371 // We issue an Exception.
372 //
373 //
374 // It means, the G4FastStep::CheckIt issues an exception only for the
375 // theMomentumDirectionChange which should be an unit vector
376 // and it corrects it because it could cause problems for the ulterior
377 // tracking.For the rest, only warning are issued.
378
379 G4bool itsOK = true;
380 G4bool exitWithError = false;
381 G4double accuracy;
382
383 // Energy should not be larger than the initial value
384 accuracy = (theEnergyChange - aTrack.GetKineticEnergy()) / MeV;
385 if (accuracy > GetAccuracyForWarning()) {
387 ed << "The energy becomes larger than the initial value, difference = " << accuracy << " MeV"
388 << G4endl;
389 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)", "FastSim006", JustWarning, ed);
390 itsOK = false;
391 if (accuracy > GetAccuracyForException()) {
392 exitWithError = true;
393 }
394 }
395
396 G4bool itsOKforMomentum = true;
397 if (theEnergyChange > 0.) {
398 accuracy = std::abs(theMomentumChange.mag2() - 1.0);
399 if (accuracy > GetAccuracyForWarning()) {
401 ed << "The Momentum Change is not a unit vector, difference = " << accuracy << G4endl;
402 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)", "FastSim007", JustWarning, ed);
403 itsOK = itsOKforMomentum = false;
404 if (accuracy > GetAccuracyForException()) {
405 exitWithError = true;
406 }
407 }
408 }
409
410 accuracy = (aTrack.GetGlobalTime() - theTimeChange) / ns;
411 if (accuracy > GetAccuracyForWarning()) {
413 ed << "The global time is getting backward, difference = " << accuracy << " ns" << G4endl;
414 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)", "FastSim008", JustWarning, ed);
415 itsOK = false;
416 }
417
418 accuracy = (aTrack.GetProperTime() - theProperTimeChange) / ns;
419 if (accuracy > GetAccuracyForWarning()) {
421 ed << "The proper time is getting backward, difference = " << accuracy << " ns" << G4endl;
422 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)", "FastSim009", JustWarning, ed);
423 itsOK = false;
424 }
425
426 if (!itsOK) {
427 G4cout << "ERROR - G4FastStep::CheckIt() " << G4endl;
428 G4cout << " Pointer : " << this << G4endl;
429 DumpInfo();
430 }
431
432 // Exit with error
433 if (exitWithError) {
435 ed << "An inaccuracy in G4FastStep is beyond tolerance." << G4endl;
436 G4Exception("G4FastStep::CheckIt(const G4Track& aTrack)", "FastSim010", FatalException, ed);
437 }
438
439 // correction for Momentum only.
440 if (!itsOKforMomentum) {
441 G4double vmag = theMomentumChange.mag();
442 theMomentumChange = (1. / vmag) * theMomentumChange;
443 }
444
445 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
446 return itsOK;
447}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double mag2() const
double mag() const
void DumpInfo() const override
G4double GetGlobalTime() const
G4double GetProperTime() const
G4double GetKineticEnergy() const
virtual G4bool CheckIt(const G4Track &)
G4double GetAccuracyForException() const
G4double GetAccuracyForWarning() const
#define ns(x)
Definition xmltok.c:1649

Referenced by UpdateStepForAtRest(), and UpdateStepForPostStep().

◆ CreateSecondaryTrack() [1/2]

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

Definition at line 199 of file G4FastStep.cc.

201{
202 // ----------------------------------------
203 // Quantities in global coordinates system.
204 //
205 // The allocated globalDynamics is deleted
206 // by the destructor of the G4Track.
207 // ----------------------------------------
208 auto globalDynamics = new G4DynamicParticle(dynamics);
209 G4ThreeVector globalPosition(position);
210
211 // -----------------------------------
212 // Convert to global system if needed:
213 // -----------------------------------
214 if (localCoordinates) {
215 // -- Momentum Direction:
216 globalDynamics->SetMomentumDirection(
218 globalDynamics->GetMomentumDirection()));
219 // -- Polarization:
220 G4ThreeVector globalPolarization;
221 globalPolarization = fFastTrack->GetInverseAffineTransformation()->TransformAxis(
222 globalDynamics->GetPolarization());
223 globalDynamics->SetPolarization(globalPolarization.x(), globalPolarization.y(),
224 globalPolarization.z());
225
226 // -- Position:
227 globalPosition = fFastTrack->GetInverseAffineTransformation()->TransformPoint(globalPosition);
228 }
229
230 //-------------------------------------
231 // Create the G4Track of the secondary:
232 //-------------------------------------
233 auto secondary = new G4Track(globalDynamics, time, globalPosition);
234
235 //-------------------------------
236 // and feed the changes:
237 //-------------------------------
238 AddSecondary(secondary);
239
240 //--------------------------------------
241 // returns the pointer on the secondary:
242 //--------------------------------------
243 return secondary;
244}
double z() const
double x() const
double y() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
const G4AffineTransform * GetInverseAffineTransformation() const
void AddSecondary(G4Track *aSecondary)

◆ CreateSecondaryTrack() [2/2]

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

Definition at line 182 of file G4FastStep.cc.

185{
186 G4DynamicParticle dummyDynamics(dynamics);
187
188 // ------------------------------------------
189 // Add the polarization to the dummyDynamics:
190 // ------------------------------------------
191 dummyDynamics.SetPolarization(polarization.x(), polarization.y(), polarization.z());
192
193 return CreateSecondaryTrack(dummyDynamics, position, time, localCoordinates);
194}
G4Track * CreateSecondaryTrack(const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)

Referenced by CreateSecondaryTrack(), and G4BaierKatkov::GeneratePhoton().

◆ DumpInfo()

void G4FastStep::DumpInfo ( ) const
overridevirtual

Reimplemented from G4VParticleChange.

Definition at line 330 of file G4FastStep.cc.

331{
332 // use base-class DumpInfo
334
335 G4cout << " Position - x (mm) : " << G4BestUnit(thePositionChange.x(), "Length")
336 << G4endl;
337 G4cout << " Position - y (mm) : " << G4BestUnit(thePositionChange.y(), "Length")
338 << G4endl;
339 G4cout << " Position - z (mm) : " << G4BestUnit(thePositionChange.z(), "Length")
340 << G4endl;
341 G4cout << " Time (ns) : " << G4BestUnit(theTimeChange, "Time") << G4endl;
342 G4cout << " Proper Time (ns) : " << G4BestUnit(theProperTimeChange, "Time") << G4endl;
343 G4long olprc = G4cout.precision(3);
344 G4cout << " Momentum Direct - x : " << std::setw(20) << theMomentumChange.x() << G4endl;
345 G4cout << " Momentum Direct - y : " << std::setw(20) << theMomentumChange.y() << G4endl;
346 G4cout << " Momentum Direct - z : " << std::setw(20) << theMomentumChange.z() << G4endl;
347 G4cout.precision(olprc);
348 G4cout << " Kinetic Energy (MeV): " << G4BestUnit(theEnergyChange, "Energy") << G4endl;
349 G4cout.precision(3);
350 G4cout << " Polarization - x : " << std::setw(20) << thePolarizationChange.x()
351 << G4endl;
352 G4cout << " Polarization - y : " << std::setw(20) << thePolarizationChange.y()
353 << G4endl;
354 G4cout << " Polarization - z : " << std::setw(20) << thePolarizationChange.z()
355 << G4endl;
356 G4cout.precision(olprc);
357}
#define G4BestUnit(a, b)
long G4long
Definition G4Types.hh:87
virtual void DumpInfo() const

Referenced by CheckIt().

◆ ForceSteppingHitInvocation()

void G4FastStep::ForceSteppingHitInvocation ( )

◆ GetNumberOfSecondaryTracks()

G4int G4FastStep::GetNumberOfSecondaryTracks ( )

◆ GetSecondaryTrack()

G4Track * G4FastStep::GetSecondaryTrack ( G4int )

◆ GetTotalEnergyDeposited()

G4double G4FastStep::GetTotalEnergyDeposited ( ) const

◆ Initialize()

void G4FastStep::Initialize ( const G4FastTrack & fastTrack)

Definition at line 53 of file G4FastStep.cc.

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

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

◆ KillPrimaryTrack()

void G4FastStep::KillPrimaryTrack ( )

Definition at line 87 of file G4FastStep.cc.

88{
91}
@ fStopAndKill
void ProposePrimaryTrackFinalKineticEnergy(G4double)
void ProposeTrackStatus(G4TrackStatus status)

◆ operator=()

G4FastStep & G4FastStep::operator= ( const G4FastStep & right)
delete

◆ ProposePrimaryTrackFinalEventBiasingWeight()

void G4FastStep::ProposePrimaryTrackFinalEventBiasingWeight ( G4double )

◆ ProposePrimaryTrackFinalKineticEnergy()

void G4FastStep::ProposePrimaryTrackFinalKineticEnergy ( G4double )

◆ ProposePrimaryTrackFinalKineticEnergyAndDirection()

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

Definition at line 138 of file G4FastStep.cc.

141{
142 // Compute global direction if needed...
143 G4ThreeVector globalDirection = direction;
144 if (localCoordinates)
145 globalDirection = fFastTrack->GetInverseAffineTransformation()->TransformAxis(direction);
146 // ...and feed the globalMomentum (ensuring unitarity)
147 SetMomentumChange(globalDirection.unit());
149}
Hep3Vector unit() const

Referenced by SetPrimaryTrackFinalKineticEnergyAndDirection().

◆ ProposePrimaryTrackFinalMomentumDirection()

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

Definition at line 117 of file G4FastStep.cc.

119{
120 // Compute the momentum in global reference
121 // system if needed ...
122 G4ThreeVector globalMomentum = momentum;
123 if (localCoordinates)
124 globalMomentum = fFastTrack->GetInverseAffineTransformation()->TransformAxis(momentum);
125 // ...and feed the globalMomentum (ensuring unitarity)
126 SetMomentumChange(globalMomentum.unit());
127}

Referenced by G4ChannelingFastSimModel::DoIt(), and SetPrimaryTrackFinalMomentum().

◆ ProposePrimaryTrackFinalPolarization()

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

Definition at line 161 of file G4FastStep.cc.

163{
164 // Compute polarization in global system if needed:
165 G4ThreeVector globalPolarization(polarization);
166 if (localCoordinates)
167 globalPolarization =
168 fFastTrack->GetInverseAffineTransformation()->TransformAxis(globalPolarization);
169 // Feed the particle globalPolarization:
170 thePolarizationChange = globalPolarization;
171}

Referenced by SetPrimaryTrackFinalPolarization().

◆ ProposePrimaryTrackFinalPosition()

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

Definition at line 96 of file G4FastStep.cc.

98{
99 // Compute the position coordinate in global
100 // reference system if needed ...
101 G4ThreeVector globalPosition = position;
102 if (localCoordinates)
103 globalPosition = fFastTrack->GetInverseAffineTransformation()->TransformPoint(position);
104 // ...and feed the globalPosition:
105 thePositionChange = globalPosition;
106}

Referenced by G4ChannelingFastSimModel::DoIt(), and SetPrimaryTrackFinalPosition().

◆ ProposePrimaryTrackFinalProperTime()

void G4FastStep::ProposePrimaryTrackFinalProperTime ( G4double )

◆ ProposePrimaryTrackFinalTime()

void G4FastStep::ProposePrimaryTrackFinalTime ( G4double )

◆ ProposePrimaryTrackPathLength()

void G4FastStep::ProposePrimaryTrackPathLength ( G4double )

◆ ProposeTotalEnergyDeposited()

void G4FastStep::ProposeTotalEnergyDeposited ( G4double anEnergyPart)

◆ SetNumberOfSecondaryTracks()

void G4FastStep::SetNumberOfSecondaryTracks ( G4int )

◆ SetPrimaryTrackFinalEventBiasingWeight()

void G4FastStep::SetPrimaryTrackFinalEventBiasingWeight ( G4double )

◆ SetPrimaryTrackFinalKineticEnergy()

void G4FastStep::SetPrimaryTrackFinalKineticEnergy ( G4double )

◆ SetPrimaryTrackFinalKineticEnergyAndDirection()

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

Definition at line 151 of file G4FastStep.cc.

154{
155 ProposePrimaryTrackFinalKineticEnergyAndDirection(kineticEnergy, direction, localCoordinates);
156}
void ProposePrimaryTrackFinalKineticEnergyAndDirection(G4double, const G4ThreeVector &, G4bool localCoordinates=true)

◆ SetPrimaryTrackFinalMomentum()

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

Definition at line 129 of file G4FastStep.cc.

131{
132 ProposePrimaryTrackFinalMomentumDirection(momentum, localCoordinates);
133}
void ProposePrimaryTrackFinalMomentumDirection(const G4ThreeVector &, G4bool localCoordinates=true)

◆ SetPrimaryTrackFinalPolarization()

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

Definition at line 173 of file G4FastStep.cc.

175{
176 ProposePrimaryTrackFinalPolarization(polarization, localCoordinates);
177}
void ProposePrimaryTrackFinalPolarization(const G4ThreeVector &, G4bool localCoordinates=true)

◆ SetPrimaryTrackFinalPosition()

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

Definition at line 108 of file G4FastStep.cc.

110{
112}
void ProposePrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition G4FastStep.cc:96

◆ SetPrimaryTrackFinalProperTime()

void G4FastStep::SetPrimaryTrackFinalProperTime ( G4double )

◆ SetPrimaryTrackFinalTime()

void G4FastStep::SetPrimaryTrackFinalTime ( G4double )

◆ SetPrimaryTrackPathLength()

void G4FastStep::SetPrimaryTrackPathLength ( G4double )

◆ SetTotalEnergyDeposited()

void G4FastStep::SetTotalEnergyDeposited ( G4double anEnergyPart)

◆ UpdateStepForAtRest()

G4Step * G4FastStep::UpdateStepForAtRest ( G4Step * Step)
overridevirtual

Reimplemented from G4VParticleChange.

Definition at line 295 of file G4FastStep.cc.

296{
297 // A physics process always calculates the final state of the particle
298
299 // G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
300 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
301 G4Track* aTrack = pStep->GetTrack();
302 // G4double mass = aTrack->GetDynamicParticle()->GetMass();
303
304 // update kinetic energy and momentum direction
305 pPostStepPoint->SetMomentumDirection(theMomentumChange);
306 pPostStepPoint->SetKineticEnergy(theEnergyChange);
307
308 // update polarization
309 pPostStepPoint->SetPolarization(thePolarizationChange);
310
311 // update position and time
312 pPostStepPoint->SetPosition(thePositionChange);
313 pPostStepPoint->SetGlobalTime(theTimeChange);
314 pPostStepPoint->AddLocalTime(theTimeChange - aTrack->GetGlobalTime());
315 pPostStepPoint->SetProperTime(theProperTimeChange);
316
317 // update weight
318 pPostStepPoint->SetWeight(theWeightChange);
319
320 if (debugFlag) CheckIt(*aTrack);
321
322 // Update the G4Step specific attributes
323 return UpdateStepInfo(pStep);
324}
G4bool CheckIt(const G4Track &) override
void SetKineticEnergy(const G4double aValue)
void SetWeight(G4double aValue)
void SetProperTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetPosition(const G4ThreeVector &aValue)
void AddLocalTime(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4Step * UpdateStepInfo(G4Step *Step)

◆ UpdateStepForPostStep()

G4Step * G4FastStep::UpdateStepForPostStep ( G4Step * Step)
overridevirtual

Reimplemented from G4VParticleChange.

Definition at line 260 of file G4FastStep.cc.

261{
262 // A physics process always calculates the final state of the particle
263
264 // Take note that the return type of GetMomentumChange is a
265 // pointer to G4ParticleMometum. Also it is a normalized
266 // momentum vector.
267
268 // G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
269 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
270 G4Track* aTrack = pStep->GetTrack();
271 // G4double mass = aTrack->GetDynamicParticle()->GetMass();
272
273 // update kinetic energy and momentum direction
274 pPostStepPoint->SetMomentumDirection(theMomentumChange);
275 pPostStepPoint->SetKineticEnergy(theEnergyChange);
276
277 // update polarization
278 pPostStepPoint->SetPolarization(thePolarizationChange);
279
280 // update position and time
281 pPostStepPoint->SetPosition(thePositionChange);
282 pPostStepPoint->SetGlobalTime(theTimeChange);
283 pPostStepPoint->AddLocalTime(theTimeChange - aTrack->GetGlobalTime());
284 pPostStepPoint->SetProperTime(theProperTimeChange);
285
286 // update weight
287 pPostStepPoint->SetWeight(theWeightChange);
288
289 if (debugFlag) CheckIt(*aTrack);
290
291 // Update the G4Step specific attributes
292 return UpdateStepInfo(pStep);
293}

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