Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes
G4FastStep Class Reference

#include <G4FastStep.hh>

Inheritance diagram for G4FastStep:
G4VParticleChange

Public Member Functions

void AddSecondary (G4Track *aSecondary)
 
G4bool CheckIt (const G4Track &)
 
void Clear ()
 
void ClearDebugFlag ()
 
G4TrackCreateSecondaryTrack (const G4DynamicParticle &, G4ThreeVector, G4double, G4bool localCoordinates=true)
 
G4TrackCreateSecondaryTrack (const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)
 
void DumpInfo () const
 
void ForceSteppingHitInvocation ()
 
 G4FastStep ()
 
G4bool GetDebugFlag () const
 
G4bool GetFirstStepInVolume () const
 
G4bool GetLastStepInVolume () const
 
G4double GetLocalEnergyDeposit () const
 
G4double GetNonIonizingEnergyDeposit () const
 
G4int GetNumberOfSecondaries () const
 
G4int GetNumberOfSecondaryTracks ()
 
G4double GetParentWeight () const
 
G4TrackGetSecondary (G4int anIndex) const
 
G4TrackGetSecondaryTrack (G4int)
 
G4SteppingControl GetSteppingControl () const
 
G4double GetTotalEnergyDeposited () const
 
G4TrackStatus GetTrackStatus () const
 
G4double GetTrueStepLength () const
 
G4int GetVerboseLevel () const
 
G4double GetWeight () const
 
void Initialize (const G4FastTrack &)
 
G4bool IsParentWeightSetByProcess () const
 
G4bool IsSecondaryWeightSetByProcess () const
 
void KillPrimaryTrack ()
 
G4bool operator!= (const G4FastStep &right) const
 
G4bool operator!= (const G4VParticleChange &right) const
 
G4bool operator== (const G4FastStep &right) const
 
G4bool operator== (const G4VParticleChange &right) const
 
void ProposeFirstStepInVolume (G4bool flag)
 
void ProposeLastStepInVolume (G4bool flag)
 
void ProposeLocalEnergyDeposit (G4double anEnergyPart)
 
void ProposeNonIonizingEnergyDeposit (G4double anEnergyPart)
 
void ProposeParentWeight (G4double finalWeight)
 
void ProposePrimaryTrackFinalEventBiasingWeight (G4double)
 
void ProposePrimaryTrackFinalKineticEnergy (G4double)
 
