Geant4 Cross Reference

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