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

G4PropagatorInField performs the navigation/propagation of a particle/track in a magnetic field. The field is in general non-uniform. For the calculation of the path, it relies on the class G4ChordFinder. It utilises an ODE solver (with the Runge-Kutta method) to evolve the particle, and drives it until the particle has traveled a set distance or it enters a new volume. More...

#include <G4PropagatorInField.hh>

Public Member Functions

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

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, const G4ThreeVector &momentumVec, G4VPhysicalVolume *physVol)
void ReportStuckParticle (G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)

Detailed Description

G4PropagatorInField performs the navigation/propagation of a particle/track in a magnetic field. The field is in general non-uniform. For the calculation of the path, it relies on the class G4ChordFinder. It utilises an ODE solver (with the Runge-Kutta method) to evolve the particle, and drives it until the particle has traveled a set distance or it enters a new volume.

Definition at line 65 of file G4PropagatorInField.hh.

Constructor & Destructor Documentation

◆ G4PropagatorInField()

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

Constructor and Destructor.

Definition at line 49 of file G4PropagatorInField.cc.

52 : fDetectorFieldMgr(detectorFieldMgr),
53 fNavigator(theNavigator),
54 fCurrentFieldMgr(detectorFieldMgr),
55 End_PointAndTangent(G4ThreeVector(0.,0.,0.),
56 G4ThreeVector(0.,0.,0.),0.0,0.0,0.0,0.0,0.0)
57{
58 fEpsilonStep = (fDetectorFieldMgr != nullptr)
59 ? fDetectorFieldMgr->GetMaximumEpsilonStep() : 1.0e-5;
60
61
62 fPreviousSftOrigin = G4ThreeVector(0.,0.,0.);
64 fZeroStepThreshold = std::max( 1.0e5 * kCarTolerance, 1.0e-1 * micrometer );
65
66 fLargestAcceptableStep = 100.0 * meter; // Reduced from 1000.0 * meter
67 fMaxStepSizeMultiplier= 0.1 ; // 0.1 in git (larger for tests.) // Reduced from 100;
68 fMinBigDistance= 100. * CLHEP::mm;
69#ifdef G4DEBUG_FIELD
70 G4cout << " PiF: Zero Step Threshold set to "
71 << fZeroStepThreshold / millimeter
72 << " mm." << G4endl;
73 G4cout << " PiF: Value of kCarTolerance = "
74 << kCarTolerance / millimeter
75 << " mm. " << G4endl;
76 fVerboseLevel = 2;
77 fVerbTracePiF = true;
78#endif
79
80 // Defining Intersection Locator and his parameters
81 if ( vLocator == nullptr )
82 {
83 fIntersectionLocator = new G4MultiLevelLocator(theNavigator);
84 fAllocatedLocator = true;
85 }
86 else
87 {
88 fIntersectionLocator = vLocator;
89 fAllocatedLocator = false;
90 }
91 RefreshIntersectionLocator(); // Copy all relevant parameters
92}
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()

◆ ~G4PropagatorInField()

G4PropagatorInField::~G4PropagatorInField ( )

Definition at line 96 of file G4PropagatorInField.cc.

97{
98 if(fAllocatedLocator) { delete fIntersectionLocator; }
99}

Member Function Documentation

◆ CheckMode()

void G4PropagatorInField::CheckMode ( G4bool mode)
inline

Enabling check mode for further diagnostics.

◆ ClearPropagatorState()

void G4PropagatorInField::ClearPropagatorState ( )

Clears the State of this class and its current associates.

Note
The current field manager & chord finder will also be called.

Definition at line 673 of file G4PropagatorInField.cc.

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

◆ ComputeStep()

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

Computes the next geometric Step.

Parameters
[in,out]pFieldTrackField track to be filled.
[in]pCurrentProposedStepLengthCurrent proposed step length.
[in,out]pNewSafetyNew safety.
[in]pPhysVolPointer to the current volume.
[in]canRelaxDeltaChordTo enable relaxing delta-chord parameter.
Returns
Step length.

Definition at line 115 of file G4PropagatorInField.cc.

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

◆ EndMomentumDir()

G4ThreeVector G4PropagatorInField::EndMomentumDir ( ) const
inline

◆ EndPosition()

G4ThreeVector G4PropagatorInField::EndPosition ( ) const
inline

Returning the state after the Step.

◆ FindAndSetFieldManager()

G4FieldManager * G4PropagatorInField::FindAndSetFieldManager ( G4VPhysicalVolume * pCurrentPhysVol)

Sets (and returns) the correct field manager (global or local), if it exists.

Note
Should be called before ComputeStep is called; Currently, ComputeStep() will call it, if it has not been called.
Parameters
[in]pCurrentPhysVolPointer to the current volume.
Returns
The pointer to the field manager.

Definition at line 697 of file G4PropagatorInField.cc.

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

Referenced by ComputeStep().

◆ GetChordFinder()

