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

#include <G4Transportation.hh>

+ Inheritance diagram for G4Transportation:

Public Member Functions

 G4Transportation (G4int verbosityLevel=1, const G4String &aName="Transportation")
 
 ~G4Transportation ()
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
 
G4bool FieldExertedForce ()
 
G4PropagatorInFieldGetPropagatorInField ()
 
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
 
G4double GetThresholdWarningEnergy () const
 
G4double GetThresholdImportantEnergy () const
 
G4int GetThresholdTrials () const
 
void SetThresholdWarningEnergy (G4double newEnWarn)
 
void SetThresholdImportantEnergy (G4double newEnImp)
 
void SetThresholdTrials (G4int newMaxTrials)
 
void SetHighLooperThresholds ()
 
void SetLowLooperThresholds ()
 
void PushThresholdsToLogger ()
 
void ReportLooperThresholds ()
 
G4double GetMaxEnergyKilled () const
 
G4double GetSumEnergyKilled () const
 
void ResetKilledStatistics (G4int report=1)
 
void EnableShortStepOptimisation (G4bool optimise=true)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
void StartTracking (G4Track *aTrack)
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void PrintStatistics (std::ostream &outStr) const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static G4bool EnableMagneticMoment (G4bool useMoment=true)
 
static G4bool EnableGravity (G4bool useGravity=true)
 
static void SetSilenceLooperWarnings (G4bool val)
 
static G4bool GetSilenceLooperWarnings ()
 
static G4bool EnableUseMagneticMoment (G4bool useMoment=true)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

void SetTouchableInformation (const G4TouchableHandle &touchable)
 
void ReportMissingLogger (const char *methodName)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4NavigatorfLinearNavigator
 
G4PropagatorInFieldfFieldPropagator
 
G4ThreeVector fTransportEndPosition = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4ThreeVector fTransportEndMomentumDir = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4double fTransportEndKineticEnergy = 0.0
 
G4ThreeVector fTransportEndSpin = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4bool fMomentumChanged = true
 
G4bool fEndGlobalTimeComputed = false
 
G4double fCandidateEndGlobalTime = 0.0
 
G4bool fAnyFieldExists = false
 
G4bool fParticleIsLooping = false
 
G4bool fNewTrack = true
 
G4bool fFirstStepInVolume = true
 
G4bool fLastStepInVolume = false
 
G4bool fGeometryLimitedStep = true
 
G4bool fFieldExertedForce = false
 
G4TouchableHandle fCurrentTouchableHandle
 
G4ThreeVector fPreviousSftOrigin
 
G4double fPreviousSafety
 
G4ParticleChangeForTransport fParticleChange
 
G4double fEndPointDistance
 
G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV
 
G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV
 
G4int fThresholdTrials = 10
 
G4int fAbandonUnstableTrials = 0
 
G4int fNoLooperTrials = 0
 
G4double fSumEnergyKilled = 0.0
 
G4double fSumEnerSqKilled = 0.0
 
G4double fMaxEnergyKilled = -1.0
 
G4int fMaxEnergyKilledPDG = 0
 
unsigned long fNumLoopersKilled = 0
 
G4double fSumEnergyKilled_NonElectron = 0.0
 
G4double fSumEnerSqKilled_NonElectron = 0.0
 
G4double fMaxEnergyKilled_NonElectron = -1.0
 
G4int fMaxEnergyKilled_NonElecPDG = 0
 
unsigned long fNumLoopersKilled_NonElectron = 0
 
G4double fSumEnergySaved = 0.0
 
G4double fMaxEnergySaved = -1.0
 
G4double fSumEnergyUnstableSaved = 0.0
 
G4bool fShortStepOptimisation
 
G4SafetyHelperfpSafetyHelper
 
G4TransportationLoggerfpLogger
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Static Protected Attributes

static G4bool fUseMagneticMoment =false
 
static G4bool fUseGravity = false
 
static G4bool fSilenceLooperWarnings = false
 

Detailed Description

Definition at line 57 of file G4Transportation.hh.

Constructor & Destructor Documentation

◆ G4Transportation()

G4Transportation::G4Transportation ( G4int verbosityLevel = 1,
const G4String & aName = "Transportation" )

Definition at line 74 of file G4Transportation.cc.

