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