void ProposePrimaryTrackFinalKineticEnergyAndDirection (G4double, const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackFinalMomentumDirection (const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackFinalPolarization (const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackFinalPosition (const G4ThreeVector &, G4bool localCoordinates=true)
 
void ProposePrimaryTrackFinalProperTime (G4double)
 
void ProposePrimaryTrackFinalTime (G4double)
 
void ProposePrimaryTrackPathLength (G4double)
 
void ProposeSteppingControl (G4SteppingControl StepControlFlag)
 
void ProposeTotalEnergyDeposited (G4double anEnergyPart)
 
void ProposeTrackStatus (G4TrackStatus status)
 
void ProposeTrueStepLength (G4double truePathLength)
 
void ProposeWeight (G4double finalWeight)
 
void SetDebugFlag ()
 
void SetNumberOfSecondaries (G4int totSecondaries)
 
void SetNumberOfSecondaryTracks (G4int)
 
void SetParentWeightByProcess (G4bool)
 
void SetPrimaryTrackFinalEventBiasingWeight (G4double)
 
void SetPrimaryTrackFinalKineticEnergy (G4double)
 
void SetPrimaryTrackFinalKineticEnergyAndDirection (G4double, const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalMomentum (const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalPolarization (const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalPosition (const G4ThreeVector &, G4bool localCoordinates=true)
 
void SetPrimaryTrackFinalProperTime (G4double)
 
void SetPrimaryTrackFinalTime (G4double)
 
void SetPrimaryTrackPathLength (G4double)
 
void SetSecondaryWeightByProcess (G4bool)
 
void SetTotalEnergyDeposited (G4double anEnergyPart)
 
void SetVerboseLevel (G4int vLevel)
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
G4StepUpdateStepForAtRest (G4Step *Step)
 
G4StepUpdateStepForPostStep (G4Step *Step)
 
virtual ~G4FastStep ()
 

Protected Member Functions

G4bool CheckSecondary (G4Track &)
 
 G4FastStep (const G4FastStep &right)
 
G4double GetAccuracyForException () const
 
G4double GetAccuracyForWarning () const
 
void InitializeLocalEnergyDeposit (const G4Track &)
 
void InitializeParentGlobalTime (const G4Track &)
 
void InitializeParentWeight (const G4Track &)
 
void InitializeSecondaries (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeStepInVolumeFlags (const G4Track &)
 
void InitializeSteppingControl (const G4Track &)
 
void InitializeTrueStepLength (const G4Track &)
 
G4FastStepoperator= (const G4FastStep &right)
 
G4StepUpdateStepInfo (G4Step *Step)
 

Protected Attributes

G4bool debugFlag = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool isParentWeightProposed = false
 
G4bool theFirstStepInVolume = false
 
G4bool theLastStepInVolume = false
 
G4TrackFastVectortheListOfSecondaries = nullptr
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4double theParentGlobalTime = 0.0
 
G4double theParentWeight = 1.0
 
G4int theSizeOftheListOfSecondaries = 0
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theTrueStepLength = 0.0
 
G4int verboseLevel = 1
 

Static Protected Attributes

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

Private Member Functions

void Initialize (const G4Track &)
 
void SetMomentumChange (const G4ThreeVector &Pfinal)
 
void SetMomentumChange (G4double Px, G4double Py, G4double Pz)
 

Private Attributes

const G4FastTrackfFastTrack = nullptr
 
G4double theEnergyChange = 0.0
 
G4ParticleMomentum theMomentumChange
 
G4ThreeVector thePolarizationChange
 
G4ThreeVector thePositionChange
 
G4double theProperTimeChange = 0.0
 
G4double theTimeChange = 0.0
 
G4double theWeightChange = 0.0
 

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

References G4cerr, G4endl, and G4VParticleChange::verboseLevel.

◆ ~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}

References G4cerr, G4endl, and G4VParticleChange::verboseLevel.

◆ 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

◆ AddSecondary()

void G4VParticleChange::AddSecondary ( G4Track aSecondary)
inherited

Definition at line 146 of file G4VParticleChange.cc.

147{
148 if(debugFlag)
149 CheckSecondary(*aTrack);
150
151 // add a secondary after size check
153 {
154 // set weight of secondary tracks
156 aTrack->SetWeight(theParentWeight);
159 }
160 else
161 {
162 delete aTrack;
163
164 G4Exception("G4VParticleChange::AddSecondary()", "TRACK101", JustWarning,
165 "Secondary buffer is full. The track is deleted!");
166 }
167}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
G4TrackFastVector * theListOfSecondaries
G4bool CheckSecondary(G4Track &)

References G4VParticleChange::CheckSecondary(), G4VParticleChange::debugFlag, G4VParticleChange::fSetSecondaryWeightByProcess, G4Exception(), JustWarning, G4FastVector< Type, N >::SetElement(), G4Track::SetWeight(), G4VParticleChange::theListOfSecondaries, G4VParticleChange::theNumberOfSecondaries, G4VParticleChange::theParentWeight, and G4VParticleChange::theSizeOftheListOfSecondaries.

Referenced by G4ParticleChangeForGamma::AddSecondary(), G4ParticleChange::AddSecondary(), G4eplusAnnihilation::AtRestDoIt(), CreateSecondaryTrack(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4VEnergyLossProcess::FillSecondariesAlongStep(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), and G4ParticleChangeForOccurenceBiasing::StealSecondaries().

◆ 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) {
556 }
557
558 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
559 return itsOK;
560}
@ FatalException
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4GLOB_DLL std::ostream G4cout
double mag2() const
double mag() const
G4double theEnergyChange
Definition: G4FastStep.hh:357
void DumpInfo() const
Definition: G4FastStep.cc:437
G4double theProperTimeChange
Definition: G4FastStep.hh:366
G4ParticleMomentum theMomentumChange
Definition: G4FastStep.hh:351
G4double theTimeChange
Definition: G4FastStep.hh:363
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

References G4VParticleChange::CheckIt(), DumpInfo(), FatalException, G4cout, G4endl, G4Exception(), G4VParticleChange::GetAccuracyForException(), G4VParticleChange::GetAccuracyForWarning(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4Track::GetProperTime(), JustWarning, CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), MeV, ns, theEnergyChange, theMomentumChange, theProperTimeChange, and theTimeChange.

Referenced by UpdateStepForAtRest(), and UpdateStepForPostStep().

◆ CheckSecondary()

G4bool G4VParticleChange::CheckSecondary ( G4Track aTrack)
protectedinherited

Definition at line 385 of file G4VParticleChange.cc.

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

References G4VParticleChange::accuracyForException, G4VParticleChange::accuracyForWarning, G4VParticleChange::DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), G4ThreadLocal, G4Track::GetDefinition(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4Track::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), m, MeV, ns, G4Track::SetKineticEnergy(), G4Track::SetMomentumDirection(), G4VParticleChange::theParentGlobalTime, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4VParticleChange::AddSecondary().

◆ Clear()

void G4VParticleChange::Clear ( )
inherited

◆ ClearDebugFlag()

void G4VParticleChange::ClearDebugFlag ( )
inherited

◆ 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}
void SetPolarization(const G4ThreeVector &)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetPolarization() const
const G4FastTrack * fFastTrack
Definition: G4FastStep.hh:369
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:236
void AddSecondary(G4Track *aSecondary)

References G4VParticleChange::AddSecondary(), fFastTrack, G4FastTrack::GetInverseAffineTransformation(), G4DynamicParticle::GetPolarization(), G4DynamicParticle::SetMomentumDirection(), G4DynamicParticle::SetPolarization(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ 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

References CreateSecondaryTrack(), G4DynamicParticle::SetPolarization(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

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)
G4ThreeVector thePolarizationChange
Definition: G4FastStep.hh:354
G4ThreeVector thePositionChange
Definition: G4FastStep.hh:360

References G4VParticleChange::DumpInfo(), G4BestUnit, G4cout, G4endl, theEnergyChange, theMomentumChange, thePolarizationChange, thePositionChange, theProperTimeChange, theTimeChange, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by CheckIt().

◆ ForceSteppingHitInvocation()

void G4FastStep::ForceSteppingHitInvocation ( )

◆ GetAccuracyForException()

G4double G4VParticleChange::GetAccuracyForException ( ) const
protectedinherited

Definition at line 507 of file G4VParticleChange.cc.

508{
510}

References G4VParticleChange::accuracyForException.

Referenced by CheckIt().

◆ GetAccuracyForWarning()

G4double G4VParticleChange::GetAccuracyForWarning ( ) const
protectedinherited

Definition at line 501 of file G4VParticleChange.cc.

502{
503 return accuracyForWarning;
504}

References G4VParticleChange::accuracyForWarning.

Referenced by CheckIt().

◆ GetDebugFlag()

G4bool G4VParticleChange::GetDebugFlag ( ) const
inherited

◆ GetFirstStepInVolume()

G4bool G4VParticleChange::GetFirstStepInVolume ( ) const
inherited

◆ GetLastStepInVolume()

G4bool G4VParticleChange::GetLastStepInVolume ( ) const
inherited

◆ GetLocalEnergyDeposit()

G4double G4VParticleChange::GetLocalEnergyDeposit ( ) const
inherited

◆ GetNonIonizingEnergyDeposit()

G4double G4VParticleChange::GetNonIonizingEnergyDeposit ( ) const
inherited

◆ GetNumberOfSecondaries()

G4int G4VParticleChange::GetNumberOfSecondaries ( ) const
inherited

◆ GetNumberOfSecondaryTracks()

G4int G4FastStep::GetNumberOfSecondaryTracks ( )

◆ GetParentWeight()

G4double G4VParticleChange::GetParentWeight ( ) const
inherited

◆ GetSecondary()

G4Track * G4VParticleChange::GetSecondary ( G4int  anIndex) const
inherited

◆ GetSecondaryTrack()

G4Track * G4FastStep::GetSecondaryTrack ( G4int  )

◆ GetSteppingControl()

G4SteppingControl G4VParticleChange::GetSteppingControl ( ) const
inherited

◆ GetTotalEnergyDeposited()

G4double G4FastStep::GetTotalEnergyDeposited ( ) const

◆ GetTrackStatus()

G4TrackStatus G4VParticleChange::GetTrackStatus ( ) const
inherited

◆ GetTrueStepLength()

G4double G4VParticleChange::GetTrueStepLength ( ) const
inherited

◆ GetVerboseLevel()

G4int G4VParticleChange::GetVerboseLevel ( ) const
inherited

◆ GetWeight()

G4double G4VParticleChange::GetWeight ( ) const
inherited

◆ Initialize() [1/2]

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();
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
G4double theWeightChange
Definition: G4FastStep.hh:372
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:206
G4double GetWeight() const
const G4DynamicParticle * GetDynamicParticle() const
G4SteppingControl theSteppingControlFlag
virtual void Initialize(const G4Track &)

References AvoidHitInvocation, fFastTrack, G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4DynamicParticle::GetPolarization(), G4Track::GetPosition(), G4FastTrack::GetPrimaryTrack(), G4DynamicParticle::GetProperTime(), G4Track::GetWeight(), G4VParticleChange::Initialize(), theEnergyChange, theMomentumChange, thePolarizationChange, thePositionChange, theProperTimeChange, G4VParticleChange::theSteppingControlFlag, theTimeChange, and theWeightChange.

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

◆ Initialize() [2/2]

void G4FastStep::Initialize ( const G4Track )
privatevirtual

Reimplemented from G4VParticleChange.

Definition at line 286 of file G4FastStep.cc.

287{
288 G4ExceptionDescription tellWhatIsWrong;
289 tellWhatIsWrong << "G4FastStep can be initialised only through G4FastTrack."
290 << G4endl;
291 G4Exception("G4FastStep::Initialize(const G4Track&)",
292 "FastSim005",
294 tellWhatIsWrong);
295}

References FatalException, G4endl, and G4Exception().

◆ InitializeLocalEnergyDeposit()

void G4VParticleChange::InitializeLocalEnergyDeposit ( const G4Track )
protectedinherited

◆ InitializeParentGlobalTime()

void G4VParticleChange::InitializeParentGlobalTime ( const G4Track )
protectedinherited

◆ InitializeParentWeight()

void G4VParticleChange::InitializeParentWeight ( const G4Track )
protectedinherited

◆ InitializeSecondaries()

void G4VParticleChange::InitializeSecondaries ( const G4Track )
protectedinherited

◆ InitializeStatusChange()

void G4VParticleChange::InitializeStatusChange ( const G4Track )
protectedinherited

◆ InitializeStepInVolumeFlags()

void G4VParticleChange::InitializeStepInVolumeFlags ( const G4Track )
protectedinherited

◆ InitializeSteppingControl()

void G4VParticleChange::InitializeSteppingControl ( const G4Track )
protectedinherited

◆ InitializeTrueStepLength()

void G4VParticleChange::InitializeTrueStepLength ( const G4Track )
protectedinherited

◆ IsParentWeightSetByProcess()

G4bool G4VParticleChange::IsParentWeightSetByProcess ( ) const
inherited

Definition at line 516 of file G4VParticleChange.cc.

516{ return true; }

◆ IsSecondaryWeightSetByProcess()

G4bool G4VParticleChange::IsSecondaryWeightSetByProcess ( ) const
inherited

◆ KillPrimaryTrack()

void G4FastStep::KillPrimaryTrack ( )

Definition at line 87 of file G4FastStep.cc.

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

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

Referenced by GFlashShowerModel::ElectronDoIt().

◆ operator!=() [1/2]

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

Definition at line 355 of file G4FastStep.cc.

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

◆ operator!=() [2/2]

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

Definition at line 140 of file G4VParticleChange.cc.

141{
142 return (this != (G4VParticleChange*) &right);
143}

◆ operator=()

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

Definition at line 322 of file G4FastStep.cc.

References fFastTrack, G4VParticleChange::operator=(), theEnergyChange, G4VParticleChange::theListOfSecondaries, G4VParticleChange::theLocalEnergyDeposit, theMomentumChange, G4VParticleChange::theNumberOfSecondaries, thePolarizationChange, thePositionChange, theProperTimeChange, G4VParticleChange::theSizeOftheListOfSecondaries, G4VParticleChange::theStatusChange, G4VParticleChange::theSteppingControlFlag, theTimeChange, G4VParticleChange::theTrueStepLength, and theWeightChange.

◆ operator==() [1/2]

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

Definition at line 350 of file G4FastStep.cc.

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

◆ operator==() [2/2]

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

Definition at line 134 of file G4VParticleChange.cc.

135{
136 return (this == (G4VParticleChange*) &right);
137}

◆ ProposeFirstStepInVolume()

void G4VParticleChange::ProposeFirstStepInVolume ( G4bool  flag)
inherited

◆ ProposeLastStepInVolume()

void G4VParticleChange::ProposeLastStepInVolume ( G4bool  flag)
inherited

◆ ProposeLocalEnergyDeposit()

void G4VParticleChange::ProposeLocalEnergyDeposit ( G4double  anEnergyPart)
inherited

Referenced by G4VEnergyLossProcess::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4OpBoundaryProcess::DoAbsorption(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4SpecialCuts::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4DNABornExcitationModel1::SampleSecondaries(), G4DNABornExcitationModel2::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNADiracRMatrixExcitationModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNAQuinnPlasmonExcitationModel::SampleSecondaries(), G4DNARelativisticIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4PAIPhotModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreIonisationModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4PenelopeBremsstrahlungModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopeIonisationModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4LEPTSAttachmentModel::SampleSecondaries(), G4LEPTSDissociationModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSExcitationModel::SampleSecondaries(), G4LEPTSIonisationModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4LEPTSRotExcitationModel::SampleSecondaries(), G4LEPTSVibExcitationModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ ProposeNonIonizingEnergyDeposit()

void G4VParticleChange::ProposeNonIonizingEnergyDeposit ( G4double  anEnergyPart)
inherited

◆ ProposeParentWeight()

void G4VParticleChange::ProposeParentWeight ( G4double  finalWeight)
inherited

◆ 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
void SetMomentumChange(G4double Px, G4double Py, G4double Pz)

References fFastTrack, G4FastTrack::GetInverseAffineTransformation(), SetMomentumChange(), SetPrimaryTrackFinalKineticEnergy(), and CLHEP::Hep3Vector::unit().

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}

References fFastTrack, G4FastTrack::GetInverseAffineTransformation(), SetMomentumChange(), and CLHEP::Hep3Vector::unit().

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}

References fFastTrack, G4FastTrack::GetInverseAffineTransformation(), and thePolarizationChange.

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

References fFastTrack, G4FastTrack::GetInverseAffineTransformation(), position, and thePositionChange.

Referenced by SetPrimaryTrackFinalPosition().

◆ ProposePrimaryTrackFinalProperTime()

void G4FastStep::ProposePrimaryTrackFinalProperTime ( G4double  )

◆ ProposePrimaryTrackFinalTime()

void G4FastStep::ProposePrimaryTrackFinalTime ( G4double  )

◆ ProposePrimaryTrackPathLength()

void G4FastStep::ProposePrimaryTrackPathLength ( G4double  )

◆ ProposeSteppingControl()

void G4VParticleChange::ProposeSteppingControl ( G4SteppingControl  StepControlFlag)
inherited

◆ ProposeTotalEnergyDeposited()

void G4FastStep::ProposeTotalEnergyDeposited ( G4double  anEnergyPart)

◆ ProposeTrackStatus()

void G4VParticleChange::ProposeTrackStatus ( G4TrackStatus  status)
inherited

Referenced by G4BiasingProcessInterface::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4OpBoundaryProcess::DoAbsorption(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), KillPrimaryTrack(), G4ImportanceProcess::KillTrack(), G4WeightWindowProcess::KillTrack(), G4DNAElectronHoleRecombination::MakeReaction(), G4SpecialCuts::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4FastSimulationManagerProcess::PostStepDoIt(), G4PhononDownconversion::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNLoss::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4BiasingProcessInterface::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4PenelopeAnnihilationModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4eeToTwoGammaModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4eplusTo2GammaOKVIModel::SampleSecondaries(), G4eplusTo3GammaOKVIModel::SampleSecondaries(), G4eeToHadronsMultiModel::SampleSecondaries(), G4LEPTSAttachmentModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), and G4BetheHeitler5DModel::SampleSecondaries().

◆ ProposeTrueStepLength()

void G4VParticleChange::ProposeTrueStepLength ( G4double  truePathLength)
inherited

◆ ProposeWeight()

void G4VParticleChange::ProposeWeight ( G4double  finalWeight)
inherited

◆ SetDebugFlag()

void G4VParticleChange::SetDebugFlag ( )
inherited

◆ SetMomentumChange() [1/2]

void G4FastStep::SetMomentumChange ( const G4ThreeVector Pfinal)
private

◆ SetMomentumChange() [2/2]

void G4FastStep::SetMomentumChange ( G4double  Px,
G4double  Py,
G4double  Pz 
)
private

◆ SetNumberOfSecondaries()

void G4VParticleChange::SetNumberOfSecondaries ( G4int  totSecondaries)
inherited

◆ SetNumberOfSecondaryTracks()

void G4FastStep::SetNumberOfSecondaryTracks ( G4int  )

◆ SetParentWeightByProcess()

void G4VParticleChange::SetParentWeightByProcess ( G4bool  )
inherited

◆ 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

References ProposePrimaryTrackFinalKineticEnergyAndDirection().

◆ 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

References ProposePrimaryTrackFinalMomentumDirection().

◆ 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

References ProposePrimaryTrackFinalPolarization().

◆ 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

References ProposePrimaryTrackFinalPosition().

◆ SetPrimaryTrackFinalProperTime()

void G4FastStep::SetPrimaryTrackFinalProperTime ( G4double  )

◆ SetPrimaryTrackFinalTime()

void G4FastStep::SetPrimaryTrackFinalTime ( G4double  )

◆ SetPrimaryTrackPathLength()

void G4FastStep::SetPrimaryTrackPathLength ( G4double  )

◆ SetSecondaryWeightByProcess()

void G4VParticleChange::SetSecondaryWeightByProcess ( G4bool  )
inherited

◆ SetTotalEnergyDeposited()

void G4FastStep::SetTotalEnergyDeposited ( G4double  anEnergyPart)

◆ SetVerboseLevel()

void G4VParticleChange::SetVerboseLevel ( G4int  vLevel)
inherited

◆ UpdateStepForAlongStep()

G4Step * G4VParticleChange::UpdateStepForAlongStep ( G4Step Step)
virtualinherited

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

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(), theEnergyChange, theMomentumChange, thePolarizationChange, thePositionChange, theProperTimeChange, theTimeChange, theWeightChange, and G4VParticleChange::UpdateStepInfo().

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

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(), theEnergyChange, theMomentumChange, thePolarizationChange, thePositionChange, theProperTimeChange, theTimeChange, theWeightChange, and G4VParticleChange::UpdateStepInfo().

◆ UpdateStepInfo()

G4Step * G4VParticleChange::UpdateStepInfo ( G4Step Step)
protectedinherited

Definition at line 170 of file G4VParticleChange.cc.

171{
172 // Update the G4Step specific attributes
173 pStep->SetStepLength(theTrueStepLength);
174 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
175 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
176 pStep->SetControlFlag(theSteppingControlFlag);
177
179 {
180 pStep->SetFirstStepFlag();
181 }
182 else
183 {
184 pStep->ClearFirstStepFlag();
185 }
187 {
188 pStep->SetLastStepFlag();
189 }
190 else
191 {
192 pStep->ClearLastStepFlag();
193 }
194
195 return pStep;
196}
G4double theNonIonizingEnergyDeposit

References G4Step::AddNonIonizingEnergyDeposit(), G4Step::AddTotalEnergyDeposit(), G4Step::ClearFirstStepFlag(), G4Step::ClearLastStepFlag(), G4Step::SetControlFlag(), G4Step::SetFirstStepFlag(), G4Step::SetLastStepFlag(), G4Step::SetStepLength(), G4VParticleChange::theFirstStepInVolume, G4VParticleChange::theLastStepInVolume, G4VParticleChange::theLocalEnergyDeposit, G4VParticleChange::theNonIonizingEnergyDeposit, G4VParticleChange::theSteppingControlFlag, and G4VParticleChange::theTrueStepLength.

Referenced by G4VParticleChange::UpdateStepForAlongStep(), UpdateStepForAtRest(), G4ParticleChangeForDecay::UpdateStepForAtRest(), G4VParticleChange::UpdateStepForAtRest(), UpdateStepForPostStep(), G4ParticleChangeForDecay::UpdateStepForPostStep(), and G4VParticleChange::UpdateStepForPostStep().

Field Documentation

◆ accuracyForException

const G4double G4VParticleChange::accuracyForException = 0.001
staticprotectedinherited

◆ accuracyForWarning

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

◆ debugFlag

G4bool G4VParticleChange::debugFlag = false
protectedinherited

◆ fFastTrack

const G4FastTrack* G4FastStep::fFastTrack = nullptr
private

◆ fSetSecondaryWeightByProcess

G4bool G4VParticleChange::fSetSecondaryWeightByProcess = false
protectedinherited

◆ isParentWeightProposed

G4bool G4VParticleChange::isParentWeightProposed = false
protectedinherited

◆ theEnergyChange

G4double G4FastStep::theEnergyChange = 0.0
private

◆ theFirstStepInVolume

G4bool G4VParticleChange::theFirstStepInVolume = false
protectedinherited

◆ theLastStepInVolume

G4bool G4VParticleChange::theLastStepInVolume = false
protectedinherited

◆ theListOfSecondaries

G4TrackFastVector* G4VParticleChange::theListOfSecondaries = nullptr
protectedinherited

◆ theLocalEnergyDeposit

G4double G4VParticleChange::theLocalEnergyDeposit = 0.0
protectedinherited

◆ theMomentumChange

G4ParticleMomentum G4FastStep::theMomentumChange
private

◆ theNonIonizingEnergyDeposit

G4double G4VParticleChange::theNonIonizingEnergyDeposit = 0.0
protectedinherited

◆ theNumberOfSecondaries

G4int G4VParticleChange::theNumberOfSecondaries = 0
protectedinherited

◆ theParentGlobalTime

G4double G4VParticleChange::theParentGlobalTime = 0.0
protectedinherited

◆ theParentWeight

G4double G4VParticleChange::theParentWeight = 1.0
protectedinherited

◆ thePolarizationChange

G4ThreeVector G4FastStep::thePolarizationChange
private

◆ thePositionChange

G4ThreeVector G4FastStep::thePositionChange
private

◆ theProperTimeChange

G4double G4FastStep::theProperTimeChange = 0.0
private

◆ theSizeOftheListOfSecondaries

G4int G4VParticleChange::theSizeOftheListOfSecondaries = 0
protectedinherited

◆ theStatusChange

G4TrackStatus G4VParticleChange::theStatusChange = fAlive
protectedinherited

◆ theSteppingControlFlag

G4SteppingControl G4VParticleChange::theSteppingControlFlag = NormalCondition
protectedinherited

◆ theTimeChange

G4double G4FastStep::theTimeChange = 0.0
private

◆ theTrueStepLength

G4double G4VParticleChange::theTrueStepLength = 0.0
protectedinherited

◆ theWeightChange

G4double G4FastStep::theWeightChange = 0.0
private

Definition at line 372 of file G4FastStep.hh.

Referenced by Initialize(), operator=(), UpdateStepForAtRest(), and UpdateStepForPostStep().

◆ verboseLevel

G4int G4VParticleChange::verboseLevel = 1
protectedinherited

Definition at line 274 of file G4VParticleChange.hh.

Referenced by G4FastStep(), G4VParticleChange::operator=(), and ~G4FastStep().


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