Geant4 10.7.0
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

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

 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 = nullptr
 
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
 
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
 

Detailed Description

Definition at line 90 of file G4FastStep.hh.

Constructor & Destructor Documentation

◆ G4FastStep() [1/2]

G4FastStep::G4FastStep ( )

Definition at line 297 of file G4FastStep.cc.

299{
300 if (verboseLevel>2)
301 {
302 G4cerr << "G4FastStep::G4FastStep()" << G4endl;
303 }
304}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57

◆ ~G4FastStep()

G4FastStep::~G4FastStep ( )
virtual

Definition at line 306 of file G4FastStep.cc.

307{
308 if (verboseLevel>2)
309 {
310 G4cerr << "G4FastStep::~G4FastStep()" << G4endl;
311 }
312}

◆ G4FastStep() [2/2]

G4FastStep::G4FastStep ( const G4FastStep right)
protected

Definition at line 315 of file G4FastStep.cc.

317{
318 *this = right;
319}

Member Function Documentation

◆ CheckIt()

G4bool G4FastStep::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 460 of file G4FastStep.cc.

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

Referenced by UpdateStepForAtRest(), and UpdateStepForPostStep().

◆ CreateSecondaryTrack() [1/2]

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

Definition at line 223 of file G4FastStep.cc.

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

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

Referenced by CreateSecondaryTrack().

◆ DumpInfo()

void G4FastStep::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 437 of file G4FastStep.cc.

438{
439// use base-class DumpInfo
441
442 G4cout << " Position - x (mm) : " << G4BestUnit( thePositionChange.x(), "Length" ) << G4endl;
443 G4cout << " Position - y (mm) : " << G4BestUnit( thePositionChange.y(), "Length" ) << G4endl;
444 G4cout << " Position - z (mm) : " << G4BestUnit( thePositionChange.z(), "Length" ) << G4endl;
445 G4cout << " Time (ns) : " << G4BestUnit( theTimeChange, "Time" ) << G4endl;
446 G4cout << " Proper Time (ns) : " << G4BestUnit( theProperTimeChange, "Time" ) << G4endl;
447 G4int olprc = G4cout.precision(3);
448 G4cout << " Momentum Direct - x : " << std::setw(20) << theMomentumChange.x() << G4endl;
449 G4cout << " Momentum Direct - y : " << std::setw(20) << theMomentumChange.y() << G4endl;
450 G4cout << " Momentum Direct - z : " << std::setw(20) << theMomentumChange.z() << G4endl;
451 G4cout.precision(olprc);
452 G4cout << " Kinetic Energy (MeV): " << G4BestUnit( theEnergyChange, "Energy" ) << G4endl;
453 G4cout.precision(3);
454 G4cout << " Polarization - x : " << std::setw(20) << thePolarizationChange.x() << G4endl;
455 G4cout << " Polarization - y : " << std::setw(20) << thePolarizationChange.y() << G4endl;
456 G4cout << " Polarization - z : " << std::setw(20) << thePolarizationChange.z() << G4endl;
457 G4cout.precision(olprc);
458}
#define G4BestUnit(a, b)
int G4int
Definition: G4Types.hh:85
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 G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:206
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 SetPrimaryTrackFinalKineticEnergy(G4double)
void ProposeTrackStatus(G4TrackStatus status)

◆ operator!=()

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

Definition at line 355 of file G4FastStep.cc.

356{
357 return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
358}

◆ operator=()

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

Definition at line 322 of file G4FastStep.cc.

323{
324 if (this != &right)
325 {
331 theMomentumChange = right.theMomentumChange;
332 thePolarizationChange = right.thePolarizationChange;
333 thePositionChange = right.thePositionChange;
334 theProperTimeChange = right.theProperTimeChange;
335 theTimeChange = right.theTimeChange;
336 theEnergyChange = right.theEnergyChange;
340 theWeightChange = right.theWeightChange;
341 fFastTrack = right.fFastTrack;
342 }
343 return *this;
344}
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4VParticleChange & operator=(const G4VParticleChange &right)

◆ operator==()

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

Definition at line 350 of file G4FastStep.cc.

351{
352 return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
353}

◆ 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 149 of file G4FastStep.cc.

153{
154 // Compute global direction if needed...
155 G4ThreeVector globalDirection = direction;
156 if (localCoordinates)
157 globalDirection =fFastTrack->GetInverseAffineTransformation()->
158 TransformAxis(direction);
159 // ...and feed the globalMomentum (ensuring unitarity)
160 SetMomentumChange(globalDirection.unit());
162}
Hep3Vector unit() const

Referenced by SetPrimaryTrackFinalKineticEnergyAndDirection().

◆ ProposePrimaryTrackFinalMomentumDirection()

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

Definition at line 123 of file G4FastStep.cc.

126{
127 // Compute the momentum in global reference
128 // system if needed ...
129 G4ThreeVector globalMomentum = momentum;
130 if (localCoordinates)
131 globalMomentum = fFastTrack->GetInverseAffineTransformation()->
132 TransformAxis(momentum);
133 // ...and feed the globalMomentum (ensuring unitarity)
134 SetMomentumChange(globalMomentum.unit());
135}

Referenced by SetPrimaryTrackFinalMomentum().

◆ ProposePrimaryTrackFinalPolarization()

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

Definition at line 177 of file G4FastStep.cc.

180{
181 // Compute polarization in global system if needed:
182 G4ThreeVector globalPolarization(polarization);
183 if (localCoordinates)
184 globalPolarization = fFastTrack->GetInverseAffineTransformation()->
185 TransformAxis(globalPolarization);
186 // Feed the particle globalPolarization:
187 thePolarizationChange = globalPolarization;
188}

Referenced by SetPrimaryTrackFinalPolarization().

◆ ProposePrimaryTrackFinalPosition()

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

Definition at line 97 of file G4FastStep.cc.

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

Referenced by 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 165 of file G4FastStep.cc.

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

◆ SetPrimaryTrackFinalMomentum()

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

Definition at line 138 of file G4FastStep.cc.

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

◆ SetPrimaryTrackFinalPolarization()

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

Definition at line 191 of file G4FastStep.cc.

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

◆ SetPrimaryTrackFinalPosition()

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

Definition at line 112 of file G4FastStep.cc.

115{
117}
void ProposePrimaryTrackFinalPosition(const G4ThreeVector &, G4bool localCoordinates=true)
Definition: G4FastStep.cc:98

◆ 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)
virtual

Reimplemented from G4VParticleChange.

Definition at line 401 of file G4FastStep.cc.

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

Reimplemented from G4VParticleChange.

Definition at line 364 of file G4FastStep.cc.

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

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