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

#include <G4PropagatorInField.hh>

Public Member Functions

void CheckMode (G4bool mode)
 
void ClearPropagatorState ()
 
G4double ComputeStep (G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=nullptr, G4bool canRelaxDeltaChord=false)
 
G4ThreeVector EndMomentumDir () const
 
G4ThreeVector EndPosition () const
 
G4FieldManagerFindAndSetFieldManager (G4VPhysicalVolume *pCurrentPhysVol)
 
 G4PropagatorInField (G4Navigator *theNavigator, G4FieldManager *detectorFieldMgr, G4VIntersectionLocator *vLocator=nullptr)
 
G4ChordFinderGetChordFinder ()
 
G4EquationOfMotionGetCurrentEquationOfMotion ()
 
G4FieldManagerGetCurrentFieldManager ()
 
G4double GetDeltaIntersection () const
 
G4double GetDeltaOneStep () const
 
G4FieldTrack GetEndState () const
 
G4double GetEpsilonStep () const
 
G4VIntersectionLocatorGetIntersectionLocator ()
 
G4int GetIterationsToIncreaseChordDistance () const
 
G4double GetLargestAcceptableStep ()
 
G4double GetMaximumEpsilonStep () const
 
G4int GetMaxLoopCount () const
 
G4double GetMinimumEpsilonStep () const
 
G4NavigatorGetNavigatorForPropagating ()
 
G4int GetThresholdNoZeroSteps (G4int i)
 
G4bool GetUseSafetyForOptimization ()
 
G4int GetVerboseLevel () const
 
G4bool GetVerboseTrace ()
 
G4double GetZeroStepThreshold ()
 
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt () const
 
