Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PropagatorInField.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4PropagatorInField
27//
28// Class description:
29//
30// This class performs the navigation/propagation of a particle/track
31// in a magnetic field. The field is in general non-uniform.
32// For the calculation of the path, it relies on the class G4ChordFinder.
33// It utilises an ODE solver (with the Runge-Kutta method) to evolve the
34// particle, and drives it until the particle has traveled a set distance
35// or it enters a new volume.
36
37// Author: John Apostolakis (CERN), 25 October 1996
38// ---------------------------------------------------------------------------
39#ifndef G4PropagatorInField_hh
40#define G4PropagatorInField_hh 1
41
42#include "G4Types.hh"
43
44#include <vector>
45
46#include "G4FieldTrack.hh"
47#include "G4FieldManager.hh"
49
50class G4ChordFinder;
51
52class G4Navigator;
55
56/**
57 * @brief G4PropagatorInField performs the navigation/propagation of a
58 * particle/track in a magnetic field. The field is in general non-uniform.
59 * For the calculation of the path, it relies on the class G4ChordFinder.
60 * It utilises an ODE solver (with the Runge-Kutta method) to evolve the
61 * particle, and drives it until the particle has traveled a set distance
62 * or it enters a new volume.
63 */
64
66{
67 public:
68
69 /**
70 * Constructor and Destructor.
71 */
72 G4PropagatorInField( G4Navigator* theNavigator,
73 G4FieldManager* detectorFieldMgr,
74 G4VIntersectionLocator* vLocator = nullptr );
76
77 /**
78 * Computes the next geometric Step.
79 * @param[in,out] pFieldTrack Field track to be filled.
80 * @param[in] pCurrentProposedStepLength Current proposed step length.
81 * @param[in,out] pNewSafety New safety.
82 * @param[in] pPhysVol Pointer to the current volume.
83 * @param[in] canRelaxDeltaChord To enable relaxing delta-chord parameter.
84 * @returns Step length.
85 */
86 G4double ComputeStep( G4FieldTrack& pFieldTrack,
87 G4double pCurrentProposedStepLength,
88 G4double& pNewSafety,
89 G4VPhysicalVolume* pPhysVol = nullptr,
90 G4bool canRelaxDeltaChord = false);
91
92 /**
93 * Returning the state after the Step.
94 */
95 inline G4ThreeVector EndPosition() const;
97 inline G4bool IsParticleLooping() const;
98
99 /**
100 * Returning the relative accuracy for the current Step.
101 */
102 inline G4double GetEpsilonStep() const;
103
104 /**
105 * Setting the relative accuracy for the current Step.
106 * The ratio DeltaOneStep()/h_current_step.
107 */
108 inline void SetEpsilonStep(G4double newEps);
109
110 /**
111 * Sets (and returns) the correct field manager (global or local),
112 * if it exists.
113 * @note Should be called before ComputeStep is called;
114 * Currently, ComputeStep() will call it, if it has not been called.
115 * @param[in] pCurrentPhysVol Pointer to the current volume.
116 * @returns The pointer to the field manager.
117 */
119
120 /**
121 * Returning the pointer to the chord finder.
122 */
124
125 /**
126 * Verbosity control.
127 */
128 G4int SetVerboseLevel( G4int verbose );
129 inline G4int GetVerboseLevel() const;
130 inline G4int Verbose() const;
131
132 /**
133 * Enabling check mode for further diagnostics.
134 */
135 inline void CheckMode(G4bool mode);
136
137 /**
138 * Accessor/modifier for tracing key parts of ComputeStep().
139 */
140 inline void SetVerboseTrace( G4bool enable );
142
143 /**
144 * Accessor/modifier for controlling the maximum for the number of
145 * substeps that a particle can take. Above this number it is signaled
146 * as 'looping'.
147 */
148 inline G4int GetMaxLoopCount() const;
149 inline void SetMaxLoopCount( G4int new_max );
150
151 /**
152 * Print method, useful mostly for debugging.
153 */
154 void printStatus( const G4FieldTrack& startFT,
155 const G4FieldTrack& currentFT,
156 G4double requestStep,
157 G4double safety,
158 G4int step,
159 G4VPhysicalVolume* startVolume);
160
161 /**
162 * Accessor for retrieving the field track.
163 */
164 inline G4FieldTrack GetEndState() const;
165
166 /**
167 * Methods to control values for the global field manager.
168 * @deprecated The four methods below are now obsolescent but *for now*
169 * will work. They are being replaced by same-name methods in
170 * G4FieldManager, allowing the specialisation in different volumes.
171 */
172 inline G4double GetMinimumEpsilonStep() const; // Min for relative accuracy
173 inline void SetMinimumEpsilonStep( G4double newEpsMin ); // of any step
175 inline void SetMaximumEpsilonStep( G4double newEpsMax );
176
177 /**
178 * Methods to obtain / change the size of the largest step the method
179 * will undertake. The Reset method uses the world volume's.
180 */
181 void SetLargestAcceptableStep( G4double newBigDist );
184
185 /**
186 * Methods to control extra Multiplier parameter for limiting long steps.
187 */
190
191 /**
192 * Methods to Control minimum 'directional' distance in case of
193 * too-large step.
194 */
196 void SetMinBigDistance(G4double val);
197
198 /**
199 * Sets the filter that examines & stores 'intermediate'
200 * curved trajectory points.
201 * @note Currently only position is stored.
202 */
204
205 /**
206 * Accesses the points which have passed by the filter.
207 * @note Responsibility for deleting the points lies with the client.
208 * This method MUST BE called exactly ONCE per step.
209 */
210 std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
211
212 /**
213 * Clears the State of this class and its current associates.
214 * @note The current field manager & chord finder will also be called.
215 */
217
218 /**
219 * Setter for global field manager. Updates the state.
220 */
221 inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
222
223 /**
224 * Toggles & views parameter for using safety to discard unneccesary calls
225 * to the navigator (thus 'optimising' performance).
226 */
229
230 /**
231 * Intersects the chord from StartPointA to EndPointB and returns
232 * whether an intersection occurred.
233 * @note Safety is changed!
234 */
235 inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
236 const G4ThreeVector& EndPointB,
237 G4double& NewSafety,
238 G4double& LinearStepLength,
239 G4ThreeVector& IntersectionPoint);
240
241 /**
242 * Returns if it is the first step in the volume.
243 */
245
246 /**
247 * Returns if it is the last step in the volume.
248 */
250
251 /**
252 * Initialises track flags.
253 */
254 inline void PrepareNewTrack();
255
256 /**
257 * Changes or gets the object which calculates the exact intersection
258 * point with the next boundary.
259 */
262
263 /**
264 * Controls the parameter which enables the temporary 'relaxation' which
265 * ensures that chord segments are short enough so that their sagitta is
266 * small than delta-chord parameter.
267 * The Set method increases the value of delta-chord temporarily, doubling
268 * it once the number of iterations substeps reach value of
269 * 'IncreaseChordDistanceThreshold'. It is also doubled again every time
270 * the iteration count reaches a multiple of this value.
271 * @note The delta-chord is reset to its original value at the end of
272 * each call to ComputeStep().
273 */
276
277 /**
278 * Accessors.
279 */
281 inline G4double GetDeltaOneStep() const;
282
283 /**
284 * Auxiliary methods.
285 * @note Their results can/will change during propagation.
286 */
289
290 /**
291 * Accessor and modifier for navigator.
292 */
293 inline void SetNavigatorForPropagating(G4Navigator* SimpleOrMultiNavigator);
295
296 /**
297 * Accessors and modifiers for no-zero steps threshold.
298 */
299 inline void SetThresholdNoZeroStep( G4int noAct,
300 G4int noHarsh,
301 G4int noAbandon );
304 inline void SetZeroStepThreshold( G4double newLength );
305
306 /**
307 * Updates the Locator with parameters from this class and from current
308 * field manager.
309 */
311
312 protected:
313
314 /**
315 * Logging methods.
316 */
317 void PrintStepLengthDiagnostic( G4double currentProposedStepLength,
318 G4double decreaseFactor,
319 G4double stepTrial,
320 const G4FieldTrack& aFieldTrack);
321 void ReportLoopingParticle( G4int count, G4double StepTaken,
322 G4double stepRequest, const char* methodName,
323 const G4ThreeVector& momentumVec,
324 G4VPhysicalVolume* physVol);
325 void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep,
326 G4double lastTriedStep, G4VPhysicalVolume* physVol);
327
328 private:
329
330 // ----------------------------------------------------------------------
331 // DATA Members
332 // ----------------------------------------------------------------------
333
334 // ==================================================================
335 // INVARIANTS - Must not change during tracking
336
337 // ** PARAMETERS -----------
338 G4int fMax_loop_count = 1000;
339 // Limit for the number of sub-steps taken in one call to ComputeStep
340 G4int fIncreaseChordDistanceThreshold = 100;
341 G4bool fUseSafetyForOptimisation = true;
342 // (false) is less sensitive to incorrect safety
343
344 // Thresholds for identifying "abnormal" cases - which cause looping
345 //
346 G4int fActionThreshold_NoZeroSteps = 2; // Threshold # - above it act
347 G4int fSevereActionThreshold_NoZeroSteps = 10; // Threshold # to act harshly
348 G4int fAbandonThreshold_NoZeroSteps = 50; // Threshold # to abandon
349 G4double fZeroStepThreshold = 0.0;
350 // Threshold *length* for counting of tiny or 'zero' steps
351
352 // Parameters related to handling of very large steps which
353 // occur typically in large volumes with vacuum or very thin gas
354 //
355 G4double fLargestAcceptableStep;
356 // Maximum size of a step - for optimization (and to avoid problems)
357 G4double fMaxStepSizeMultiplier = 3;
358 // Multiplier for directional exit distance used as extra long-step limit
359 G4double fMinBigDistance= 100. ; // * CLHEP::mm
360 // Minimum distance added to directional exit distance
361 // ** End of PARAMETERS -----
362
363 G4double kCarTolerance;
364 // Geometrical tolerance defining surface thickness
365
366 G4bool fAllocatedLocator; // Book-keeping
367
368 // --------------------------------------------------------
369 // ** Dependent Objects - to which work is delegated
370
371 G4FieldManager* fDetectorFieldMgr;
372 // The Field Manager of the whole Detector. (default)
373
374 G4VIntersectionLocator* fIntersectionLocator;
375 // Refines candidate intersection
376
377 G4VCurvedTrajectoryFilter* fpTrajectoryFilter = nullptr;
378 // The filter encapsulates the algorithm which selects which
379 // intermediate points should be stored in a trajectory.
380 // When it is NULL, no intermediate points will be stored.
381 // Else PIF::ComputeStep must submit (all) intermediate
382 // points it calculates, to this filter. (jacek 04/11/2002)
383
384 G4Navigator* fNavigator;
385 // Set externally - only by tracking / run manager
386 //
387 // ** End of Dependent Objects ----------------------------
388
389 // End of INVARIANTS
390 // ==================================================================
391
392 // STATE information
393 // -----------------
394 G4FieldManager* fCurrentFieldMgr;
395 // The Field Manager of the current volume (may be the global)
396 G4bool fSetFieldMgr = false; // Has it been set for the current step?
397
398 // Parameters of current step
399 //
400 G4double fEpsilonStep; // Relative accuracy of current Step
401 G4FieldTrack End_PointAndTangent; // End point storage
402 G4bool fParticleIsLooping = false;
403 G4int fNoZeroStep = 0; // Count of zero Steps
404
405 // State used for Optimisation
406 //
407 G4double fFull_CurveLen_of_LastAttempt = -1;
408 G4double fLast_ProposedStepLength = -1;
409 // Previous step information -- for use in adjust step size
410 G4ThreeVector fPreviousSftOrigin;
411 G4double fPreviousSafety = 0.0;
412 // Last safety origin & value: for optimisation
413
414 G4int fVerboseLevel = 0;
415 G4bool fVerbTracePiF = false;
416 G4bool fCheck = false;
417 // For debugging purposes
418
419 G4bool fFirstStepInVolume = true;
420 G4bool fLastStepInVolume = true;
421 G4bool fNewTrack = true;
422};
423
424// Inline methods
425//
426#include "G4PropagatorInField.icc"
427
428#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4ChordFinder is a class that provides Runge-Kutta integration of motion ODE and also has a method th...
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
G4FieldManager is a manager (store) for a pointer to the Field subclass that describes the field of a...
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
G4Navigator is a class for use by the tracking management, able to obtain/calculate dynamic tracking ...
void SetVerboseTrace(G4bool enable)
G4ThreeVector EndPosition() const
void SetThresholdNoZeroStep(G4int noAct, G4int noHarsh, G4int noAbandon)
G4int GetIterationsToIncreaseChordDistance() const
G4int GetThresholdNoZeroSteps(G4int i)
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4double GetDeltaOneStep() const
G4FieldManager * GetCurrentFieldManager()
void CheckMode(G4bool mode)
G4int SetVerboseLevel(G4int verbose)
void ResetLargestAcceptableStep()
G4bool IsLastStepInVolume()
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
G4PropagatorInField(G4Navigator *theNavigator, G4FieldManager *detectorFieldMgr, G4VIntersectionLocator *vLocator=nullptr)
G4double GetZeroStepThreshold()
void SetMaxLoopCount(G4int new_max)
void SetNavigatorForPropagating(G4Navigator *SimpleOrMultiNavigator)
G4int Verbose() const
G4bool GetUseSafetyForOptimization()
G4bool IsParticleLooping() const
void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)
void SetZeroStepThreshold(G4double newLength)
void SetDetectorFieldManager(G4FieldManager *newGlobalFieldManager)
G4double GetMaximumEpsilonStep() const
void SetMinBigDistance(G4double val)
G4ChordFinder * GetChordFinder()
void SetIterationsToIncreaseChordDistance(G4int numIters)
G4EquationOfMotion * GetCurrentEquationOfMotion()
void SetUseSafetyForOptimization(G4bool)
G4double GetDeltaIntersection() const
G4bool IsFirstStepInVolume()
G4double GetEpsilonStep() const
void SetMaximumEpsilonStep(G4double newEpsMax)
G4Navigator * GetNavigatorForPropagating()
void SetMaxStepSizeMultiplier(G4double vm)
G4FieldTrack GetEndState() const
G4ThreeVector EndMomentumDir() const
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
G4VIntersectionLocator * GetIntersectionLocator()
G4int GetVerboseLevel() const
void SetMinimumEpsilonStep(G4double newEpsMin)
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
void ReportLoopingParticle(G4int count, G4double StepTaken, G4double stepRequest, const char *methodName, const G4ThreeVector &momentumVec, G4VPhysicalVolume *physVol)
G4int GetMaxLoopCount() const
void SetIntersectionLocator(G4VIntersectionLocator *pLocator)
G4bool GetVerboseTrace()
void SetTrajectoryFilter(G4VCurvedTrajectoryFilter *filter)
G4double GetMinimumEpsilonStep() const
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=nullptr, G4bool canRelaxDeltaChord=false)
void SetLargestAcceptableStep(G4double newBigDist)
void SetEpsilonStep(G4double newEps)
G4VCurvedTrajectoryFilter defines a filter for deciding which intermediate points on a curved traject...
G4VIntersectionLocator is a base class for the calculation of the intersection point with a boundary ...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....