Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/navigation/include/G4PathFinder.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 G4PathFinder 
 27 //
 28 // Class description:
 29 // 
 30 // This class directs the lock-stepped propagation of a track in the 
 31 // 'mass' and other parallel geometries.  It ensures that tracking 
 32 // in a magnetic field sees these parallel geometries at each trial step, 
 33 // and that the earliest boundary limits the step.
 34 // 
 35 // For the movement in field, it relies on the class G4PropagatorInField
 36 
 37 // History:
 38 // -------
 39 //  7.10.05 John Apostolakis,  Draft design 
 40 // 26.04.06 John Apostolakis,  Revised design and first implementation 
 41 // ---------------------------------------------------------------------------
 42 #ifndef G4PATHFINDER_HH 
 43 #define G4PATHFINDER_HH  1
 44 
 45 #include <vector>
 46 #include "G4Types.hh"
 47 
 48 #include "G4FieldTrack.hh"
 49 
 50 class G4TransportationManager; 
 51 class G4Navigator;
 52 
 53 #include "G4MultiNavigator.hh"
 54 #include "G4TouchableHandle.hh"
 55 
 56 class G4PropagatorInField;
 57 
 58 class G4PathFinder
 59 {
 60 
 61  public:  // with description
 62 
 63    static G4PathFinder* GetInstance();
 64      // Retrieve singleton instance and create it if not existing.
 65 
 66    static G4PathFinder* GetInstanceIfExist();
 67      // Retrieve singleton instance pointer.
 68 
 69    G4double ComputeStep( const G4FieldTrack& pFieldTrack,
 70                          G4double  pCurrentProposedStepLength,
 71                          G4int     navigatorId, // Identifies the geometry
 72                          G4int     stepNo,      // See next step/check 
 73                          G4double& pNewSafety,  // Only for this geometry
 74                          ELimited& limitedStep,      
 75                          G4FieldTrack& EndState, 
 76                          G4VPhysicalVolume* currentVolume );
 77      // Compute the next geometric Step -- curved or linear
 78      // If it is called with a larger 'stepNo' it will execute a new step;
 79      // if 'stepNo' is same as last call, then the results for 
 80      // the geometry with Id. number 'navigatorId' will be returned. 
 81 
 82    void Locate( const G4ThreeVector& position, 
 83                 const G4ThreeVector& direction,
 84                       G4bool relativeSearch = true); 
 85      // Make primary relocation of global point in all navigators,
 86      // and update them.
 87 
 88    void ReLocate( const G4ThreeVector& position ); 
 89      // Make secondary relocation of global point (within safety only) 
 90      // in all navigators, and update them.
 91 
 92    void PrepareNewTrack( const G4ThreeVector& position,
 93                          const G4ThreeVector& direction,
 94                                G4VPhysicalVolume* massStartVol = nullptr); 
 95      // Check and cache set of active navigators.
 96 
 97    void EndTrack();
 98      // Signal end of tracking of current track.  
 99      //   Reset internal state
100      //   Inform TransportationManager to use 'ordinary' Navigator
101 
102     G4TouchableHandle CreateTouchableHandle( G4int navId ) const;
103     inline G4VPhysicalVolume* GetLocatedVolume( G4int navId ) const;
104 
105    // -----------------------------------------------------------------
106   
107    inline G4bool IsParticleLooping() const;
108 
109    inline G4double GetCurrentSafety() const;
110      // Minimum value of safety after last ComputeStep
111    inline G4double GetMinimumStep() const;      
112      // Get the minimum step size from the last ComputeStep call
113      //   - in case full step is taken, this is kInfinity
114    inline unsigned int  GetNumberGeometriesLimitingStep() const; 
115 
116    G4double ComputeSafety( const G4ThreeVector& globalPoint); 
117      // Recompute safety for the relevant point the endpoint of the last step!!
118      // Maintain vector of individual safety values (for next method)
119 
120    G4double ObtainSafety( G4int navId, G4ThreeVector& globalCenterPoint );
121      // Obtain safety for navigator/geometry navId for last point 'computed'
122      //   --> last point for which ComputeSafety was called
123      //   Returns the point (center) for which this safety is valid
124 
125    void EnableParallelNavigation( G4bool enableChoice = true ); 
126      // Must call it to ensure that PathFinder is prepared,  
127      // especially for curved tracks. If true it switches PropagatorInField
128      // to use MultiNavigator. Must call it with false to undo (=PiF use
129      // Navigator for tracking!)
130 
131    inline G4int  SetVerboseLevel(G4int lev = -1);
132 
133  public:  // with description
134 
135    inline G4int GetMaxLoopCount() const;
136    inline void  SetMaxLoopCount( G4int new_max );
137      // A maximum for the number of steps that a (looping) particle can take
138 
139  public:  // without description
140 
141    inline void MovePoint();
142      // Signal that location will be moved -- internal use primarily
143 
144    // To provide best compatibility between Coupled and Old Transportation
145    // the next two methods are provided:
146 
147    G4double LastPreSafety( G4int navId, G4ThreeVector& globalCenterPoint,
148                            G4double& minSafety ); 
149      // Obtain last safety needed in ComputeStep (for geometry navId)
150      // --> last point at which ComputeStep recalculated safety
151      //     Returns the point (center) for which this safety is valid
152      //     and also the minimum safety over all navigators (i.e. full)
153 
154    void PushPostSafetyToPreSafety(); 
155      // Tell PathFinder to copy PostStep Safety to PreSafety
156      // (for use at next step)
157 
158    G4String& LimitedString( ELimited lim );
159      // Convert ELimited to string
160 
161    ~G4PathFinder();
162      // Destructor
163 
164  protected:  // without description
165 
166    G4double DoNextLinearStep( const G4FieldTrack& FieldTrack,
167                                     G4double      proposedStepLength); 
168 
169    G4double DoNextCurvedStep( const G4FieldTrack& FieldTrack,
170                                     G4double      proposedStepLength,
171                                     G4VPhysicalVolume* pCurrentPhysVolume); 
172 
173   void WhichLimited();
174   void PrintLimited();
175     // Print key details out for debugging
176 
177   inline G4bool UseSafetyForOptimization( G4bool );
178     // Whether use safety to discard unneccesary calls to navigator
179 
180   void ReportMove( const G4ThreeVector& OldV,
181                    const G4ThreeVector& NewV,
182                    const G4String& Quantity ) const;
183     // Helper method to report movement (likely of initial point)
184 
185  protected:
186 
187    G4PathFinder();  //  Singleton 
188 
189    inline G4Navigator* GetNavigator(G4int n) const;
190 
191  private:
192 
193    // ----------------------------------------------------------------------
194    //  DATA Members
195    // ----------------------------------------------------------------------
196 
197    G4MultiNavigator* fpMultiNavigator; 
198      // Object that enables G4PropagatorInField to see many geometries
199 
200    G4int fNoActiveNavigators = 0; 
201    G4bool fNewTrack = false; // Flag a new track (ensure first step)
202 
203    static const G4int fMaxNav = 16;
204 
205    // Global state (retained during stepping for one track)
206 
207    G4Navigator*  fpNavigator[fMaxNav];
208 
209    // State changed in a step computation
210 
211    ELimited      fLimitedStep[fMaxNav];
212    G4bool        fLimitTruth[fMaxNav];
213    G4double      fCurrentStepSize[fMaxNav]; 
214    G4int fNoGeometriesLimiting = 0;  // How many processes contribute to limit
215 
216    G4ThreeVector fPreSafetyLocation;
217      // last initial position for which safety evaluated
218    G4double      fPreSafetyMinValue = -1.0;
219      // - corresponding value of full safety
220    G4double      fPreSafetyValues[ fMaxNav ];
221      // Safeties for the above point
222 
223    // This part of the state can be retained for several calls --> CARE
224 
225    G4ThreeVector fPreStepLocation;
226      // point where last ComputeStep called
227    G4double fMinSafety_PreStepPt = -1.0;
228      // - corresponding value of full safety
229    G4double fCurrentPreStepSafety[ fMaxNav ];
230      // Safeties for the above point.
231      // This changes at each step, so it can differ when steps
232      // inside min-safety are made
233 
234    G4bool fPreStepCenterRenewed = false;
235      // Whether PreSafety coincides with PreStep point 
236 
237    G4double fMinStep = -1.0;  // As reported by Navigators -- can be kInfinity
238    G4double fTrueMinStep = -1.0;  // Corrected in case >= proposed
239 
240    // State after calling 'locate'
241    //
242    G4VPhysicalVolume* fLocatedVolume[fMaxNav];
243    G4ThreeVector      fLastLocatedPosition; 
244 
245    // State after calling 'ComputeStep'
246    // (others member variables will be affected)
247    //
248    G4FieldTrack fEndState;  // Point, velocity, ... at proposed step end
249    G4bool fFieldExertedForce = false; // In current proposed step
250 
251    G4bool fRelocatedPoint = false; // Signals that point was or is being moved 
252                                    // from the position of the last location or
253                                    // the endpoint resulting from ComputeStep()
254                                    // -- invalidates fEndState
255 
256    // State for 'ComputeSafety' and related methods
257    //
258    G4ThreeVector fSafetyLocation;
259      // point where ComputeSafety is called
260    G4double      fMinSafety_atSafLocation = -1.0;
261      // - corresponding value of safety
262    G4double      fNewSafetyComputed[ fMaxNav ];
263      // Safeties for last ComputeSafety
264 
265    // State for Step numbers
266    //
267    G4int fLastStepNo = -1, fCurrentStepNo = -1; 
268 
269    G4int fVerboseLevel = 0;  // For debugging purposes
270 
271    G4TransportationManager* fpTransportManager; // Cache for frequent use
272    G4PropagatorInField* fpFieldPropagator;
273 
274    G4double kCarTolerance;
275 
276    static G4ThreadLocal G4PathFinder* fpPathFinder;
277 };
278 
279 // ********************************************************************
280 // Inline methods.
281 // ********************************************************************
282 
283 inline G4VPhysicalVolume* G4PathFinder::GetLocatedVolume( G4int navId ) const
284 {
285   G4VPhysicalVolume* vol = nullptr;  
286   if( (navId < fMaxNav) && (navId >= 0) ) { vol= fLocatedVolume[navId]; }
287   return vol; 
288 }
289 
290 inline G4int G4PathFinder::SetVerboseLevel(G4int newLevel)
291 {
292   G4int old = fVerboseLevel;
293   fVerboseLevel = newLevel;
294   return old;
295 }
296 
297 inline G4double G4PathFinder::GetMinimumStep() const
298 { 
299   return fMinStep; 
300 } 
301 
302 inline unsigned int G4PathFinder::GetNumberGeometriesLimitingStep() const
303 {
304   unsigned int noGeometries = fNoGeometriesLimiting;
305   return noGeometries;
306 }
307 
308 inline G4double G4PathFinder::GetCurrentSafety() const
309 {
310   return fMinSafety_PreStepPt;
311 }
312 
313 inline void G4PathFinder::MovePoint()
314 {
315   fRelocatedPoint = true;
316 }
317 
318 inline G4Navigator* G4PathFinder::GetNavigator(G4int n) const
319 { 
320   if( (n>fNoActiveNavigators) || (n<0) ) { n=0; }
321   return fpNavigator[n];
322 }
323 
324 inline G4double
325 G4PathFinder::ObtainSafety( G4int navId, G4ThreeVector& globalCenterPoint )
326 {
327   globalCenterPoint = fSafetyLocation; 
328   return fNewSafetyComputed[ navId ];
329 }
330 
331 inline G4double
332 G4PathFinder::LastPreSafety( G4int navId, G4ThreeVector& globalCenterPoint, 
333                              G4double& minSafety )
334 {
335   globalCenterPoint = fPreSafetyLocation;
336   minSafety =         fPreSafetyMinValue;
337   return  fPreSafetyValues[ navId ];
338 }
339 
340 #endif 
341