75 : G4VProcess( aName, fTransportation ),
76 fFieldExertedForce( false ),
77 fPreviousSftOrigin( 0.,0.,0. ),
78 fPreviousSafety( 0.0 ),
79 fEndPointDistance( -1.0 ),
80 fShortStepOptimisation( false ) // Old default: true (=fast short steps)
81{
83 pParticleChange= &fParticleChange; // Required to conform to G4VProcess
84 SetVerboseLevel(verbosity);
85
86 G4TransportationManager* transportMgr ;
87
89
90 fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
91
92 fFieldPropagator = transportMgr->GetPropagatorInField() ;
93
94 fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
95
96 fpLogger = new G4TransportationLogger("G4Transportation", verbosity);
97
99 {
101
102 SetThresholdWarningEnergy( trParams->GetWarningEnergy() );
103 SetThresholdImportantEnergy( trParams->GetImportantEnergy() );
104 SetThresholdTrials( trParams->GetNumberOfTrials() );
105 G4Transportation::fSilenceLooperWarnings= trParams->GetSilenceAllLooperWarnings();
106 }
107 else {
109 // Use the old defaults: Warning = 100 MeV, Important = 250 MeV, No Trials = 10;
110 }
111
113 // Should be done by Set methods in SetHighLooperThresholds -- making sure
114
115 static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
116 if ( !pNullTouchableHandle)
117 {
118 pNullTouchableHandle = new G4TouchableHandle;
119 }
120 fCurrentTouchableHandle = *pNullTouchableHandle;
121 // Points to (G4VTouchable*) 0
122
123
124#ifdef G4VERBOSE
125 if( verboseLevel > 0)
126 {
127 G4cout << " G4Transportation constructor> set fShortStepOptimisation to ";
128 if ( fShortStepOptimisation ) { G4cout << "true" << G4endl; }
129 else { G4cout << "false" << G4endl; }
130 }
131#endif
132}
@ fTransportation
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4SafetyHelper * GetSafetyHelper() const
G4Navigator * GetNavigatorForTracking() const
static G4TransportationParameters * Instance()
void SetThresholdImportantEnergy(G4double newEnImp)
G4ThreeVector fPreviousSftOrigin
void PushThresholdsToLogger()
static G4bool fSilenceLooperWarnings
G4PropagatorInField * fFieldPropagator
void SetThresholdTrials(G4int newMaxTrials)
G4ParticleChangeForTransport fParticleChange
G4TouchableHandle fCurrentTouchableHandle
G4Navigator * fLinearNavigator
G4TransportationLogger * fpLogger
void SetThresholdWarningEnergy(G4double newEnWarn)
G4SafetyHelper * fpSafetyHelper
void SetVerboseLevel(G4int value)
G4int verboseLevel
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
#define G4ThreadLocal
Definition tls.hh:77

◆ ~G4Transportation()

G4Transportation::~G4Transportation ( )

Definition at line 136 of file G4Transportation.cc.

137{
138 if( fSumEnergyKilled > 0.0 )
139 {
141 }
142 delete fpLogger;
143}
void PrintStatistics(std::ostream &outStr) const

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4Transportation::AlongStepDoIt ( const G4Track & track,
const G4Step & stepData )
virtual

Implements G4VProcess.

Definition at line 496 of file G4Transportation.cc.

498{
499#if defined(G4VERBOSE) || defined(G4DEBUG_TRANSPORT)
501 noCallsASDI++;
502#else
503 #define noCallsASDI 0
504#endif
505
507 {
509 }
510
512
513 // Code for specific process
514 //
519
521
522 G4double deltaTime = 0.0 ;
523
524 // Calculate Lab Time of Flight (ONLY if field Equations used it!)
525 // G4double endTime = fCandidateEndGlobalTime;
526 // G4double delta_time = endTime - startTime;
527
528 G4double startTime = track.GetGlobalTime() ;
529
531 {
532 // The time was not integrated .. make the best estimate possible
533 //
534 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
535 G4double stepLength = track.GetStepLength();
536
537 deltaTime= 0.0; // in case initialVelocity = 0
538 if ( initialVelocity > 0.0 ) { deltaTime = stepLength/initialVelocity; }
539
540 fCandidateEndGlobalTime = startTime + deltaTime ;
541 fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
542 }
543 else
544 {
545 deltaTime = fCandidateEndGlobalTime - startTime ;
547 }
548
549
550 // Now Correct by Lorentz factor to get delta "proper" Time
551
552 G4double restMass = track.GetDynamicParticle()->GetMass() ;
553 G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
554
555 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
556 //fParticleChange.ProposeTrueStepLength( track.GetStepLength() ) ;
557
558 // If the particle is caught looping or is stuck (in very difficult
559 // boundaries) in a magnetic field (doing many steps) THEN can kill it ...
560 //
561 if ( fParticleIsLooping )
562 {
564 fNoLooperTrials ++;
565 auto particleType= track.GetDynamicParticle()->GetParticleDefinition();
566
567 G4bool stable = particleType->GetPDGStable();
568 G4bool candidateForEnd = (endEnergy < fThreshold_Important_Energy)
570 G4bool unstableAndKillable = !stable && ( fAbandonUnstableTrials != 0);
571 G4bool unstableForEnd = (endEnergy < fThreshold_Important_Energy)
573 if( (candidateForEnd && stable) || (unstableAndKillable && unstableForEnd) )
574 {
575 // Kill the looping particle
576 //
578 G4int particlePDG= particleType->GetPDGEncoding();
579 const G4int electronPDG= 11; // G4Electron::G4Electron()->GetPDGEncoding();
580
581 // Simple statistics
582 fSumEnergyKilled += endEnergy;
583 fSumEnerSqKilled += endEnergy * endEnergy;
585
586 if( endEnergy > fMaxEnergyKilled ) {
587 fMaxEnergyKilled = endEnergy;
588 fMaxEnergyKilledPDG = particlePDG;
589 }
590 if( particleType->GetPDGEncoding() != electronPDG )
591 {
592 fSumEnergyKilled_NonElectron += endEnergy;
593 fSumEnerSqKilled_NonElectron += endEnergy * endEnergy;
595
596 if( endEnergy > fMaxEnergyKilled_NonElectron )
597 {
599 fMaxEnergyKilled_NonElecPDG = particlePDG;
600 }
601 }
602
604 {
606 noCallsASDI, __func__ );
607 }
609 }
610 else
611 {
612 fMaxEnergySaved = std::max( endEnergy, fMaxEnergySaved);
613 if( fNoLooperTrials == 1 ) {
614 fSumEnergySaved += endEnergy;
615 if ( !stable )
616 fSumEnergyUnstableSaved += endEnergy;
617 }
618#ifdef G4VERBOSE
620 {
621 G4cout << " " << __func__
622 << " Particle is looping but is saved ..." << G4endl
623 << " Number of trials = " << fNoLooperTrials << G4endl
624 << " No of calls to = " << noCallsASDI << G4endl;
625 }
626#endif
627 }
628 }
629 else
630 {
632 }
633
634 // Another (sometimes better way) is to use a user-limit maximum Step size
635 // to alleviate this problem ..
636
637 // Introduce smooth curved trajectories to particle-change
638 //
641
642 return &fParticleChange ;
643}
@ fGeomBoundary
@ fStopAndKill
#define noCallsASDI
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
G4double GetMass() const
const G4ParticleDefinition * GetParticleDefinition() const
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *vec)
void Initialize(const G4Track &) final
void SetMomentumChanged(G4bool b)
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void ProposePosition(G4double x, G4double y, G4double z)
void ProposeLocalTime(G4double t)
void ProposeProperTime(G4double finalProperTime)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeGlobalTime(G4double t)
G4bool GetPDGStable() const
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4double GetVelocity() const
void SetStepStatus(const G4StepStatus aValue)
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetStepLength() const
G4double GetTotalEnergy() const
void ReportLoopingTrack(const G4Track &track, const G4Step &stepInfo, G4int numTrials, G4long noCalls, const char *methodName) const
G4double fThreshold_Important_Energy
G4double fMaxEnergyKilled_NonElectron
G4ThreeVector fTransportEndSpin
G4double fTransportEndKineticEnergy
unsigned long fNumLoopersKilled_NonElectron
G4ThreeVector fTransportEndPosition
G4double fThreshold_Warning_Energy
G4double fSumEnerSqKilled_NonElectron
G4double fSumEnergyUnstableSaved
G4double fSumEnergyKilled_NonElectron
G4ThreeVector fTransportEndMomentumDir
unsigned long fNumLoopersKilled
G4double fCandidateEndGlobalTime
void ProposeTrackStatus(G4TrackStatus status)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4Transportation::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety,
G4GPILSelection * selection )
virtual