G4bool IntersectChord (const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
 
G4bool IsFirstStepInVolume ()
 
G4bool IsLastStepInVolume ()
 
G4bool IsParticleLooping () const
 
void PrepareNewTrack ()
 
void printStatus (const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
 
void RefreshIntersectionLocator ()
 
void SetDetectorFieldManager (G4FieldManager *newGlobalFieldManager)
 
void SetEpsilonStep (G4double newEps)
 
void SetIntersectionLocator (G4VIntersectionLocator *pLocator)
 
void SetIterationsToIncreaseChordDistance (G4int numIters)
 
void SetLargestAcceptableStep (G4double newBigDist)
 
void SetMaximumEpsilonStep (G4double newEpsMax)
 
void SetMaxLoopCount (G4int new_max)
 
void SetMinimumEpsilonStep (G4double newEpsMin)
 
void SetNavigatorForPropagating (G4Navigator *SimpleOrMultiNavigator)
 
void SetThresholdNoZeroStep (G4int noAct, G4int noHarsh, G4int noAbandon)
 
void SetTrajectoryFilter (G4VCurvedTrajectoryFilter *filter)
 
void SetUseSafetyForOptimization (G4bool)
 
G4int SetVerboseLevel (G4int verbose)
 
void SetVerboseTrace (G4bool enable)
 
void SetZeroStepThreshold (G4double newLength)
 
G4int Verbose () const
 
 ~G4PropagatorInField ()
 

Protected Member Functions

void PrintStepLengthDiagnostic (G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
 
void ReportLoopingParticle (G4int count, G4double StepTaken, G4double stepRequest, const char *methodName, G4ThreeVector momentumVec, G4VPhysicalVolume *physVol)
 
void ReportStuckParticle (G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)
 

Private Attributes

G4FieldTrack End_PointAndTangent
 
G4int fAbandonThreshold_NoZeroSteps = 50
 
G4int fActionThreshold_NoZeroSteps = 2
 
G4bool fAllocatedLocator
 
G4bool fCheck = false
 
G4FieldManagerfCurrentFieldMgr
 
G4FieldManagerfDetectorFieldMgr
 
G4double fEpsilonStep
 
G4bool fFirstStepInVolume = true
 
G4double fFull_CurveLen_of_LastAttempt = -1
 
G4int fIncreaseChordDistanceThreshold = 100
 
G4VIntersectionLocatorfIntersectionLocator
 
G4double fLargestAcceptableStep
 
G4double fLast_ProposedStepLength = -1
 
G4bool fLastStepInVolume = true
 
G4int fMax_loop_count = 1000
 
G4NavigatorfNavigator
 
G4bool fNewTrack = true
 
G4int fNoZeroStep = 0
 
G4bool fParticleIsLooping = false
 
G4double fPreviousSafety = 0.0
 
G4ThreeVector fPreviousSftOrigin
 
G4VCurvedTrajectoryFilterfpTrajectoryFilter = nullptr
 
G4bool fSetFieldMgr = false
 
G4int fSevereActionThreshold_NoZeroSteps = 10
 
G4bool fUseSafetyForOptimisation = true
 
G4int fVerboseLevel = 0
 
G4bool fVerbTracePiF = false
 
G4double fZeroStepThreshold = 0.0
 
G4double kCarTolerance
 

Detailed Description

Definition at line 57 of file G4PropagatorInField.hh.

Constructor & Destructor Documentation

◆ G4PropagatorInField()

G4PropagatorInField::G4PropagatorInField ( G4Navigator theNavigator,
G4FieldManager detectorFieldMgr,
G4VIntersectionLocator vLocator = nullptr 
)

Definition at line 55 of file G4PropagatorInField.cc.

58 : fDetectorFieldMgr(detectorFieldMgr),
59 fNavigator(theNavigator),
60 fCurrentFieldMgr(detectorFieldMgr),
62 G4ThreeVector(0.,0.,0.),0.0,0.0,0.0,0.0,0.0)
63{
64 fEpsilonStep = (fDetectorFieldMgr != nullptr)
67
71
72#ifdef G4DEBUG_FIELD
73 G4cout << " PiF: Zero Step Threshold set to "
75 << " mm." << G4endl;
76 G4cout << " PiF: Value of kCarTolerance = "
78 << " mm. " << G4endl;
79 fVerboseLevel = 2;
80 fVerbTracePiF = true;
81#endif
82
83 // Defining Intersection Locator and his parameters
84 if ( vLocator == nullptr )
85 {
86 fIntersectionLocator = new G4MultiLevelLocator(theNavigator);
87 fAllocatedLocator = true;
88 }
89 else
90 {
91 fIntersectionLocator = vLocator;
92 fAllocatedLocator = false;
93 }
94 RefreshIntersectionLocator(); // Copy all relevant parameters
95}
static constexpr double micrometer
Definition: G4SIunits.hh:80
static constexpr double millimeter
Definition: G4SIunits.hh:66
static constexpr double meter
Definition: G4SIunits.hh:62
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetMaximumEpsilonStep() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4FieldManager * fCurrentFieldMgr
G4ThreeVector fPreviousSftOrigin
G4VIntersectionLocator * fIntersectionLocator
G4FieldTrack End_PointAndTangent
G4FieldManager * fDetectorFieldMgr
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References fAllocatedLocator, fDetectorFieldMgr, fEpsilonStep, fIntersectionLocator, fLargestAcceptableStep, fPreviousSftOrigin, fVerboseLevel, fVerbTracePiF, fZeroStepThreshold, G4cout, G4endl, G4GeometryTolerance::GetInstance(), G4FieldManager::GetMaximumEpsilonStep(), G4GeometryTolerance::GetSurfaceTolerance(), kCarTolerance, G4INCL::Math::max(), meter, micrometer, millimeter, and RefreshIntersectionLocator().

◆ ~G4PropagatorInField()

G4PropagatorInField::~G4PropagatorInField ( )

Definition at line 99 of file G4PropagatorInField.cc.

100{
102}

References fAllocatedLocator, and fIntersectionLocator.

Member Function Documentation

◆ CheckMode()

void G4PropagatorInField::CheckMode ( G4bool  mode)
inline

◆ ClearPropagatorState()

void G4PropagatorInField::ClearPropagatorState ( )

Definition at line 674 of file G4PropagatorInField.cc.

675{
676 // Goal: Clear all memory of previous steps, cached information
677
678 fParticleIsLooping = false;
679 fNoZeroStep = 0;
680
681 fSetFieldMgr = false; // Has field-manager been set for the current step?
682 fEpsilonStep= 1.0e-5; // Relative accuracy of current Step
683
685 G4ThreeVector(0.,0.,0.),
686 0.0,0.0,0.0,0.0,0.0);
689
691 fPreviousSafety= 0.0;
692
693 fNewTrack = true;
694}

References End_PointAndTangent, fEpsilonStep, fFull_CurveLen_of_LastAttempt, fLast_ProposedStepLength, fNewTrack, fNoZeroStep, fParticleIsLooping, fPreviousSafety, fPreviousSftOrigin, and fSetFieldMgr.

Referenced by G4ITTransportation::StartTracking(), G4CoupledTransportation::StartTracking(), and G4Transportation::StartTracking().

◆ ComputeStep()

G4double G4PropagatorInField::ComputeStep ( G4FieldTrack pFieldTrack,
G4double  pCurrentProposedStepLength,
G4double pNewSafety,
G4VPhysicalVolume pPhysVol = nullptr,
G4bool  canRelaxDeltaChord = false 
)

Definition at line 118 of file G4PropagatorInField.cc.

124{
126 const G4double deltaChord = GetChordFinder()->GetDeltaChord();
127
128 // If CurrentProposedStepLength is too small for finding Chords
129 // then return with no action (for now - TODO: some action)
130 //
131 const char* methodName = "G4PropagatorInField::ComputeStep";
132 if (CurrentProposedStepLength<kCarTolerance)
133 {
134 return kInfinity;
135 }
136
137 // Introducing smooth trajectory display (jacek 01/11/2002)
138 //
140 {
142 }
143
145 fLastStepInVolume = false;
146 fNewTrack = false;
147
148 if( fVerboseLevel > 2 )
149 {
150 G4cout << methodName << " called" << G4endl;
151 G4cout << " Starting FT: " << pFieldTrack;
152 G4cout << " Requested length = " << CurrentProposedStepLength << G4endl;
153 G4cout << " PhysVol = ";
154 if( pPhysVol )
155 G4cout << pPhysVol->GetName() << G4endl;
156 else
157 G4cout << " N/A ";
158 G4cout << G4endl;
159 }
160
161 // Parameters for adaptive Runge-Kutta integration
162
163 G4double h_TrialStepSize; // 1st Step Size
164 G4double TruePathLength = CurrentProposedStepLength;
165 G4double StepTaken = 0.0;
166 G4double s_length_taken, epsilon;
167 G4bool intersects;
168 G4bool first_substep = true;
169
170 G4double NewSafety;
171 fParticleIsLooping = false;
172
173 // If not yet done,
174 // Set the field manager to the local one if the volume has one,
175 // or to the global one if not
176 //
177 if( !fSetFieldMgr )
178 {
180 }
181 fSetFieldMgr = false; // For next call, the field manager must be set again
182
183 G4FieldTrack CurrentState(pFieldTrack);
184 G4FieldTrack OriginalState = CurrentState;
185
186 // If the Step length is "infinite", then an approximate-maximum Step
187 // length (used to calculate the relative accuracy) must be guessed
188 //
189 if( CurrentProposedStepLength >= fLargestAcceptableStep )
190 {
191 G4ThreeVector StartPointA, VelocityUnit;
192 StartPointA = pFieldTrack.GetPosition();
193 VelocityUnit = pFieldTrack.GetMomentumDir();
194
195 G4double trialProposedStep = 1.e2 * ( 10.0 * cm +
197 GetSolid()->DistanceToOut(StartPointA, VelocityUnit) );
198 CurrentProposedStepLength = std::min( trialProposedStep,
200 }
201 epsilon = fCurrentFieldMgr->GetDeltaOneStep() / CurrentProposedStepLength;
204 if( epsilon < epsilonMin ) { epsilon = epsilonMin; }
205 if( epsilon > epsilonMax ) { epsilon = epsilonMax; }
207
208 // Values for Intersection Locator has to be updated on each call for the
209 // case that CurrentFieldManager has changed from the one of previous step
210 //
212
213 // Shorten the proposed step in case of earlier problems (zero steps)
214 //
216 {
217 G4double stepTrial;
218
220 if( (stepTrial <= 0.0) && (fLast_ProposedStepLength > 0.0) )
221 {
222 stepTrial = fLast_ProposedStepLength;
223 }
224
225 G4double decreaseFactor = 0.9; // Unused default
227 && (stepTrial > 100.0*fZeroStepThreshold) )
228 {
229 // Attempt quick convergence
230 //
231 decreaseFactor= 0.25;
232 }
233 else
234 {
235 // We are in significant difficulties, probably at a boundary that
236 // is either geometrically sharp or between very different materials.
237 // Careful decreases to cope with tolerance are required
238 //
239 if( stepTrial > 100.0*fZeroStepThreshold )
240 decreaseFactor = 0.35; // Try decreasing slower
241 else if( stepTrial > 30.0*fZeroStepThreshold )
242 decreaseFactor= 0.5; // Try yet slower decrease
243 else if( stepTrial > 10.0*fZeroStepThreshold )
244 decreaseFactor= 0.75; // Try even slower decreases
245 else
246 decreaseFactor= 0.9; // Try very slow decreases
247 }
248 stepTrial *= decreaseFactor;
249
250#ifdef G4DEBUG_FIELD
251 if( fVerboseLevel > 2
253 {
254 G4cerr << " " << methodName
255 << " Decreasing step after " << fNoZeroStep << " zero steps "
256 << " - in volume " << pPhysVol;
257 if( pPhysVol )
258 G4cerr << " with name " << pPhysVol->GetName();
259 else
260 G4cerr << " i.e. *unknown* volume.";
261 G4cerr << G4endl;
262 PrintStepLengthDiagnostic(CurrentProposedStepLength, decreaseFactor,
263 stepTrial, pFieldTrack);
264 }
265#endif
266 if( stepTrial == 0.0 ) // Change to make it < 0.1 * kCarTolerance ??
267 {
268 std::ostringstream message;
269 message << "Particle abandoned due to lack of progress in field."
270 << G4endl
271 << " Properties : " << pFieldTrack << G4endl
272 << " Attempting a zero step = " << stepTrial << G4endl
273 << " while attempting to progress after " << fNoZeroStep
274 << " trial steps. Will abandon step.";
275 G4Exception(methodName, "GeomNav1002", JustWarning, message);
276 fParticleIsLooping = true;
277 return 0; // = stepTrial;
278 }
279 if( stepTrial < CurrentProposedStepLength )
280 {
281 CurrentProposedStepLength = stepTrial;
282 }
283 }
284 fLast_ProposedStepLength = CurrentProposedStepLength;
285
286 G4int do_loop_count = 0;
287 do // Loop checking, 07.10.2016, JA
288 {
289 G4FieldTrack SubStepStartState = CurrentState;
290 G4ThreeVector SubStartPoint = CurrentState.GetPosition();
291
292 if(!first_substep)
293 {
294 if( fVerboseLevel > 4 )
295 {
296 G4cout << " PiF: Calling Nav/Locate Global Point within-Volume "
297 << G4endl;
298 }
300 }
301
302 // How far to attempt to move the particle !
303 //
304 h_TrialStepSize = CurrentProposedStepLength - StepTaken;
305
306 if (canRelaxDeltaChord &&
308 do_loop_count > fIncreaseChordDistanceThreshold &&
309 do_loop_count % fIncreaseChordDistanceThreshold == 0)
310 {
312 GetChordFinder()->GetDeltaChord() * 2.0
313 );
314 }
315
316 // Integrate as far as "chord miss" rule allows.
317 //
318 s_length_taken = GetChordFinder()->AdvanceChordLimited(
319 CurrentState, // Position & velocity
320 h_TrialStepSize,
324 // CurrentState is now updated with the final position and velocity
325
326 fFull_CurveLen_of_LastAttempt = s_length_taken;
327
328 G4ThreeVector EndPointB = CurrentState.GetPosition();
329 G4ThreeVector InterSectionPointE;
330 G4double LinearStepLength;
331
332 // Intersect chord AB with geometry
333 //
334 intersects= IntersectChord( SubStartPoint, EndPointB,
335 NewSafety, LinearStepLength,
336 InterSectionPointE );
337 // E <- Intersection Point of chord AB and either volume A's surface
338 // or a daughter volume's surface ..
339
340 if( first_substep )
341 {
342 currentSafety = NewSafety;
343 } // Updating safety in other steps is potential future extention
344
345 if( intersects )
346 {
347 G4FieldTrack IntersectPointVelct_G(CurrentState); // FT-Def-Construct
348
349 // Find the intersection point of AB true path with the surface
350 // of vol(A), if it exists. Start with point E as first "estimate".
351 G4bool recalculatedEndPt = false;
352
353 G4bool found_intersection = fIntersectionLocator->
354 EstimateIntersectionPoint( SubStepStartState, CurrentState,
355 InterSectionPointE, IntersectPointVelct_G,
356 recalculatedEndPt, fPreviousSafety,
358 intersects = found_intersection;
359 if( found_intersection )
360 {
361 End_PointAndTangent= IntersectPointVelct_G; // G is our EndPoint ...
362 StepTaken = TruePathLength = IntersectPointVelct_G.GetCurveLength()
363 - OriginalState.GetCurveLength();
364 }
365 else
366 {
367 // Either "minor" chords do not intersect
368 // or else stopped (due to too many steps)
369 //
370 if( recalculatedEndPt )
371 {
372 G4double endAchieved = IntersectPointVelct_G.GetCurveLength();
373 G4double endExpected = CurrentState.GetCurveLength();
374
375 // Detect failure - due to too many steps
376 G4bool shortEnd = endAchieved
377 < (endExpected*(1.0-CLHEP::perMillion));
378
379 G4double stepAchieved = endAchieved
380 - SubStepStartState.GetCurveLength();
381
382 // Update remaining state - must work for 'full' step or
383 // abandonned intersection
384 //
385 CurrentState = IntersectPointVelct_G;
386 s_length_taken = stepAchieved;
387 if( shortEnd )
388 {
389 fParticleIsLooping = true;
390 }
391 }
392 }
393 }
394 if( !intersects )
395 {
396 StepTaken += s_length_taken;
397
398 if (fpTrajectoryFilter) // For smooth trajectory display (jacek 1/11/2002)
399 {
400 fpTrajectoryFilter->TakeIntermediatePoint(CurrentState.GetPosition());
401 }
402 }
403 first_substep = false;
404
405#ifdef G4DEBUG_FIELD
407 {
409 G4cout << " Above 'Severe Action' threshold -- for Zero steps. ";
410 else
411 G4cout << " Above 'action' threshold -- for Zero steps. ";
412 G4cout << " Number of zero steps = " << fNoZeroStep << G4endl;
413 printStatus( SubStepStartState, // or OriginalState,
414 CurrentState, CurrentProposedStepLength,
415 NewSafety, do_loop_count, pPhysVol );
416 }
417 if( (fVerboseLevel > 1) && (do_loop_count > fMax_loop_count-10 ))
418 {
419 if( do_loop_count == fMax_loop_count-9 )
420 {
421 G4cout << " G4PropagatorInField::ComputeStep(): " << G4endl
422 << " Difficult track - taking many sub steps." << G4endl;
423 printStatus( SubStepStartState, SubStepStartState, CurrentProposedStepLength,
424 NewSafety, 0, pPhysVol );
425 }
426 printStatus( SubStepStartState, CurrentState, CurrentProposedStepLength,
427 NewSafety, do_loop_count, pPhysVol );
428 }
429#endif
430
431 ++do_loop_count;
432
433 } while( (!intersects )
435 && (StepTaken + kCarTolerance < CurrentProposedStepLength)
436 && ( do_loop_count < fMax_loop_count ) );
437
438 if( do_loop_count >= fMax_loop_count
439 && (StepTaken + kCarTolerance < CurrentProposedStepLength) )
440 {
441 fParticleIsLooping = true;
442 }
443 if ( ( fParticleIsLooping ) && (fVerboseLevel > 0) )
444 {
445 ReportLoopingParticle( do_loop_count, StepTaken,
446 CurrentProposedStepLength, methodName,
447 CurrentState.GetMomentum(), pPhysVol );
448 }
449
450 if( !intersects )
451 {
452 // Chord AB or "minor chords" do not intersect
453 // B is the endpoint Step of the current Step.
454 //
455 End_PointAndTangent = CurrentState;
456 TruePathLength = StepTaken; // Original code
457
458 // Tried the following to avoid potential issue with round-off error
459 // - but has issues... Suppressing this change JA 2015/05/02
460 // TruePathLength = CurrentProposedStepLength;
461 }
462 fLastStepInVolume = intersects;
463
464 // Set pFieldTrack to the return value
465 //
466 pFieldTrack = End_PointAndTangent;
467
468#ifdef G4VERBOSE
469 // Check that "s" is correct
470 //
471 if( std::fabs(OriginalState.GetCurveLength() + TruePathLength
472 - End_PointAndTangent.GetCurveLength()) > 3.e-4 * TruePathLength )
473 {
474 std::ostringstream message;
475 message << "Curve length mis-match between original state "
476 << "and proposed endpoint of propagation." << G4endl
477 << " The curve length of the endpoint should be: "
478 << OriginalState.GetCurveLength() + TruePathLength << G4endl
479 << " and it is instead: "
481 << " A difference of: "
482 << OriginalState.GetCurveLength() + TruePathLength
484 << " Original state = " << OriginalState << G4endl
485 << " Proposed state = " << End_PointAndTangent;
486 G4Exception(methodName, "GeomNav0003", FatalException, message);
487 }
488#endif
489
490 if( TruePathLength+kCarTolerance >= CurrentProposedStepLength )
491 {
492 fNoZeroStep = 0;
493 }
494 else
495 {
496 // In particular anomalous cases, we can get repeated zero steps
497 // We identify these cases and take corrective action when they occur.
498 //
499 if( TruePathLength < std::max( fZeroStepThreshold, 0.5*kCarTolerance ) )
500 {
501 ++fNoZeroStep;
502 }
503 else
504 {
505 fNoZeroStep = 0;
506 }
507 }
509 {
510 fParticleIsLooping = true;
511 ReportStuckParticle( fNoZeroStep, CurrentProposedStepLength,
513 fNoZeroStep = 0;
514 }
515
516 GetChordFinder()->SetDeltaChord(deltaChord);
517 return TruePathLength;
518}
G4double epsilon(G4double density, G4double temperature)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double cm
Definition: G4SIunits.hh:99
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
void OnComputeStep()
G4double GetDeltaChord() const
void SetDeltaChord(G4double newval)
G4double AdvanceChordLimited(G4FieldTrack &yCurrent, G4double stepInitial, G4double epsStep_Relative, const G4ThreeVector &latestSafetyOrigin, G4double lasestSafetyRadius)
G4double GetMinimumEpsilonStep() const
G4double GetDeltaOneStep() const
G4double GetCurveLength() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:614
G4VPhysicalVolume * GetWorldVolume() const
G4VCurvedTrajectoryFilter * fpTrajectoryFilter
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)
G4ChordFinder * GetChordFinder()
void ReportLoopingParticle(G4int count, G4double StepTaken, G4double stepRequest, const char *methodName, G4ThreeVector momentumVec, G4VPhysicalVolume *physVol)
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
void SetEpsilonStep(G4double newEps)
virtual void TakeIntermediatePoint(G4ThreeVector newPoint)=0
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
static const G4double kInfinity
Definition: geomdefs.hh:41
static constexpr double perMillion
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References G4ChordFinder::AdvanceChordLimited(), cm, G4VCurvedTrajectoryFilter::CreateNewTrajectorySegment(), End_PointAndTangent, epsilon(), fAbandonThreshold_NoZeroSteps, fActionThreshold_NoZeroSteps, FatalException, fCurrentFieldMgr, fEpsilonStep, fFirstStepInVolume, fFull_CurveLen_of_LastAttempt, fIncreaseChordDistanceThreshold, FindAndSetFieldManager(), fIntersectionLocator, fLargestAcceptableStep, fLast_ProposedStepLength, fLastStepInVolume, fMax_loop_count, fNavigator, fNewTrack, fNoZeroStep, fParticleIsLooping, fPreviousSafety, fPreviousSftOrigin, fpTrajectoryFilter, fSetFieldMgr, fSevereActionThreshold_NoZeroSteps, fVerboseLevel, fZeroStepThreshold, G4cerr, G4cout, G4endl, G4Exception(), GetChordFinder(), G4FieldTrack::GetCurveLength(), G4ChordFinder::GetDeltaChord(), G4FieldManager::GetDeltaOneStep(), G4VPhysicalVolume::GetLogicalVolume(), G4FieldManager::GetMaximumEpsilonStep(), G4FieldManager::GetMinimumEpsilonStep(), G4FieldTrack::GetMomentum(), G4FieldTrack::GetMomentumDir(), G4VPhysicalVolume::GetName(), G4FieldTrack::GetPosition(), G4Navigator::GetWorldVolume(), IntersectChord(), JustWarning, kCarTolerance, kInfinity, G4Navigator::LocateGlobalPointWithinVolume(), G4INCL::Math::max(), G4INCL::Math::min(), G4ChordFinder::OnComputeStep(), CLHEP::perMillion, printStatus(), PrintStepLengthDiagnostic(), RefreshIntersectionLocator(), ReportLoopingParticle(), ReportStuckParticle(), G4ChordFinder::SetDeltaChord(), SetEpsilonStep(), and G4VCurvedTrajectoryFilter::TakeIntermediatePoint().

Referenced by G4Transportation::AlongStepGetPhysicalInteractionLength(), and G4PathFinder::DoNextCurvedStep().

◆ EndMomentumDir()

G4ThreeVector G4PropagatorInField::EndMomentumDir ( ) const
inline

◆ EndPosition()

G4ThreeVector G4PropagatorInField::EndPosition ( ) const
inline

◆ FindAndSetFieldManager()

G4FieldManager * G4PropagatorInField::FindAndSetFieldManager ( G4VPhysicalVolume pCurrentPhysVol)

Definition at line 698 of file G4PropagatorInField.cc.

700{
701 G4FieldManager* currentFieldMgr;
702
703 currentFieldMgr = fDetectorFieldMgr;
704 if( pCurrentPhysicalVolume != nullptr )
705 {
706 G4FieldManager *pRegionFieldMgr = nullptr, *localFieldMgr = nullptr;
707 G4LogicalVolume* pLogicalVol = pCurrentPhysicalVolume->GetLogicalVolume();
708
709 if( pLogicalVol != nullptr )
710 {
711 // Value for Region, if any, overrides
712 //
713 G4Region* pRegion = pLogicalVol->GetRegion();
714 if( pRegion != nullptr )
715 {
716 pRegionFieldMgr = pRegion->GetFieldManager();
717 if( pRegionFieldMgr != nullptr )
718 {
719 currentFieldMgr= pRegionFieldMgr;
720 }
721 }
722
723 // 'Local' Value from logical volume, if any, overrides
724 //
725 localFieldMgr = pLogicalVol->GetFieldManager();
726 if ( localFieldMgr != nullptr )
727 {
728 currentFieldMgr = localFieldMgr;
729 }
730 }
731 }
732 fCurrentFieldMgr = currentFieldMgr;
733
734 // Flag that field manager has been set
735 //
736 fSetFieldMgr = true;
737
738 return currentFieldMgr;
739}
G4Region * GetRegion() const
G4FieldManager * GetFieldManager() const
G4FieldManager * GetFieldManager() const

References fCurrentFieldMgr, fDetectorFieldMgr, fSetFieldMgr, G4LogicalVolume::GetFieldManager(), G4Region::GetFieldManager(), G4VPhysicalVolume::GetLogicalVolume(), and G4LogicalVolume::GetRegion().

Referenced by G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4Transportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4PathFinder::ComputeStep(), ComputeStep(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4SynchrotronRadiation::PostStepDoIt(), and G4SynchrotronRadiationInMat::PostStepDoIt().

◆ GetChordFinder()

G4ChordFinder * G4PropagatorInField::GetChordFinder ( )
inline

◆ GetCurrentEquationOfMotion()

G4EquationOfMotion * G4PropagatorInField::GetCurrentEquationOfMotion ( )
inline

◆ GetCurrentFieldManager()

G4FieldManager * G4PropagatorInField::GetCurrentFieldManager ( )
inline

◆ GetDeltaIntersection()

G4double G4PropagatorInField::GetDeltaIntersection ( ) const
inline

◆ GetDeltaOneStep()

G4double G4PropagatorInField::GetDeltaOneStep ( ) const
inline

◆ GetEndState()

G4FieldTrack G4PropagatorInField::GetEndState ( ) const
inline

◆ GetEpsilonStep()

G4double G4PropagatorInField::GetEpsilonStep ( ) const
inline

◆ GetIntersectionLocator()

G4VIntersectionLocator * G4PropagatorInField::GetIntersectionLocator ( )
inline

◆ GetIterationsToIncreaseChordDistance()

G4int G4PropagatorInField::GetIterationsToIncreaseChordDistance ( ) const
inline

◆ GetLargestAcceptableStep()

G4double G4PropagatorInField::GetLargestAcceptableStep ( )
inline

◆ GetMaximumEpsilonStep()

G4double G4PropagatorInField::GetMaximumEpsilonStep ( ) const
inline

◆ GetMaxLoopCount()

G4int G4PropagatorInField::GetMaxLoopCount ( ) const
inline

◆ GetMinimumEpsilonStep()

G4double G4PropagatorInField::GetMinimumEpsilonStep ( ) const
inline

◆ GetNavigatorForPropagating()

G4Navigator * G4PropagatorInField::GetNavigatorForPropagating ( )
inline

◆ GetThresholdNoZeroSteps()

G4int G4PropagatorInField::GetThresholdNoZeroSteps ( G4int  i)
inline

◆ GetUseSafetyForOptimization()

G4bool G4PropagatorInField::GetUseSafetyForOptimization ( )
inline

◆ GetVerboseLevel()

G4int G4PropagatorInField::GetVerboseLevel ( ) const
inline

◆ GetVerboseTrace()

G4bool G4PropagatorInField::GetVerboseTrace ( )
inline

◆ GetZeroStepThreshold()

G4double G4PropagatorInField::GetZeroStepThreshold ( )
inline

◆ GimmeTrajectoryVectorAndForgetIt()

std::vector< G4ThreeVector > * G4PropagatorInField::GimmeTrajectoryVectorAndForgetIt ( ) const

Definition at line 649 of file G4PropagatorInField.cc.

650{
651 // NB, GimmeThePointsAndForgetThem really forgets them, so it can
652 // only be called (exactly) once for each step.
653
654 if (fpTrajectoryFilter != nullptr)
655 {
657 }
658 else
659 {
660 return nullptr;
661 }
662}
std::vector< G4ThreeVector > * GimmeThePointsAndForgetThem()

References fpTrajectoryFilter, and G4VCurvedTrajectoryFilter::GimmeThePointsAndForgetThem().

Referenced by G4ITTransportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepDoIt(), and G4Transportation::AlongStepDoIt().

◆ IntersectChord()

G4bool G4PropagatorInField::IntersectChord ( const G4ThreeVector StartPointA,
const G4ThreeVector EndPointB,
G4double NewSafety,
G4double LinearStepLength,
G4ThreeVector IntersectionPoint 
)
inline

Referenced by ComputeStep().

◆ IsFirstStepInVolume()

G4bool G4PropagatorInField::IsFirstStepInVolume ( )
inline

◆ IsLastStepInVolume()

G4bool G4PropagatorInField::IsLastStepInVolume ( )
inline

◆ IsParticleLooping()

G4bool G4PropagatorInField::IsParticleLooping ( ) const
inline

◆ PrepareNewTrack()

void G4PropagatorInField::PrepareNewTrack ( )
inline

◆ printStatus()

void G4PropagatorInField::printStatus ( const G4FieldTrack startFT,
const G4FieldTrack currentFT,
G4double  requestStep,
G4double  safety,
G4int  step,
G4VPhysicalVolume startVolume 
)

Definition at line 524 of file G4PropagatorInField.cc.

530{
531 const G4int verboseLevel = fVerboseLevel;
532 const G4ThreeVector StartPosition = StartFT.GetPosition();
533 const G4ThreeVector StartUnitVelocity = StartFT.GetMomentumDir();
534 const G4ThreeVector CurrentPosition = CurrentFT.GetPosition();
535 const G4ThreeVector CurrentUnitVelocity = CurrentFT.GetMomentumDir();
536
537 G4double step_len = CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
538
539 G4int oldprec; // cout/cerr precision settings
540
541 if( ((stepNo == 0) && (verboseLevel <3)) || (verboseLevel >= 3) )
542 {
543 oldprec = G4cout.precision(4);
544 G4cout << std::setw( 5) << "Step#"
545 << std::setw(10) << " s " << " "
546 << std::setw(10) << "X(mm)" << " "
547 << std::setw(10) << "Y(mm)" << " "
548 << std::setw(10) << "Z(mm)" << " "
549 << std::setw( 7) << " N_x " << " "
550 << std::setw( 7) << " N_y " << " "
551 << std::setw( 7) << " N_z " << " " ;
552 G4cout << std::setw( 7) << " Delta|N|" << " "
553 << std::setw( 9) << "StepLen" << " "
554 << std::setw(12) << "StartSafety" << " "
555 << std::setw( 9) << "PhsStep" << " ";
556 if( startVolume != nullptr )
557 { G4cout << std::setw(18) << "NextVolume" << " "; }
558 G4cout.precision(oldprec);
559 G4cout << G4endl;
560 }
561 if((stepNo == 0) && (verboseLevel <=3))
562 {
563 // Recurse to print the start values
564 //
565 printStatus( StartFT, StartFT, -1.0, safety, -1, startVolume);
566 }
567 if( verboseLevel <= 3 )
568 {
569 if( stepNo >= 0)
570 { G4cout << std::setw( 4) << stepNo << " "; }
571 else
572 { G4cout << std::setw( 5) << "Start" ; }
573 oldprec = G4cout.precision(8);
574 G4cout << std::setw(10) << CurrentFT.GetCurveLength() << " ";
575 G4cout.precision(8);
576 G4cout << std::setw(10) << CurrentPosition.x() << " "
577 << std::setw(10) << CurrentPosition.y() << " "
578 << std::setw(10) << CurrentPosition.z() << " ";
579 G4cout.precision(4);
580 G4cout << std::setw( 7) << CurrentUnitVelocity.x() << " "
581 << std::setw( 7) << CurrentUnitVelocity.y() << " "
582 << std::setw( 7) << CurrentUnitVelocity.z() << " ";
583 G4cout.precision(3);
584 G4cout << std::setw( 7)
585 << CurrentFT.GetMomentum().mag()-StartFT.GetMomentum().mag() << " ";
586 G4cout << std::setw( 9) << step_len << " ";
587 G4cout << std::setw(12) << safety << " ";
588 if( requestStep != -1.0 )
589 { G4cout << std::setw( 9) << requestStep << " "; }
590 else
591 { G4cout << std::setw( 9) << "Init/NotKnown" << " "; }
592 if( startVolume != 0)
593 { G4cout << std::setw(12) << startVolume->GetName() << " "; }
594 G4cout.precision(oldprec);
595 G4cout << G4endl;
596 }
597 else // if( verboseLevel > 3 )
598 {
599 // Multi-line output
600
601 G4cout << "Step taken was " << step_len
602 << " out of PhysicalStep = " << requestStep << G4endl;
603 G4cout << "Final safety is: " << safety << G4endl;
604 G4cout << "Chord length = " << (CurrentPosition-StartPosition).mag()
605 << G4endl;
606 G4cout << G4endl;
607 }
608}
double z() const
double x() const
double y() const

References fVerboseLevel, G4cout, G4endl, G4FieldTrack::GetCurveLength(), G4FieldTrack::GetMomentum(), G4FieldTrack::GetMomentumDir(), G4VPhysicalVolume::GetName(), G4FieldTrack::GetPosition(), CLHEP::Hep3Vector::mag(), printStatus(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by ComputeStep(), and printStatus().

◆ PrintStepLengthDiagnostic()

void G4PropagatorInField::PrintStepLengthDiagnostic ( G4double  currentProposedStepLength,
G4double  decreaseFactor,
G4double  stepTrial,
const G4FieldTrack aFieldTrack 
)
protected

Definition at line 614 of file G4PropagatorInField.cc.

619{
620 G4int iprec= G4cout.precision(8);
621 G4cout << " " << std::setw(12) << " PiF: NoZeroStep "
622 << " " << std::setw(20) << " CurrentProposed len "
623 << " " << std::setw(18) << " Full_curvelen_last"
624 << " " << std::setw(18) << " last proposed len "
625 << " " << std::setw(18) << " decrease factor "
626 << " " << std::setw(15) << " step trial "
627 << G4endl;
628
629 G4cout << " " << std::setw(10) << fNoZeroStep << " "
630 << " " << std::setw(20) << CurrentProposedStepLength
631 << " " << std::setw(18) << fFull_CurveLen_of_LastAttempt
632 << " " << std::setw(18) << fLast_ProposedStepLength
633 << " " << std::setw(18) << decreaseFactor
634 << " " << std::setw(15) << stepTrial
635 << G4endl;
636 G4cout.precision( iprec );
637}

References fFull_CurveLen_of_LastAttempt, fLast_ProposedStepLength, fNoZeroStep, G4cout, and G4endl.

Referenced by ComputeStep().

◆ RefreshIntersectionLocator()

void G4PropagatorInField::RefreshIntersectionLocator ( )

◆ ReportLoopingParticle()

void G4PropagatorInField::ReportLoopingParticle ( G4int  count,
G4double  StepTaken,
G4double  stepRequest,
const char *  methodName,
G4ThreeVector  momentumVec,
G4VPhysicalVolume physVol 
)
protected

Definition at line 760 of file G4PropagatorInField.cc.

766{
767 std::ostringstream message;
768 G4double fraction = StepTaken / StepRequested;
769 message << " Unfinished integration of track (likely looping particle) "
770 << " of momentum " << momentumVec << " ( magnitude = "
771 << momentumVec.mag() << " ) " << G4endl
772 << " after " << count << " field substeps "
773 << " totaling " << std::setprecision(12) << StepTaken / mm << " mm "
774 << " out of requested step " << std::setprecision(12)
775 << StepRequested / mm << " mm ";
776 message << " a fraction of ";
777 G4int prec = 4;
778 if( fraction > 0.99 )
779 {
780 prec = 7;
781 }
782 else
783 {
784 if (fraction > 0.97 ) { prec = 5; }
785 }
786 message << std::setprecision(prec)
787 << 100. * StepTaken / StepRequested << " % " << G4endl ;
788 if( pPhysVol )
789 {
790 message << " in volume " << pPhysVol->GetName() ;
791 auto material = pPhysVol->GetLogicalVolume()->GetMaterial();
792 if( material != nullptr )
793 message << " with material " << material->GetName()
794 << " ( density = "
795 << material->GetDensity() / ( g/(cm*cm*cm) ) << " g / cm^3 ) ";
796 }
797 else
798 {
799 message << " in unknown (null) volume. " ;
800 }
801 G4Exception(methodName, "GeomNav1002", JustWarning, message);
802}
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double g
Definition: G4SIunits.hh:168
double mag() const
static const double prec
Definition: RanecuEngine.cc:61
string material
Definition: eplot.py:19

References cm, g, G4endl, G4Exception(), G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetMaterial(), G4VPhysicalVolume::GetName(), JustWarning, CLHEP::Hep3Vector::mag(), eplot::material, mm, and CLHEP::prec.

Referenced by ComputeStep().

◆ ReportStuckParticle()

void G4PropagatorInField::ReportStuckParticle ( G4int  noZeroSteps,
G4double  proposedStep,
G4double  lastTriedStep,
G4VPhysicalVolume physVol 
)
protected

Definition at line 806 of file G4PropagatorInField.cc.

810{
811 std::ostringstream message;
812 message << "Particle is stuck; it will be killed." << G4endl
813 << " Zero progress for " << noZeroSteps << " attempted steps."
814 << G4endl
815 << " Proposed Step is " << proposedStep
816 << " but Step Taken is "<< lastTriedStep << G4endl;
817 if( physVol != nullptr )
818 message << " in volume " << physVol->GetName() ;
819 else
820 message << " in unknown or null volume. " ;
821 G4Exception("G4PropagatorInField::ComputeStep()",
822 "GeomNav1002", JustWarning, message);
823}

References G4endl, G4Exception(), G4VPhysicalVolume::GetName(), and JustWarning.

Referenced by ComputeStep().

◆ SetDetectorFieldManager()

void G4PropagatorInField::SetDetectorFieldManager ( G4FieldManager newGlobalFieldManager)
inline

◆ SetEpsilonStep()

void G4PropagatorInField::SetEpsilonStep ( G4double  newEps)
inline

Referenced by ComputeStep().

◆ SetIntersectionLocator()

void G4PropagatorInField::SetIntersectionLocator ( G4VIntersectionLocator pLocator)
inline

◆ SetIterationsToIncreaseChordDistance()

void G4PropagatorInField::SetIterationsToIncreaseChordDistance ( G4int  numIters)
inline

◆ SetLargestAcceptableStep()

void G4PropagatorInField::SetLargestAcceptableStep ( G4double  newBigDist)
inline

◆ SetMaximumEpsilonStep()

void G4PropagatorInField::SetMaximumEpsilonStep ( G4double  newEpsMax)
inline

◆ SetMaxLoopCount()

void G4PropagatorInField::SetMaxLoopCount ( G4int  new_max)
inline

◆ SetMinimumEpsilonStep()

void G4PropagatorInField::SetMinimumEpsilonStep ( G4double  newEpsMin)
inline

◆ SetNavigatorForPropagating()

void G4PropagatorInField::SetNavigatorForPropagating ( G4Navigator SimpleOrMultiNavigator)
inline

◆ SetThresholdNoZeroStep()

void G4PropagatorInField::SetThresholdNoZeroStep ( G4int  noAct,
G4int  noHarsh,
G4int  noAbandon 
)
inline

◆ SetTrajectoryFilter()

void G4PropagatorInField::SetTrajectoryFilter ( G4VCurvedTrajectoryFilter filter)

Definition at line 667 of file G4PropagatorInField.cc.

668{
669 fpTrajectoryFilter = filter;
670}

References fpTrajectoryFilter.

Referenced by G4TrackingMessenger::SetNewValue().

◆ SetUseSafetyForOptimization()

void G4PropagatorInField::SetUseSafetyForOptimization ( G4bool  )
inline

◆ SetVerboseLevel()

G4int G4PropagatorInField::SetVerboseLevel ( G4int  verbose)

Definition at line 743 of file G4PropagatorInField.cc.

744{
745 G4int oldval = fVerboseLevel;
746 fVerboseLevel = level;
747
748 // Forward the verbose level 'reduced' to ChordFinder,
749 // MagIntegratorDriver ... ?
750 //
751 auto integrDriver = GetChordFinder()->GetIntegrationDriver();
752 integrDriver->SetVerboseLevel( fVerboseLevel - 2 );
753 G4cout << "Set Driver verbosity to " << fVerboseLevel - 2 << G4endl;
754
755 return oldval;
756}
G4VIntegrationDriver * GetIntegrationDriver()
virtual void SetVerboseLevel(G4int level)=0

References fVerboseLevel, G4cout, G4endl, GetChordFinder(), G4ChordFinder::GetIntegrationDriver(), and G4VIntegrationDriver::SetVerboseLevel().

◆ SetVerboseTrace()

void G4PropagatorInField::SetVerboseTrace ( G4bool  enable)
inline

◆ SetZeroStepThreshold()

void G4PropagatorInField::SetZeroStepThreshold ( G4double  newLength)
inline

◆ Verbose()

G4int G4PropagatorInField::Verbose ( ) const
inline

Field Documentation

◆ End_PointAndTangent

G4FieldTrack G4PropagatorInField::End_PointAndTangent
private

Definition at line 283 of file G4PropagatorInField.hh.

Referenced by ClearPropagatorState(), and ComputeStep().

◆ fAbandonThreshold_NoZeroSteps

G4int G4PropagatorInField::fAbandonThreshold_NoZeroSteps = 50
private

Definition at line 238 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fActionThreshold_NoZeroSteps

G4int G4PropagatorInField::fActionThreshold_NoZeroSteps = 2
private

Definition at line 236 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fAllocatedLocator

G4bool G4PropagatorInField::fAllocatedLocator
private

Definition at line 249 of file G4PropagatorInField.hh.

Referenced by G4PropagatorInField(), and ~G4PropagatorInField().

◆ fCheck

G4bool G4PropagatorInField::fCheck = false
private

Definition at line 297 of file G4PropagatorInField.hh.

◆ fCurrentFieldMgr

G4FieldManager* G4PropagatorInField::fCurrentFieldMgr
private

◆ fDetectorFieldMgr

G4FieldManager* G4PropagatorInField::fDetectorFieldMgr
private

Definition at line 254 of file G4PropagatorInField.hh.

Referenced by FindAndSetFieldManager(), and G4PropagatorInField().

◆ fEpsilonStep

G4double G4PropagatorInField::fEpsilonStep
private

◆ fFirstStepInVolume

G4bool G4PropagatorInField::fFirstStepInVolume = true
private

Definition at line 300 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fFull_CurveLen_of_LastAttempt

G4double G4PropagatorInField::fFull_CurveLen_of_LastAttempt = -1
private

◆ fIncreaseChordDistanceThreshold

G4int G4PropagatorInField::fIncreaseChordDistanceThreshold = 100
private

Definition at line 230 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fIntersectionLocator

G4VIntersectionLocator* G4PropagatorInField::fIntersectionLocator
private

◆ fLargestAcceptableStep

G4double G4PropagatorInField::fLargestAcceptableStep
private

Definition at line 242 of file G4PropagatorInField.hh.

Referenced by ComputeStep(), and G4PropagatorInField().

◆ fLast_ProposedStepLength

G4double G4PropagatorInField::fLast_ProposedStepLength = -1
private

◆ fLastStepInVolume

G4bool G4PropagatorInField::fLastStepInVolume = true
private

Definition at line 301 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fMax_loop_count

G4int G4PropagatorInField::fMax_loop_count = 1000
private

Definition at line 228 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fNavigator

G4Navigator* G4PropagatorInField::fNavigator
private

Definition at line 267 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fNewTrack

G4bool G4PropagatorInField::fNewTrack = true
private

Definition at line 302 of file G4PropagatorInField.hh.

Referenced by ClearPropagatorState(), and ComputeStep().

◆ fNoZeroStep

G4int G4PropagatorInField::fNoZeroStep = 0
private

◆ fParticleIsLooping

G4bool G4PropagatorInField::fParticleIsLooping = false
private

Definition at line 284 of file G4PropagatorInField.hh.

Referenced by ClearPropagatorState(), and ComputeStep().

◆ fPreviousSafety

G4double G4PropagatorInField::fPreviousSafety = 0.0
private

Definition at line 292 of file G4PropagatorInField.hh.

Referenced by ClearPropagatorState(), and ComputeStep().

◆ fPreviousSftOrigin

G4ThreeVector G4PropagatorInField::fPreviousSftOrigin
private

Definition at line 291 of file G4PropagatorInField.hh.

Referenced by ClearPropagatorState(), ComputeStep(), and G4PropagatorInField().

◆ fpTrajectoryFilter

G4VCurvedTrajectoryFilter* G4PropagatorInField::fpTrajectoryFilter = nullptr
private

◆ fSetFieldMgr

G4bool G4PropagatorInField::fSetFieldMgr = false
private

◆ fSevereActionThreshold_NoZeroSteps

G4int G4PropagatorInField::fSevereActionThreshold_NoZeroSteps = 10
private

Definition at line 237 of file G4PropagatorInField.hh.

Referenced by ComputeStep().

◆ fUseSafetyForOptimisation

G4bool G4PropagatorInField::fUseSafetyForOptimisation = true
private

Definition at line 231 of file G4PropagatorInField.hh.

Referenced by RefreshIntersectionLocator().

◆ fVerboseLevel

G4int G4PropagatorInField::fVerboseLevel = 0
private

◆ fVerbTracePiF

G4bool G4PropagatorInField::fVerbTracePiF = false
private

Definition at line 296 of file G4PropagatorInField.hh.

Referenced by G4PropagatorInField().

◆ fZeroStepThreshold

G4double G4PropagatorInField::fZeroStepThreshold = 0.0
private

Definition at line 239 of file G4PropagatorInField.hh.

Referenced by ComputeStep(), and G4PropagatorInField().

◆ kCarTolerance

G4double G4PropagatorInField::kCarTolerance
private

Definition at line 246 of file G4PropagatorInField.hh.

Referenced by ComputeStep(), and G4PropagatorInField().


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