G4ChordFinder * G4PropagatorInField::GetChordFinder ( )
inline

Returning the pointer to the chord finder.

Referenced by ComputeStep(), RefreshIntersectionLocator(), and SetVerboseLevel().

◆ GetCurrentEquationOfMotion()

G4EquationOfMotion * G4PropagatorInField::GetCurrentEquationOfMotion ( )
inline

◆ GetCurrentFieldManager()

G4FieldManager * G4PropagatorInField::GetCurrentFieldManager ( )
inline

Auxiliary methods.

Note
Their results can/will change during propagation.

Referenced by G4DecayWithSpin::AtRestDoIt().

◆ GetDeltaIntersection()

G4double G4PropagatorInField::GetDeltaIntersection ( ) const
inline

Accessors.

◆ GetDeltaOneStep()

G4double G4PropagatorInField::GetDeltaOneStep ( ) const
inline

◆ GetEndState()

G4FieldTrack G4PropagatorInField::GetEndState ( ) const
inline

Accessor for retrieving the field track.

◆ GetEpsilonStep()

G4double G4PropagatorInField::GetEpsilonStep ( ) const
inline

Returning the relative accuracy for the current Step.

◆ GetIntersectionLocator()

G4VIntersectionLocator * G4PropagatorInField::GetIntersectionLocator ( )
inline

Changes or gets the object which calculates the exact intersection point with the next boundary.

◆ GetIterationsToIncreaseChordDistance()

G4int G4PropagatorInField::GetIterationsToIncreaseChordDistance ( ) const
inline

Controls the parameter which enables the temporary 'relaxation' which ensures that chord segments are short enough so that their sagitta is small than delta-chord parameter. The Set method increases the value of delta-chord temporarily, doubling it once the number of iterations substeps reach value of 'IncreaseChordDistanceThreshold'. It is also doubled again every time the iteration count reaches a multiple of this value.

Note
The delta-chord is reset to its original value at the end of each call to ComputeStep().

◆ GetLargestAcceptableStep()

G4double G4PropagatorInField::GetLargestAcceptableStep ( )

Definition at line 837 of file G4PropagatorInField.cc.

838{
839 return fLargestAcceptableStep;
840}

◆ GetMaximumEpsilonStep()

G4double G4PropagatorInField::GetMaximumEpsilonStep ( ) const
inline

◆ GetMaxLoopCount()

G4int G4PropagatorInField::GetMaxLoopCount ( ) const
inline

Accessor/modifier for controlling the maximum for the number of substeps that a particle can take. Above this number it is signaled as 'looping'.

◆ GetMaxStepSizeMultiplier()

G4double G4PropagatorInField::GetMaxStepSizeMultiplier ( )

Methods to control extra Multiplier parameter for limiting long steps.

Definition at line 854 of file G4PropagatorInField.cc.

855{
856 return fMaxStepSizeMultiplier;
857}

◆ GetMinBigDistance()

G4double G4PropagatorInField::GetMinBigDistance ( )

Methods to Control minimum 'directional' distance in case of too-large step.

Definition at line 868 of file G4PropagatorInField.cc.

869{
870 return fMinBigDistance;
871}

◆ GetMinimumEpsilonStep()

G4double G4PropagatorInField::GetMinimumEpsilonStep ( ) const
inline

Methods to control values for the global field manager.

Deprecated
The four methods below are now obsolescent but for now will work. They are being replaced by same-name methods in G4FieldManager, allowing the specialisation in different volumes.

◆ 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

Accesses the points which have passed by the filter.

Note
Responsibility for deleting the points lies with the client. This method MUST BE called exactly ONCE per step.

Definition at line 651 of file G4PropagatorInField.cc.

652{
653 // NB, GimmeThePointsAndForgetThem really forgets them, so it can
654 // only be called (exactly) once for each step.
655
656 if (fpTrajectoryFilter != nullptr)
657 {
658 return fpTrajectoryFilter->GimmeThePointsAndForgetThem();
659 }
660 return nullptr;
661}

◆ IntersectChord()

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

Intersects the chord from StartPointA to EndPointB and returns whether an intersection occurred.

Note
Safety is changed!

Referenced by ComputeStep().

◆ IsFirstStepInVolume()

G4bool G4PropagatorInField::IsFirstStepInVolume ( )
inline

Returns if it is the first step in the volume.

◆ IsLastStepInVolume()

G4bool G4PropagatorInField::IsLastStepInVolume ( )
inline

Returns if it is the last step in the volume.

◆ IsParticleLooping()

G4bool G4PropagatorInField::IsParticleLooping ( ) const
inline

◆ PrepareNewTrack()

void G4PropagatorInField::PrepareNewTrack ( )
inline

Initialises track flags.

◆ printStatus()

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

Print method, useful mostly for debugging.

Definition at line 526 of file G4PropagatorInField.cc.

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

Referenced by ComputeStep(), and printStatus().