Implements G4VProcess.

Reimplemented in G4TransportationWithMsc.

Definition at line 190 of file G4Transportation.cc.

195{
196 // Initial actions moved to StartTrack()
197 // --------------------------------------
198 // Note: in case another process changes touchable handle
199 // it will be necessary to add here (for all steps)
200 // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
201
202 // GPILSelection is set to defaule value of CandidateForSelection
203 // It is a return value
204 //
205 *selection = CandidateForSelection;
206
207 // Get initial Energy/Momentum of the track
208 //
209 const G4ThreeVector startPosition = track.GetPosition();
210 const G4ThreeVector startMomentumDir = track.GetMomentumDirection();
211
212 // The Step Point safety can be limited by other geometries and/or the
213 // assumptions of any process - it's not always the geometrical safety.
214 // We calculate the starting point's isotropic safety here.
215 {
216 const G4double MagSqShift = (startPosition - fPreviousSftOrigin).mag2();
217
218 if(MagSqShift >= sqr(fPreviousSafety))
219 currentSafety = 0.0;
220 else
221 currentSafety = fPreviousSafety - std::sqrt(MagSqShift);
222 }
223
224 // Is the particle charged or has it a magnetic moment?
225 //
226 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
227
228 const G4double particleMass = pParticle->GetMass();
229 const G4double particleCharge = pParticle->GetCharge();
230 const G4double kineticEnergy = pParticle->GetKineticEnergy();
231
232 const G4double magneticMoment = pParticle->GetMagneticMoment();
233 const G4ThreeVector particleSpin = pParticle->GetPolarization();
234
235 // There is no need to locate the current volume. It is Done elsewhere:
236 // On track construction
237 // By the tracking, after all AlongStepDoIts, in "Relocation"
238
239 // Check if the particle has a force, EM or gravitational, exerted on it
240 //
241
242 G4bool eligibleEM =
243 (particleCharge != 0.0) || ((magneticMoment != 0.0) && fUseMagneticMoment);
244 G4bool eligibleGrav = (particleMass != 0.0) && fUseGravity;
245
246 fFieldExertedForce = false;
247
248 if(eligibleEM || eligibleGrav)
249 {
250 if(G4FieldManager* fieldMgr =
252 {
253 // User can configure the field Manager for this track
254 fieldMgr->ConfigureForTrack(&track);
255 // Called here to allow a transition from no-field pointer
256 // to finite field (non-zero pointer).
257
258 // If the field manager has no field ptr, the field is zero
259 // by definition ( = there is no field ! )
260 if(const G4Field* ptrField = fieldMgr->GetDetectorField())
262 eligibleEM || (eligibleGrav && ptrField->IsGravityActive());
263 }
264 }
265
266 G4double geometryStepLength = currentMinimumStep;
267
268 if(currentMinimumStep == 0.0)
269 {
270 fEndPointDistance = 0.0;
271 fGeometryLimitedStep = false; // Old code: = (currentSafety == 0.0);
272 // Changed to avoid problems when setting the step status (now done in AlongStepDoIt)
273 fMomentumChanged = false;
274 fParticleIsLooping = false;
276 fTransportEndPosition = startPosition;
277 fTransportEndMomentumDir = startMomentumDir;
278 fTransportEndKineticEnergy = kineticEnergy;
279 fTransportEndSpin = particleSpin;
280 }
281 else if(!fFieldExertedForce)
282 {
283 fGeometryLimitedStep = false;
284 if(geometryStepLength > currentSafety || !fShortStepOptimisation)
285 {
286 const G4double linearStepLength = fLinearNavigator->ComputeStep(
287 startPosition, startMomentumDir, currentMinimumStep, currentSafety);
288
289 if(linearStepLength <= currentMinimumStep)
290 {
291 geometryStepLength = linearStepLength;
293 }
294 // Remember last safety origin & value.
295 //
296 fPreviousSftOrigin = startPosition;
297 fPreviousSafety = currentSafety;
298 fpSafetyHelper->SetCurrentSafety(currentSafety, startPosition);
299 }
300
301 fEndPointDistance = geometryStepLength;
302
303 fMomentumChanged = false;
304 fParticleIsLooping = false;
307 startPosition + geometryStepLength * startMomentumDir;
308 fTransportEndMomentumDir = startMomentumDir;
309 fTransportEndKineticEnergy = kineticEnergy;
310 fTransportEndSpin = particleSpin;
311 }
312 else // A field exerts force
313 {
314 const auto pParticleDef = pParticle->GetDefinition();
315 const auto particlePDGSpin = pParticleDef->GetPDGSpin();
316 const auto particlePDGMagM = pParticleDef->GetPDGMagneticMoment();
317
318 auto equationOfMotion = fFieldPropagator->GetCurrentEquationOfMotion();
319
320 // The charge can change (dynamic), therefore the use of G4ChargeState
321 //
322 equationOfMotion->SetChargeMomentumMass(
323 G4ChargeState(particleCharge, magneticMoment, particlePDGSpin),
324 pParticle->GetTotalMomentum(), particleMass);
325
326 G4FieldTrack aFieldTrack(startPosition,
327 track.GetGlobalTime(), // Lab.
328 startMomentumDir, kineticEnergy, particleMass,
329 particleCharge, particleSpin, particlePDGMagM,
330 0.0, // Length along track
331 particlePDGSpin);
332
333 // Do the Transport in the field (non recti-linear)
334 //
335 const G4double lengthAlongCurve = fFieldPropagator->ComputeStep(
336 aFieldTrack, currentMinimumStep, currentSafety, track.GetVolume(),
337 kineticEnergy < fThreshold_Important_Energy);
338
339 if(lengthAlongCurve < geometryStepLength)
340 geometryStepLength = lengthAlongCurve;
341
342 // Remember last safety origin & value.
343 //
344 fPreviousSftOrigin = startPosition;
345 fPreviousSafety = currentSafety;
346 fpSafetyHelper->SetCurrentSafety(currentSafety, startPosition);
347
349 //
350 // It is possible that step was reduced in PropagatorInField due to
351 // previous zero steps. To cope with case that reduced step is taken
352 // in full, we must rely on PiF to obtain this value
353
354 G4bool changesEnergy =
356
357 fMomentumChanged = true;
359
360 fEndGlobalTimeComputed = changesEnergy;
361 fTransportEndPosition = aFieldTrack.GetPosition();
362 fTransportEndMomentumDir = aFieldTrack.GetMomentumDir();
363
364 // G4cout << " G4Transport: End of step pMag = " << aFieldTrack.GetMomentum().mag() << G4endl;
365
366 fEndPointDistance = (fTransportEndPosition - startPosition).mag();
367
368 // Ignore change in energy for fields that conserve energy
369 // This hides the integration error, but gives a better physical answer
371 changesEnergy ? aFieldTrack.GetKineticEnergy() : kineticEnergy;
372 fTransportEndSpin = aFieldTrack.GetSpin();
373
375 {
376 // If the field can change energy, then the time must be integrated
377 // - so this should have been updated
378 //
379 fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
380
381 // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
382 // a cleaner way is to have FieldTrack knowing whether time is updated.
383 }
384#if defined(G4VERBOSE) || defined(G4DEBUG_TRANSPORT)
385 else
386 {
387 // The energy should be unchanged by field transport,
388 // - so the time changed will be calculated elsewhere
389 //
390 // Check that the integration preserved the energy
391 // - and if not correct this!
392 G4double startEnergy = kineticEnergy;
394
395 static G4ThreadLocal G4int no_large_ediff;
396 if(verboseLevel > 1)
397 {
398 if(std::fabs(startEnergy - endEnergy) > perThousand * endEnergy)
399 {
400 static G4ThreadLocal G4int no_warnings = 0, warnModulo = 1,
401 moduloFactor = 10;
402 no_large_ediff++;
403 if((no_large_ediff % warnModulo) == 0)
404 {
405 no_warnings++;
406 std::ostringstream message;
407 message << "Energy change in Step is above 1^-3 relative value. "
408 << G4endl << " Relative change in 'tracking' step = "
409 << std::setw(15) << (endEnergy - startEnergy) / startEnergy
410 << G4endl << " Starting E= " << std::setw(12)
411 << startEnergy / MeV << " MeV " << G4endl
412 << " Ending E= " << std::setw(12) << endEnergy / MeV
413 << " MeV " << G4endl
414 << "Energy has been corrected -- however, review"
415 << " field propagation parameters for accuracy." << G4endl;
416 if((verboseLevel > 2) || (no_warnings < 4) ||
417 (no_large_ediff == warnModulo * moduloFactor))
418 {
419 message << "These include EpsilonStepMax(/Min) in G4FieldManager "
420 << G4endl
421 << "which determine fractional error per step for "
422 "integrated quantities. "
423 << G4endl
424 << "Note also the influence of the permitted number of "
425 "integration steps."
426 << G4endl;
427 }
428 message << "Bad 'endpoint'. Energy change detected and corrected."
429 << G4endl << "Has occurred already " << no_large_ediff
430 << " times.";
431 G4Exception("G4Transportation::AlongStepGetPIL()", "EnergyChange",
432 JustWarning, message);
433 if(no_large_ediff == warnModulo * moduloFactor)
434 {
435 warnModulo *= moduloFactor;
436 }
437 }
438 }
439 } // end of if (verboseLevel)
440 }
441#endif
442 }
443
444 // Update the safety starting from the end-point,
445 // if it will become negative at the end-point.
446 //
447 if(currentSafety < fEndPointDistance)
448 {
449 if(particleCharge != 0.0)
450 {
451 G4double endSafety =
453 currentSafety = endSafety;
455 fPreviousSafety = currentSafety;
457
458 // Because the Stepping Manager assumes it is from the start point,
459 // add the StepLength
460 //
461 currentSafety += fEndPointDistance;
462
463#ifdef G4DEBUG_TRANSPORT
464 G4cout.precision(12);
465 G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl;
466 G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
467 << " and it returned safety= " << endSafety << G4endl;
468 G4cout << " Adding endpoint distance " << fEndPointDistance
469 << " to obtain pseudo-safety= " << currentSafety << G4endl;
470 }
471 else
472 {
473 G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl;
474 G4cout << " Avoiding call to ComputeSafety : " << G4endl;
475 G4cout << " charge = " << particleCharge << G4endl;
476 G4cout << " mag moment = " << magneticMoment << G4endl;
477#endif
478 }
479 }
480
482 fLastStepInVolume = false;
483 fNewTrack = false;
484
486 fParticleChange.ProposeTrueStepLength(geometryStepLength);
487
488 return geometryStepLength;
489}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ CandidateForSelection
G4double GetCharge() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetMagneticMoment() const
G4double GetTotalMomentum() const
const G4ThreeVector & GetPolarization() const
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
G4bool DoesFieldChangeEnergy() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4FieldManager * GetCurrentFieldManager()
G4bool IsLastStepInVolume()
G4bool IsParticleLooping() const
G4EquationOfMotion * GetCurrentEquationOfMotion()
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=nullptr, G4bool canRelaxDeltaChord=false)
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
static G4bool fUseGravity
static G4bool fUseMagneticMoment
void ProposeFirstStepInVolume(G4bool flag)
void ProposeTrueStepLength(G4double truePathLength)
T sqr(const T &x)
Definition templates.hh:128

