Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/navigation/include/G4PropagatorInField.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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"
 49 #include "G4VIntersectionLocator.hh"
 50 
 51 class G4ChordFinder; 
 52 
 53 class G4Navigator;
 54 class G4VPhysicalVolume;
 55 class G4VCurvedTrajectoryFilter;
 56 
 57 class G4PropagatorInField
 58 {
 59 
 60  public:  // with description
 61 
 62    G4PropagatorInField( G4Navigator* theNavigator, 
 63                         G4FieldManager* detectorFieldMgr,
 64                         G4VIntersectionLocator* vLocator = nullptr );
 65   ~G4PropagatorInField();
 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;       
 75    inline G4ThreeVector EndMomentumDir() 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 
 84    G4FieldManager* FindAndSetFieldManager(G4VPhysicalVolume* pCurrentPhysVol);
 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  
 90    inline G4ChordFinder* GetChordFinder();
 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 );
 98    inline G4bool GetVerboseTrace();
 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
118    inline G4double GetMaximumEpsilonStep() const;
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    void     SetLargestAcceptableStep( G4double newBigDist );
127    G4double GetLargestAcceptableStep();
128    void     ResetLargestAcceptableStep();
129      // Obtain / change the size of the largest step the method will undertake
130      // Reset method uses the world volume's 
131 
132    G4double GetMaxStepSizeMultiplier();
133    void     SetMaxStepSizeMultiplier(G4double vm);
134      // Control extra Multiplier parameter for limiting long steps.
135    G4double GetMinBigDistance();
136    void     SetMinBigDistance(G4double val);
137      // Control minimum 'directional' distance in case of too-large step
138 
139    void SetTrajectoryFilter(G4VCurvedTrajectoryFilter* filter);
140      // Set the filter that examines & stores 'intermediate' 
141      // curved trajectory points.  Currently only position is stored.
142 
143    std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
144      // Access the points which have passed by the filter.
145      // Responsibility for deleting the points lies with the client.
146      // This method MUST BE called exactly ONCE per step. 
147 
148    void ClearPropagatorState();
149      // Clear all the State of this class and its current associates
150      // --> the current field manager & chord finder will also be called
151 
152    inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
153      // Update this (dangerous) state -- for the time being
154   
155    inline void   SetUseSafetyForOptimization( G4bool );
156    inline G4bool GetUseSafetyForOptimization();
157      // Toggle & view parameter for using safety to discard 
158      // unneccesary calls to navigator (thus 'optimising' performance)
159    inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
160                                  const G4ThreeVector& EndPointB,
161                                        G4double&      NewSafety,
162                                        G4double&      LinearStepLength,
163                                        G4ThreeVector& IntersectionPoint);
164      // Intersect the chord from StartPointA to EndPointB
165      // and return whether an intersection occurred
166      // NOTE: Safety is changed!
167 
168    inline G4bool IsFirstStepInVolume();
169    inline G4bool IsLastStepInVolume();
170    inline void PrepareNewTrack();
171       
172    inline G4VIntersectionLocator* GetIntersectionLocator();
173    inline void SetIntersectionLocator(G4VIntersectionLocator* pLocator );
174      // Change or get the object which calculates the exact 
175      // intersection point with the next boundary
176 
177    inline G4int GetIterationsToIncreaseChordDistance() const;
178    inline void  SetIterationsToIncreaseChordDistance(G4int numIters);
179      // Control the parameter which enables the temporary 'relaxation'
180      //   which ensures that chord segments are short enough so that
181      //   their sagitta is small than delta-chord parameter.
182      // The Set method increases the value of delta-chord temporarily,
183      //   doubling it once the number of iterations substeps reach
184      //   value of 'IncreaseChordDistanceThreshold'.  It is also doubled
185      //   again every time the iteration count reaches a multiple of this
186      //   value.
187      // Note: delta-chord is reset to its original value at the end of
188      //   each call to ComputeStep.
189 
190  public:  // without description
191 
192    inline G4double GetDeltaIntersection() const;
193    inline G4double GetDeltaOneStep() const;
194 
195    inline G4FieldManager* GetCurrentFieldManager();
196    inline G4EquationOfMotion* GetCurrentEquationOfMotion();
197       // Auxiliary methods - their results can/will change during propagation
198 
199    inline void SetNavigatorForPropagating(G4Navigator* SimpleOrMultiNavigator); 
200    inline G4Navigator* GetNavigatorForPropagating();
201 
202    inline void SetThresholdNoZeroStep( G4int noAct,
203                                        G4int noHarsh,
204                                        G4int noAbandon );
205    inline G4int GetThresholdNoZeroSteps( G4int i ); 
206 
207    inline G4double GetZeroStepThreshold(); 
208    inline void     SetZeroStepThreshold( G4double newLength ); 
209    
210    void RefreshIntersectionLocator(); 
211      // Update the Locator with parameters from this class
212      // and from current field manager
213 
214  protected:  // without description
215 
216    void PrintStepLengthDiagnostic( G4double      currentProposedStepLength,
217                                    G4double      decreaseFactor,
218                                    G4double      stepTrial,
219                              const G4FieldTrack& aFieldTrack);
220 
221    void ReportLoopingParticle( G4int count,  G4double StepTaken,
222                                G4double stepRequest, const char* methodName,
223                                const G4ThreeVector&      momentumVec,
224                                G4VPhysicalVolume* physVol);
225    void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep,
226                             G4double lastTriedStep, G4VPhysicalVolume* physVol);
227 
228  private:
229 
230    // ----------------------------------------------------------------------
231    //  DATA Members
232    // ----------------------------------------------------------------------
233 
234    //  ==================================================================
235    //  INVARIANTS - Must not change during tracking
236 
237    //  ** PARAMETERS -----------
238    G4int fMax_loop_count = 1000;
239      // Limit for the number of sub-steps taken in one call to ComputeStep
240    G4int fIncreaseChordDistanceThreshold = 100;
241    G4bool fUseSafetyForOptimisation = true;
242      // (false) is less sensitive to incorrect safety
243 
244    //  Thresholds for identifying "abnormal" cases - which cause looping
245    //
246    G4int fActionThreshold_NoZeroSteps = 2;        // Threshold # - above it act
247    G4int fSevereActionThreshold_NoZeroSteps = 10; // Threshold # to act harshly
248    G4int fAbandonThreshold_NoZeroSteps = 50;      // Threshold # to abandon
249    G4double fZeroStepThreshold = 0.0; 
250      // Threshold *length* for counting of tiny or 'zero' steps 
251 
252    // Parameters related to handling of very large steps which
253    //   occur typically in large volumes with vacuum or very thin gas
254    G4double fLargestAcceptableStep;
255      // Maximum size of a step - for optimization (and to avoid problems)
256    G4double fMaxStepSizeMultiplier = 3;
257      // Multiplier for directional exit distance used as extra long-step limit 
258    G4double fMinBigDistance= 100. ; // * CLHEP::mm
259      // Minimum distance added to directional exit distance
260    //  ** End of PARAMETERS -----
261 
262    G4double kCarTolerance;
263        // Geometrical tolerance defining surface thickness
264 
265    G4bool fAllocatedLocator;                    //  Book-keeping
266 
267    //  --------------------------------------------------------
268    //  ** Dependent Objects - to which work is delegated 
269 
270    G4FieldManager* fDetectorFieldMgr; 
271        // The  Field Manager of the whole Detector.  (default)
272 
273    G4VIntersectionLocator* fIntersectionLocator;
274      // Refines candidate intersection
275 
276    G4VCurvedTrajectoryFilter* fpTrajectoryFilter = nullptr;
277      // The filter encapsulates the algorithm which selects which
278      // intermediate points should be stored in a trajectory. 
279      // When it is NULL, no intermediate points will be stored.
280      // Else PIF::ComputeStep must submit (all) intermediate
281      // points it calculates, to this filter.  (jacek 04/11/2002)
282 
283    G4Navigator* fNavigator;
284      // Set externally - only by tracking / run manager
285    //
286    //  ** End of Dependent Objects ----------------------------
287 
288    //  End of INVARIANTS 
289    //  ==================================================================
290 
291    //  STATE information
292    //  -----------------
293    G4FieldManager* fCurrentFieldMgr;
294      // The  Field Manager of the current volume (may be the global)
295    G4bool fSetFieldMgr = false;  // Has it been set for the current step?
296 
297    // Parameters of current step
298    G4double fEpsilonStep;            // Relative accuracy of current Step
299    G4FieldTrack End_PointAndTangent; // End point storage
300    G4bool fParticleIsLooping = false;
301    G4int fNoZeroStep = 0;            // Count of zero Steps
302 
303    // State used for Optimisation
304    G4double fFull_CurveLen_of_LastAttempt = -1; 
305    G4double fLast_ProposedStepLength = -1; 
306      // Previous step information -- for use in adjust step size
307    G4ThreeVector fPreviousSftOrigin;
308    G4double fPreviousSafety = 0.0; 
309      // Last safety origin & value: for optimisation
310 
311    G4int fVerboseLevel = 0;
312    G4bool fVerbTracePiF = false;
313    G4bool fCheck = false;
314      // For debugging purposes
315 
316    G4bool fFirstStepInVolume = true; 
317    G4bool fLastStepInVolume = true; 
318    G4bool fNewTrack = true;
319 };
320 
321 // Inline methods
322 //
323 #include "G4PropagatorInField.icc"
324 
325 #endif 
326