◆ PrintStepLengthDiagnostic()

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

Logging methods.

Definition at line 616 of file G4PropagatorInField.cc.

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

Referenced by ComputeStep().

◆ RefreshIntersectionLocator()

void G4PropagatorInField::RefreshIntersectionLocator ( )

Updates the Locator with parameters from this class and from current field manager.

Definition at line 104 of file G4PropagatorInField.cc.

105{
106 fIntersectionLocator->SetEpsilonStepFor(fEpsilonStep);
107 fIntersectionLocator->SetDeltaIntersectionFor(fCurrentFieldMgr->GetDeltaIntersection());
108 fIntersectionLocator->SetChordFinderFor(GetChordFinder());
109 fIntersectionLocator->SetSafetyParametersFor( fUseSafetyForOptimisation);
110}

Referenced by ComputeStep(), and G4PropagatorInField().

◆ ReportLoopingParticle()

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

Definition at line 759 of file G4PropagatorInField.cc.

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

Referenced by ComputeStep().

◆ ReportStuckParticle()

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

Definition at line 807 of file G4PropagatorInField.cc.

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

Referenced by ComputeStep().

◆ ResetLargestAcceptableStep()

void G4PropagatorInField::ResetLargestAcceptableStep ( )

◆ SetDetectorFieldManager()

void G4PropagatorInField::SetDetectorFieldManager ( G4FieldManager * newGlobalFieldManager)
inline

Setter for global field manager. Updates the state.

◆ SetEpsilonStep()

void G4PropagatorInField::SetEpsilonStep ( G4double newEps)
inline

Setting the relative accuracy for the current Step. The ratio DeltaOneStep()/h_current_step.

Referenced by ComputeStep().

◆ SetIntersectionLocator()

void G4PropagatorInField::SetIntersectionLocator ( G4VIntersectionLocator * pLocator)
inline

◆ SetIterationsToIncreaseChordDistance()

void G4PropagatorInField::SetIterationsToIncreaseChordDistance ( G4int numIters)
inline

◆ SetLargestAcceptableStep()

void G4PropagatorInField::SetLargestAcceptableStep ( G4double newBigDist)

Methods to obtain / change the size of the largest step the method will undertake. The Reset method uses the world volume's.

Definition at line 844 of file G4PropagatorInField.cc.

845{
846 if( fLargestAcceptableStep>0.0 )
847 {
848 fLargestAcceptableStep = newBigDist;
849 }
850}

◆ SetMaximumEpsilonStep()

void G4PropagatorInField::SetMaximumEpsilonStep ( G4double newEpsMax)
inline

◆ SetMaxLoopCount()

void G4PropagatorInField::SetMaxLoopCount ( G4int new_max)
inline

◆ SetMaxStepSizeMultiplier()

void G4PropagatorInField::SetMaxStepSizeMultiplier ( G4double vm)

Definition at line 861 of file G4PropagatorInField.cc.

862{
863 fMaxStepSizeMultiplier=vm;
864}

◆ SetMinBigDistance()

void G4PropagatorInField::SetMinBigDistance ( G4double val)

Definition at line 875 of file G4PropagatorInField.cc.

876{
877 fMinBigDistance= val;
878}

◆ SetMinimumEpsilonStep()

void G4PropagatorInField::SetMinimumEpsilonStep ( G4double newEpsMin)
inline

◆ SetNavigatorForPropagating()

void G4PropagatorInField::SetNavigatorForPropagating ( G4Navigator * SimpleOrMultiNavigator)
inline

Accessor and modifier for navigator.

◆ SetThresholdNoZeroStep()

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

Accessors and modifiers for no-zero steps threshold.

◆ SetTrajectoryFilter()

void G4PropagatorInField::SetTrajectoryFilter ( G4VCurvedTrajectoryFilter * filter)

Sets the filter that examines & stores 'intermediate' curved trajectory points.

Note
Currently only position is stored.

Definition at line 666 of file G4PropagatorInField.cc.

667{
668 fpTrajectoryFilter = filter;
669}

Referenced by G4TrackingMessenger::SetNewValue().

◆ SetUseSafetyForOptimization()

void G4PropagatorInField::SetUseSafetyForOptimization ( G4bool )
inline

Toggles & views parameter for using safety to discard unneccesary calls to the navigator (thus 'optimising' performance).

◆ SetVerboseLevel()

G4int G4PropagatorInField::SetVerboseLevel ( G4int verbose)

Verbosity control.

Definition at line 742 of file G4PropagatorInField.cc.

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

◆ SetVerboseTrace()

void G4PropagatorInField::SetVerboseTrace ( G4bool enable)
inline

Accessor/modifier for tracing key parts of ComputeStep().

◆ SetZeroStepThreshold()

void G4PropagatorInField::SetZeroStepThreshold ( G4double newLength)
inline

◆ Verbose()

G4int G4PropagatorInField::Verbose ( ) const
inline

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