Referenced by G4TransportationWithMsc::AlongStepGetPhysicalInteractionLength().

◆ AtRestDoIt()

G4VParticleChange * G4Transportation::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inlinevirtual

Implements G4VProcess.

Definition at line 144 of file G4Transportation.hh.

145 { return 0; } // No operation in AtRestDoIt

◆ AtRestGetPhysicalInteractionLength()

G4double G4Transportation::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inlinevirtual

Implements G4VProcess.

Definition at line 140 of file G4Transportation.hh.

142 { return -1.0; } // No operation in AtRestGPIL

◆ EnableGravity()

G4bool G4Transportation::EnableGravity ( G4bool useGravity = true)
static

Definition at line 845 of file G4Transportation.cc.

846{
847 G4bool lastValue= fUseGravity;
848 fUseGravity= useGravity;
849 return lastValue;
850}

◆ EnableMagneticMoment()

G4bool G4Transportation::EnableMagneticMoment ( G4bool useMoment = true)
static

Definition at line 835 of file G4Transportation.cc.

836{
837 G4bool lastValue= fUseMagneticMoment;
838 fUseMagneticMoment= useMoment;
839 return lastValue;
840}

Referenced by G4CoupledTransportation::EnableUseMagneticMoment(), and EnableUseMagneticMoment().

