Geant4 10.7.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// Class 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
34// History:
35// -------
36// 25.10.96 John Apostolakis, design and implementation
37// 25.03.97 John Apostolakis, adaptation for G4Transportation and cleanup
38// 8.11.02 John Apostolakis, changes to enable use of safety in intersecting
39// ---------------------------------------------------------------------------
40#ifndef G4PropagatorInField_hh
41#define G4PropagatorInField_hh 1
42
43#include "G4Types.hh"
44
45#include <vector>
46
47#include "G4FieldTrack.hh"
48#include "G4FieldManager.hh"
50
51class G4ChordFinder;
52
53class G4Navigator;
56
58{
59
60 public: // with description
61
62 G4PropagatorInField( G4Navigator* theNavigator,
63 G4FieldManager* detectorFieldMgr,
64 G4VIntersectionLocator* vLocator = nullptr );
66
67 G4double ComputeStep( G4FieldTrack& pFieldTrack,
68 G4double pCurrentProposedStepLength,
69 G4double& pNewSafety,
70 G4VPhysicalVolume* pPhysVol = nullptr,
71 G4bool canRelaxDeltaChord = false);
72 // Compute the next geometric Step
73
74 inline G4ThreeVector EndPosition() const;
76 inline G4bool IsParticleLooping() const;
77 // Return the state after the Step
78
79 inline G4double GetEpsilonStep() const;
80 // Relative accuracy for current Step (Calc.)
81 inline void SetEpsilonStep(G4double newEps);
82 // The ratio DeltaOneStep()/h_current_step
83
85 // Set (and return) the correct field manager (global or local),
86 // if it exists.
87 // Should be called before ComputeStep is called;
88 // Currently, ComputeStep will call it, if it has not been called.
89
91
92 G4int SetVerboseLevel( G4int verbose );
93 inline G4int GetVerboseLevel() const;
94 inline G4int Verbose() const;
95 inline void CheckMode(G4bool mode);
96
97 inline void SetVerboseTrace( G4bool enable );
99 // Tracing key parts of Compute Step
100
101 inline G4int GetMaxLoopCount() const;
102 inline void SetMaxLoopCount( G4int new_max );
103 // A maximum for the number of substeps that a particle can take.
104 // Above this number it is signaled as 'looping'.
105
106 void printStatus( const G4FieldTrack& startFT,
107 const G4FieldTrack& currentFT,
108 G4double requestStep,
109 G4double safety,
110 G4int step,
111 G4VPhysicalVolume* startVolume);
112 // Print Method - useful mostly for debugging.
113
114 inline G4FieldTrack GetEndState() const;
115
116 inline G4double GetMinimumEpsilonStep() const; // Min for relative accuracy
117 inline void SetMinimumEpsilonStep( G4double newEpsMin ); // of any step
119 inline void SetMaximumEpsilonStep( G4double newEpsMax );
120 // The 4 above methods are now obsolescent but *for now* will work
121 // They are being replaced by same-name methods in G4FieldManager,
122 // allowing the specialisation in different volumes.
123 // Their new behaviour is to change the values for the global field
124 // manager
125
126 inline void SetLargestAcceptableStep( G4double newBigDist );
128
130 // Set the filter that examines & stores 'intermediate'
131 // curved trajectory points. Currently only position is stored.
132
133 std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
134 // Access the points which have passed by the filter.
135 // Responsibility for deleting the points lies with the client.
136 // This method MUST BE called exactly ONCE per step.
137
139 // Clear all the State of this class and its current associates
140 // --> the current field manager & chord finder will also be called
141
142 inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
143 // Update this (dangerous) state -- for the time being
144
147 // Toggle & view parameter for using safety to discard
148 // unneccesary calls to navigator (thus 'optimising' performance)
149 inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
150 const G4ThreeVector& EndPointB,
151 G4double& NewSafety,
152 G4double& LinearStepLength,
153 G4ThreeVector& IntersectionPoint);
154 // Intersect the chord from StartPointA to EndPointB
155 // and return whether an intersection occurred
156 // NOTE: Safety is changed!
157
160 inline void PrepareNewTrack();
161
164 // Change or get the object which calculates the exact
165 // intersection point with the next boundary
166
169 // Control the parameter which enables the temporary 'relaxation'
170 // which ensures that chord segments are short enough so that
171 // their sagitta is small than delta-chord parameter.
172 // The Set method increases the value of delta-chord temporarily,
173 // doubling it once the number of iterations substeps reach
174 // value of 'IncreaseChordDistanceThreshold'. It is also doubled
175 // again every time the iteration count reaches a multiple of this
176 // value.
177 // Note: delta-chord is reset to its original value at the end of
178 // each call to ComputeStep.
179
180 public: // without description
181
183 inline G4double GetDeltaOneStep() const;
184
187 // Auxiliary methods - their results can/will change during propagation
188
189 inline void SetNavigatorForPropagating(G4Navigator* SimpleOrMultiNavigator);
191
192 inline void SetThresholdNoZeroStep( G4int noAct,
193 G4int noHarsh,
194 G4int noAbandon );
196
198 inline void SetZeroStepThreshold( G4double newLength );
199
201 // Update the Locator with parameters from this class
202 // and from current field manager
203
204 protected: // without description
205
206 void PrintStepLengthDiagnostic( G4double currentProposedStepLength,
207 G4double decreaseFactor,
208 G4double stepTrial,
209 const G4FieldTrack& aFieldTrack);
210
211 void ReportLoopingParticle( G4int count, G4double StepTaken,
212 G4double stepRequest, const char* methodName,
213 G4ThreeVector momentumVec,
214 G4VPhysicalVolume* physVol);
215 void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep,
216 G4double lastTriedStep, G4VPhysicalVolume* physVol);
217
218 private:
219
220 // ----------------------------------------------------------------------
221 // DATA Members
222 // ----------------------------------------------------------------------
223
224 // ==================================================================
225 // INVARIANTS - Must not change during tracking
226
227 // ** PARAMETERS -----------
228 G4int fMax_loop_count = 1000;
229 // Limit for the number of sub-steps taken in one call to ComputeStep
230 G4int fIncreaseChordDistanceThreshold = 100;
231 G4bool fUseSafetyForOptimisation = true;
232 // (false) is less sensitive to incorrect safety
233
234 // Thresholds for identifying "abnormal" cases - which cause looping
235 //
236 G4int fActionThreshold_NoZeroSteps = 2; // Threshold # - above it act
237 G4int fSevereActionThreshold_NoZeroSteps = 10; // Threshold # to act harshly
238 G4int fAbandonThreshold_NoZeroSteps = 50; // Threshold # to abandon
239 G4double fZeroStepThreshold = 0.0;
240 // Threshold *length* for counting of tiny or 'zero' steps
241
242 G4double fLargestAcceptableStep;
243 // Maximum size of a step - for optimization (and to avoid problems)
244 // ** End of PARAMETERS -----
245
246 G4double kCarTolerance;
247 // Geometrical tolerance defining surface thickness
248
249 G4bool fAllocatedLocator; // Book-keeping
250
251 // --------------------------------------------------------
252 // ** Dependent Objects - to which work is delegated
253
254 G4FieldManager* fDetectorFieldMgr;
255 // The Field Manager of the whole Detector. (default)
256
257 G4VIntersectionLocator* fIntersectionLocator;
258 // Refines candidate intersection
259
260 G4VCurvedTrajectoryFilter* fpTrajectoryFilter = nullptr;
261 // The filter encapsulates the algorithm which selects which
262 // intermediate points should be stored in a trajectory.
263 // When it is NULL, no intermediate points will be stored.
264 // Else PIF::ComputeStep must submit (all) intermediate
265 // points it calculates, to this filter. (jacek 04/11/2002)
266
267 G4Navigator* fNavigator;
268 // Set externally - only by tracking / run manager
269 //
270 // ** End of Dependent Objects ----------------------------
271
272 // End of INVARIANTS
273 // ==================================================================
274
275 // STATE information
276 // -----------------
277 G4FieldManager* fCurrentFieldMgr;
278 // The Field Manager of the current volume (may be the global)
279 G4bool fSetFieldMgr = false; // Has it been set for the current step?
280
281 // Parameters of current step
282 G4double fEpsilonStep; // Relative accuracy of current Step
283 G4FieldTrack End_PointAndTangent; // End point storage
284 G4bool fParticleIsLooping = false;
285 G4int fNoZeroStep = 0; // Count of zero Steps
286
287 // State used for Optimisation
288 G4double fFull_CurveLen_of_LastAttempt = -1;
289 G4double fLast_ProposedStepLength = -1;
290 // Previous step information -- for use in adjust step size
291 G4ThreeVector fPreviousSftOrigin;
292 G4double fPreviousSafety = 0.0;
293 // Last safety origin & value: for optimisation
294
295 G4int fVerboseLevel = 0;
296 G4bool fVerbTracePiF = false;
297 G4bool fCheck = false;
298 // For debugging purposes
299
300 G4bool fFirstStepInVolume = true;
301 G4bool fLastStepInVolume = true;
302 G4bool fNewTrack = true;
303};
304
305// Inline methods
306//
307#include "G4PropagatorInField.icc"
308
309#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
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)
G4bool IsLastStepInVolume()
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
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
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()
G4FieldTrack GetEndState() const
G4ThreeVector EndMomentumDir() const
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)
G4VIntersectionLocator * GetIntersectionLocator()
G4int GetVerboseLevel() const
void SetMinimumEpsilonStep(G4double newEpsMin)
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4int GetMaxLoopCount() const
void SetIntersectionLocator(G4VIntersectionLocator *pLocator)
G4bool GetVerboseTrace()
G4double GetLargestAcceptableStep()
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)