Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 87869 2015-01-16 08:24:36Z 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 
102  void EndTrack();
103  // Signal end of tracking of current track.
104  // Reset internal state
105  // Inform TransportationManager to use 'ordinary' Navigator
106 
108  inline G4VPhysicalVolume* GetLocatedVolume( G4int navId ) const;
109 
111  const G4ThreeVector &pGlobalPoint,
112  const G4ThreeVector &pDirection,
113  const G4double pCurrentProposedStepLength,
114  G4double *prDistance,
115  G4double *prNewSafety= 0)const;
116  // Trial method for checking potential displacement for MS
117 
118  // -----------------------------------------------------------------
119 
120  inline G4bool IsParticleLooping() const;
121 
122  inline G4double GetCurrentSafety() const;
123  // Minimum value of safety after last ComputeStep
124  inline G4double GetMinimumStep() const;
125  // Get the minimum step size from the last ComputeStep call
126  // - in case full step is taken, this is kInfinity
127  inline unsigned int GetNumberGeometriesLimitingStep() const;
128 
129  G4double ComputeSafety( const G4ThreeVector& globalPoint);
130  // Recompute safety for the relevant point the endpoint of the last step!!
131  // Maintain vector of individual safety values (for next method)
132 
133  G4double ObtainSafety( G4int navId, G4ThreeVector& globalCenterPoint );
134  // Obtain safety for navigator/geometry navId for last point 'computed'
135  // --> last point for which ComputeSafety was called
136  // Returns the point (center) for which this safety is valid
137 
138  void EnableParallelNavigation( G4bool enableChoice=true );
139  //
140  // Must call it to ensure that PathFinder is prepared,
141  // especially for curved tracks. If true it switches PropagatorInField
142  // to use MultiNavigator. Must call it with false to undo (=PiF use
143  // Navigator for tracking!)
144 
145  inline G4int SetVerboseLevel(G4int lev=-1);
146 
147  public: // with description
148 
149  inline G4int GetMaxLoopCount() const;
150  inline void SetMaxLoopCount( G4int new_max );
151  //
152  // A maximum for the number of steps that a (looping) particle can take.
153 
154  public: // without description
155 
156  inline void MovePoint();
157  //
158  // Signal that location will be moved -- internal use primarily
159 
160  // To provide best compatibility between Coupled and Old Transportation
161  // the next two methods are provided:
162  G4double LastPreSafety( G4int navId, G4ThreeVector& globalCenterPoint, G4double& minSafety );
163  // Obtain last safety needed in ComputeStep (for geometry navId)
164  // --> last point at which ComputeStep recalculated safety
165  // Returns the point (center) for which this safety is valid
166  // and also the minimum safety over all navigators (ie full)
167 
169  // Tell PathFinder to copy PostStep Safety to PreSafety (for use at next step)
170 
172  // Convert ELimited to string
173 
174  protected: // without description
175 
176  G4double DoNextLinearStep( const G4FieldTrack &FieldTrack,
177  G4double proposedStepLength);
178 
179  G4double DoNextCurvedStep( const G4FieldTrack &FieldTrack,
180  G4double proposedStepLength,
181  G4VPhysicalVolume* pCurrentPhysVolume);
182 
183  void WhichLimited();
184  void PrintLimited();
185  //
186  // Print key details out - for debugging
187 
188  // void ClearState();
189  //
190  // Clear all the State of this class and its current associates
191 
193  //
194  // Whether use safety to discard unneccesary calls to navigator
195 
196  void ReportMove( const G4ThreeVector& OldV, const G4ThreeVector& NewV, const G4String& Quantity ) const;
197  // Helper method to report movement (likely of initial point)
198 
199  protected:
200 
201  G4PathFinder(); // Singleton
202  ~G4PathFinder();
203 
204  inline G4Navigator* GetNavigator(G4int n) const;
205 
206  private:
207 
208  // ----------------------------------------------------------------------
209  // DATA Members
210  // ----------------------------------------------------------------------
211 
212  G4MultiNavigator *fpMultiNavigator;
213  //
214  // Object that enables G4PropagatorInField to see many geometries
215 
216  G4int fNoActiveNavigators;
217  G4bool fNewTrack; // Flag a new track (ensure first step)
218 
219  static const G4int fMaxNav = 16; // rename to kMaxNoNav ??
220 
221  // Global state (retained during stepping for one track)
222 
223  G4Navigator* fpNavigator[fMaxNav];
224 
225  // State changed in a step computation
226 
227  ELimited fLimitedStep[fMaxNav];
228  G4bool fLimitTruth[fMaxNav];
229  G4double fCurrentStepSize[fMaxNav];
230  G4int fNoGeometriesLimiting; // How many processes contribute to limit
231 
232  G4ThreeVector fPreSafetyLocation; // last initial position for which safety evaluated
233  G4double fPreSafetyMinValue; // /\ corresponding value of full safety
234  G4double fPreSafetyValues[ fMaxNav ]; // Safeties for the above point
235  // This part of the state can be retained for severall calls --> CARE
236 
237  G4ThreeVector fPreStepLocation; // point where last ComputeStep called
238  G4double fMinSafety_PreStepPt; // /\ corresponding value of full safety
239  G4double fCurrentPreStepSafety[ fMaxNav ]; // Safeties for the above point
240  // This changes at each step,
241  // so it can differ when steps inside min-safety are made
242 
243  G4bool fPreStepCenterRenewed; // Whether PreSafety coincides with PreStep point
244 
245  G4double fMinStep; // As reported by Navigators -- can be kInfinity
246  G4double fTrueMinStep; // Corrected in case >= proposed
247 
248  // State after calling 'locate'
249 
250  G4VPhysicalVolume* fLocatedVolume[fMaxNav];
251  G4ThreeVector fLastLocatedPosition;
252 
253  // State after calling 'ComputeStep' (others member variables will be affected)
254  G4FieldTrack fEndState; // Point, velocity, ... at proposed step end
255  G4bool fFieldExertedForce; // In current proposed step
256 
257  G4bool fRelocatedPoint; // Signals that point was or is being moved
258  // from the position of the last location
259  // or the endpoint resulting from ComputeStep
260  // -- invalidates fEndState
261 
262  // State for 'ComputeSafety' and related methods
263  G4ThreeVector fSafetyLocation; // point where ComputeSafety is called
264  G4double fMinSafety_atSafLocation; // /\ corresponding value of safety
265  G4double fNewSafetyComputed[ fMaxNav ]; // Safeties for last ComputeSafety
266 
267  // State for Step numbers
268  G4int fLastStepNo, fCurrentStepNo;
269 
270  G4int fVerboseLevel; // For debuging purposes
271 
272  G4TransportationManager* fpTransportManager; // Cache for frequent use
273  G4PropagatorInField* fpFieldPropagator;
274 
275  G4double kCarTolerance;
276 
277  static G4ThreadLocal G4PathFinder* fpPathFinder;
278 };
279 
280 // ********************************************************************
281 // Inline methods.
282 // ********************************************************************
283 
285 {
286  G4VPhysicalVolume* vol=0;
287  if( (navId < fMaxNav) && (navId >=0) ) { vol= fLocatedVolume[navId]; }
288  return vol;
289 }
290 
292 {
293  G4int old= fVerboseLevel; fVerboseLevel= newLevel; return old;
294 }
295 
297 {
298  return fMinStep;
299 }
300 
302 {
303  unsigned int noGeometries=fNoGeometriesLimiting;
304  return noGeometries;
305 }
306 
308 {
309  return fMinSafety_PreStepPt;
310 }
311 
313 {
314  fRelocatedPoint= true;
315 }
316 
318 {
319  if( (n>fNoActiveNavigators)||(n<0)) { n=0; }
320  return fpNavigator[n];
321 }
322 
323 inline G4double G4PathFinder::ObtainSafety( G4int navId, G4ThreeVector& globalCenterPoint )
324 {
325  globalCenterPoint= fSafetyLocation;
326  // navId = std::min( navId, fMaxNav-1 );
327  return fNewSafetyComputed[ navId ];
328 }
329 
331  G4ThreeVector& globalCenterPoint,
332  G4double& minSafety )
333 {
334  globalCenterPoint= fPreSafetyLocation;
335  minSafety= fPreSafetyMinValue;
336  // navId = std::min( navId, fMaxNav-1 );
337  return fPreSafetyValues[ navId ];
338 }
339 #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:89
int G4int
Definition: G4Types.hh:78
void MovePoint()
G4double DoNextCurvedStep(const G4FieldTrack &FieldTrack, G4double proposedStepLength, G4VPhysicalVolume *pCurrentPhysVolume)
unsigned int GetNumberGeometriesLimitingStep() const
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)
const G4int n
void PushPostSafetyToPreSafety()
G4bool RecheckDistanceToCurrentBoundary(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=0) const
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