◆ EnableShortStepOptimisation()

void G4Transportation::EnableShortStepOptimisation ( G4bool optimise = true)
inline

◆ EnableUseMagneticMoment()

static G4bool G4Transportation::EnableUseMagneticMoment ( G4bool useMoment = true)
inlinestatic

Definition at line 135 of file G4Transportation.hh.

136 { return EnableMagneticMoment(useMoment); } // Old name - will be deprecated
static G4bool EnableMagneticMoment(G4bool useMoment=true)

◆ FieldExertedForce()

G4bool G4Transportation::FieldExertedForce ( )
inline

Definition at line 93 of file G4Transportation.hh.

93{ return fFieldExertedForce; }

◆ GetMaxEnergyKilled()

G4double G4Transportation::GetMaxEnergyKilled ( ) const
inline

◆ GetPropagatorInField()

G4PropagatorInField * G4Transportation::GetPropagatorInField ( )

◆ GetSilenceLooperWarnings()

G4bool G4Transportation::GetSilenceLooperWarnings ( )
static

Definition at line 863 of file G4Transportation.cc.

864{
866}

◆ GetSumEnergyKilled()

G4double G4Transportation::GetSumEnergyKilled ( ) const
inline

◆ GetThresholdImportantEnergy()

G4double G4Transportation::GetThresholdImportantEnergy ( ) const
inline

◆ GetThresholdTrials()

G4int G4Transportation::GetThresholdTrials ( ) const
inline

◆ GetThresholdWarningEnergy()

G4double G4Transportation::GetThresholdWarningEnergy ( ) const
inline

◆ PostStepDoIt()

G4VParticleChange * G4Transportation::PostStepDoIt ( const G4Track & track,
const G4Step & stepData )
virtual

Implements G4VProcess.

Definition at line 708 of file G4Transportation.cc.

710{
711 G4TouchableHandle retCurrentTouchable ; // The one to return
712 G4bool isLastStep= false;
713
714 // Initialize ParticleChange (by setting all its members equal
715 // to corresponding members in G4Track)
716 // fParticleChange.Initialize(track) ; // To initialise TouchableChange
717
719
720 // If the Step was determined by the volume boundary,
721 // logically relocate the particle
722
724 {
725 // fCurrentTouchable will now become the previous touchable,
726 // and what was the previous will be freed.
727 // (Needed because the preStepPoint can point to the previous touchable)
728
731 LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
732 track.GetMomentumDirection(),
734 true ) ;
735 // Check whether the particle is out of the world volume
736 // If so it has exited and must be killed.
737 //
739 {
741 }
742 retCurrentTouchable = fCurrentTouchableHandle ;
744
745 // Update the Step flag which identifies the Last Step in a volume
746 if( !fFieldExertedForce )
747 isLastStep = fLinearNavigator->ExitedMotherVolume()
749 else
750 isLastStep = fFieldPropagator->IsLastStepInVolume();
751 }
752 else // fGeometryLimitedStep is false
753 {
754 // This serves only to move the Navigator's location
755 //
757
758 // The value of the track's current Touchable is retained.
759 // (and it must be correct because we must use it below to
760 // overwrite the (unset) one in particle change)
761 // It must be fCurrentTouchable too ??
762 //
764 retCurrentTouchable = track.GetTouchableHandle() ;
765
766 isLastStep= false;
767 } // endif ( fGeometryLimitedStep )
768 fLastStepInVolume= isLastStep;
769
772
773 SetTouchableInformation(retCurrentTouchable);
774
775 return &fParticleChange ;
776}
void SetGeometricallyLimitedStep()
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
G4bool EnteredDaughterVolume() const
G4bool ExitedMotherVolume() const
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
G4TrackStatus GetTrackStatus() const
const G4TouchableHandle & GetTouchableHandle() const
void SetTouchableInformation(const G4TouchableHandle &touchable)
void ProposeLastStepInVolume(G4bool flag)

◆ PostStepGetPhysicalInteractionLength()

G4double G4Transportation::PostStepGetPhysicalInteractionLength ( const G4Track & ,
G4double previousStepSize,
G4ForceCondition * pForceCond )
virtual

Implements G4VProcess.

Definition at line 651 of file G4Transportation.cc.

655{
656 fFieldExertedForce = false; // Not known
657 *pForceCond = Forced ;
658 return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
659}
@ Forced
#define DBL_MAX
Definition templates.hh:62

◆ PrintStatistics()

void G4Transportation::PrintStatistics ( std::ostream & outStr) const

Definition at line 148 of file G4Transportation.cc.

149{
150 outStr << " G4Transportation: Statistics for looping particles " << G4endl;
151 if( fSumEnergyKilled > 0.0 || fNumLoopersKilled > 0 )
152 {
153 outStr << " Sum of energy of looping tracks killed: "
154 << fSumEnergyKilled / CLHEP::MeV << " MeV "
155 << " from " << fNumLoopersKilled << " tracks " << G4endl
156 << " Sum of energy of non-electrons : "
157 << fSumEnergyKilled_NonElectron / CLHEP::MeV << " MeV "
158 << " from " << fNumLoopersKilled_NonElectron << " tracks "
159 << G4endl;
160 outStr << " Max energy of *any type* looper killed: " << fMaxEnergyKilled
161 << " its PDG was " << fMaxEnergyKilledPDG << G4endl;
163 {
164 outStr << " Max energy of non-electron looper killed: "
166 << " its PDG was " << fMaxEnergyKilled_NonElecPDG << G4endl;
167 }
168 if( fMaxEnergySaved > 0.0 )
169 {
170 outStr << " Max energy of loopers 'saved': " << fMaxEnergySaved << G4endl;
171 outStr << " Sum of energy of loopers 'saved': "
173 outStr << " Sum of energy of unstable loopers 'saved': "
175 }
176 }
177 else
178 {
179 outStr << " No looping tracks found or killed. " << G4endl;
180 }
181}

Referenced by ~G4Transportation().

◆ ProcessDescription()

void G4Transportation::ProcessDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4VProcess.

Definition at line 922 of file G4Transportation.cc.

926{
927 G4String indent = " "; // : "");
928 G4long oldPrec= outStr.precision(6);
929 // outStr << std::setprecision(6);
930 outStr << G4endl << indent << GetProcessName() << ": ";
931
932 outStr << " Parameters for looping particles: " << G4endl
933 << " warning-E = " << fThreshold_Warning_Energy / CLHEP::MeV << " MeV " << G4endl
934 << " important E = " << fThreshold_Important_Energy / CLHEP::MeV << " MeV " << G4endl
935 << " thresholdTrials " << fThresholdTrials << G4endl;
936 outStr.precision(oldPrec);
937}
const G4String & GetProcessName() const

◆ PushThresholdsToLogger()

void G4Transportation::PushThresholdsToLogger ( )

◆ ReportLooperThresholds()

void G4Transportation::ReportLooperThresholds ( )

Definition at line 914 of file G4Transportation.cc.

915{
916 PushThresholdsToLogger(); // To be absolutely certain they are in sync
917 fpLogger->ReportLooperThresholds("G4Transportation");
918}
void ReportLooperThresholds(const char *className)

Referenced by SetHighLooperThresholds(), and SetLowLooperThresholds().

◆ ReportMissingLogger()

void G4Transportation::ReportMissingLogger ( const char * methodName)
protected

Definition at line 903 of file G4Transportation.cc.

904{
905 const char* message= "Logger object missing from G4Transportation object";
906 G4String classAndMethod= G4String("G4Transportation") + G4String( methodName );
907 G4Exception(classAndMethod, "Missing Logger", JustWarning, message);
908}

◆ ResetKilledStatistics()

void G4Transportation::ResetKilledStatistics ( G4int report = 1)
inline

◆ SetHighLooperThresholds()

void G4Transportation::SetHighLooperThresholds ( )

Definition at line 871 of file G4Transportation.cc.

872{
873 // Restores the old high values -- potentially appropriate for energy-frontier
874 // HEP experiments.
875 // Caution: All tracks with E < 100 MeV that are found to loop are
876 SetThresholdWarningEnergy( 100.0 * CLHEP::MeV ); // Warn above this energy
877 SetThresholdImportantEnergy( 250.0 * CLHEP::MeV ); // Extra trial above this En
878
879 G4int maxTrials = 10;
880 SetThresholdTrials( maxTrials );
881
882 PushThresholdsToLogger(); // Again, to be sure
884}

Referenced by G4Transportation().

◆ SetLowLooperThresholds()

void G4Transportation::SetLowLooperThresholds ( )

Definition at line 887 of file G4Transportation.cc.

888{
889 // These values were the default in Geant4 10.5 - beta
890 SetThresholdWarningEnergy( 1.0 * CLHEP::keV ); // Warn above this En
891 SetThresholdImportantEnergy( 1.0 * CLHEP::MeV ); // Extra trials above it
892
893 G4int maxTrials = 30; // A new value - was 10
894 SetThresholdTrials( maxTrials );
895
896 PushThresholdsToLogger(); // Again, to be sure
898}

◆ SetPropagatorInField()

void G4Transportation::SetPropagatorInField ( G4PropagatorInField * pFieldPropagator)

◆ SetSilenceLooperWarnings()

void G4Transportation::SetSilenceLooperWarnings ( G4bool val)
static

Definition at line 856 of file G4Transportation.cc.

857{
858 fSilenceLooperWarnings= val; // Flag to *Supress* all 'looper' warnings
859}

◆ SetThresholdImportantEnergy()

void G4Transportation::SetThresholdImportantEnergy ( G4double newEnImp)
inline

◆ SetThresholdTrials()

void G4Transportation::SetThresholdTrials ( G4int newMaxTrials)
inline

◆ SetThresholdWarningEnergy()

void G4Transportation::SetThresholdWarningEnergy ( G4double newEnWarn)
inline

◆ SetTouchableInformation()

void G4Transportation::SetTouchableInformation ( const G4TouchableHandle & touchable)
protected

Definition at line 663 of file G4Transportation.cc.

664{
665 const G4VPhysicalVolume* pNewVol = touchable->GetVolume() ;
666 const G4Material* pNewMaterial = 0 ;
667 G4VSensitiveDetector* pNewSensitiveDetector = 0;
668
669 if( pNewVol != 0 )
670 {
671 pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
672 pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
673 }
674
676 fParticleChange.SetSensitiveDetectorInTouchable( pNewSensitiveDetector ) ;
677 // temporarily until Get/Set Material of ParticleChange,
678 // and StepPoint can be made const.
679
680 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
681 if( pNewVol != 0 )
682 {
683 pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
684 }
685
686 if ( pNewVol!=0 && pNewMaterialCutsCouple!=0
687 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
688 {
689 // for parametrized volume
690 //
691 pNewMaterialCutsCouple =
693 ->GetMaterialCutsCouple(pNewMaterial,
694 pNewMaterialCutsCouple->GetProductionCuts());
695 }
696 fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
697
698 // Set the touchable in ParticleChange
699 // this must always be done because the particle change always
700 // uses this value to overwrite the current touchable pointer.
701 //
703}
G4VSensitiveDetector * GetSensitiveDetector() const
G4Material * GetMaterial() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
void SetMaterialInTouchable(G4Material *fMaterial)
void SetMaterialCutsCoupleInTouchable(const G4MaterialCutsCouple *fMaterialCutsCouple)
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4LogicalVolume * GetLogicalVolume() const

Referenced by G4CoupledTransportation::PostStepDoIt(), and PostStepDoIt().

◆ StartTracking()

void G4Transportation::StartTracking ( G4Track * aTrack)
virtual

Reimplemented from G4VProcess.

Reimplemented in G4TransportationWithMsc.

Definition at line 785 of file G4Transportation.cc.

786{
788 fNewTrack= true;
789 fFirstStepInVolume= true;
790 fLastStepInVolume= false;
791
792 // The actions here are those that were taken in AlongStepGPIL
793 // when track.GetCurrentStepNumber()==1
794
795 // Whether field exists should be determined at run level -- TODO
797 fAnyFieldExists = fieldMgrStore->size() > 0;
798
799 // reset safety value and center
800 //
801 fPreviousSafety = 0.0 ;
803
804 // reset looping counter -- for motion in field
805 fNoLooperTrials= 0;
806 // Must clear this state .. else it depends on last track's value
807 // --> a better solution would set this from state of suspended track TODO ?
808 // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
809
810 // ChordFinder reset internal state
811 //
813 {
815 // Resets all state of field propagator class (ONLY) including safety
816 // values (in case of overlaps and to wipe for first track).
817 }
818
819 // Make sure to clear the chord finders of all fields (i.e. managers)
820 //
821 fieldMgrStore->ClearAllChordFindersState();
822
823 // Update the current touchable handle (from the track's)
824 //
826
827 // Inform field propagator of new track
828 //
830}
CLHEP::Hep3Vector G4ThreeVector
static G4FieldManagerStore * GetInstance()
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87

Referenced by G4CoupledTransportation::StartTracking(), and G4TransportationWithMsc::StartTracking().

Member Data Documentation

◆ fAbandonUnstableTrials

G4int G4Transportation::fAbandonUnstableTrials = 0
protected

Definition at line 206 of file G4Transportation.hh.

Referenced by AlongStepDoIt().

◆ fAnyFieldExists

G4bool G4Transportation::fAnyFieldExists = false
protected

Definition at line 176 of file G4Transportation.hh.

Referenced by StartTracking().

◆ fCandidateEndGlobalTime

G4double G4Transportation::fCandidateEndGlobalTime = 0.0
protected

◆ fCurrentTouchableHandle

G4TouchableHandle G4Transportation::fCurrentTouchableHandle
protected

◆ fEndGlobalTimeComputed

G4bool G4Transportation::fEndGlobalTimeComputed = false
protected

◆ fEndPointDistance

◆ fFieldExertedForce

G4bool G4Transportation::fFieldExertedForce = false
protected

◆ fFieldPropagator

◆ fFirstStepInVolume

◆ fGeometryLimitedStep

◆ fLastStepInVolume

G4bool G4Transportation::fLastStepInVolume = false
protected

◆ fLinearNavigator

◆ fMaxEnergyKilled

G4double G4Transportation::fMaxEnergyKilled = -1.0
protected

Definition at line 216 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fMaxEnergyKilled_NonElecPDG

G4int G4Transportation::fMaxEnergyKilled_NonElecPDG = 0
protected

Definition at line 222 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fMaxEnergyKilled_NonElectron

G4double G4Transportation::fMaxEnergyKilled_NonElectron = -1.0
protected

Definition at line 221 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fMaxEnergyKilledPDG

G4int G4Transportation::fMaxEnergyKilledPDG = 0
protected

Definition at line 217 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fMaxEnergySaved

G4double G4Transportation::fMaxEnergySaved = -1.0
protected

Definition at line 225 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fMomentumChanged

◆ fNewTrack

G4bool G4Transportation::fNewTrack = true
protected

◆ fNoLooperTrials

G4int G4Transportation::fNoLooperTrials = 0
protected

Definition at line 210 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and StartTracking().

◆ fNumLoopersKilled

unsigned long G4Transportation::fNumLoopersKilled = 0
protected

Definition at line 218 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fNumLoopersKilled_NonElectron

unsigned long G4Transportation::fNumLoopersKilled_NonElectron = 0
protected

Definition at line 223 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fParticleChange

◆ fParticleIsLooping

G4bool G4Transportation::fParticleIsLooping = false
protected

◆ fpLogger

G4TransportationLogger* G4Transportation::fpLogger
protected

◆ fPreviousSafety

G4double G4Transportation::fPreviousSafety
protected

Definition at line 191 of file G4Transportation.hh.

Referenced by AlongStepGetPhysicalInteractionLength(), and StartTracking().

◆ fPreviousSftOrigin

◆ fpSafetyHelper

◆ fShortStepOptimisation

G4bool G4Transportation::fShortStepOptimisation
protected

Definition at line 230 of file G4Transportation.hh.

Referenced by AlongStepGetPhysicalInteractionLength(), and G4Transportation().

◆ fSilenceLooperWarnings

G4bool G4Transportation::fSilenceLooperWarnings = false
staticprotected

◆ fSumEnergyKilled

G4double G4Transportation::fSumEnergyKilled = 0.0
protected

Definition at line 214 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), PrintStatistics(), and ~G4Transportation().

◆ fSumEnergyKilled_NonElectron

G4double G4Transportation::fSumEnergyKilled_NonElectron = 0.0
protected

Definition at line 219 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fSumEnergySaved

G4double G4Transportation::fSumEnergySaved = 0.0
protected

Definition at line 224 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fSumEnergyUnstableSaved

G4double G4Transportation::fSumEnergyUnstableSaved = 0.0
protected

Definition at line 226 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and PrintStatistics().

◆ fSumEnerSqKilled

G4double G4Transportation::fSumEnerSqKilled = 0.0
protected

Definition at line 215 of file G4Transportation.hh.

Referenced by AlongStepDoIt().

◆ fSumEnerSqKilled_NonElectron

G4double G4Transportation::fSumEnerSqKilled_NonElectron = 0.0
protected

Definition at line 220 of file G4Transportation.hh.

Referenced by AlongStepDoIt().

◆ fThreshold_Important_Energy

G4double G4Transportation::fThreshold_Important_Energy = 1.0 * CLHEP::MeV
protected

Definition at line 202 of file G4Transportation.hh.

Referenced by AlongStepDoIt(), and AlongStepGetPhysicalInteractionLength().

◆ fThreshold_Warning_Energy

G4double G4Transportation::fThreshold_Warning_Energy = 1.0 * CLHEP::keV
protected

Definition at line 201 of file G4Transportation.hh.

Referenced by AlongStepDoIt().

◆ fThresholdTrials

G4int G4Transportation::fThresholdTrials = 10
protected

Definition at line 203 of file G4Transportation.hh.

Referenced by AlongStepDoIt().

◆ fTransportEndKineticEnergy

◆ fTransportEndMomentumDir

◆ fTransportEndPosition

◆ fTransportEndSpin

G4ThreeVector G4Transportation::fTransportEndSpin = G4ThreeVector( 0.0, 0.0, 0.0 )
protected

◆ fUseGravity

G4bool G4Transportation::fUseGravity = false
staticprotected

◆ fUseMagneticMoment

G4bool G4Transportation::fUseMagneticMoment =false
staticprotected

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