Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// 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
51class G4Navigator;
52
53#include "G4MultiNavigator.hh"
54#include "G4TouchableHandle.hh"
55
57
59{
60
61 public: // with description
62
63 static G4PathFinder* GetInstance();
64 // Retrieve singleton instance and create it if not existing.
65
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
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
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
155 // Tell PathFinder to copy PostStep Safety to PreSafety
156 // (for use at next step)
157
159 // Convert ELimited to string
160
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
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
284{
285 G4VPhysicalVolume* vol = nullptr;
286 if( (navId < fMaxNav) && (navId >= 0) ) { vol= fLocatedVolume[navId]; }
287 return vol;
288}
289
291{
292 G4int old = fVerboseLevel;
293 fVerboseLevel = newLevel;
294 return old;
295}
296
298{
299 return fMinStep;
300}
301
303{
304 unsigned int noGeometries = fNoGeometriesLimiting;
305 return noGeometries;
306}
307
309{
310 return fMinSafety_PreStepPt;
311}
312
314{
315 fRelocatedPoint = true;
316}
317
319{
320 if( (n>fNoActiveNavigators) || (n<0) ) { n=0; }
321 return fpNavigator[n];
322}
323
324inline G4double
326{
327 globalCenterPoint = fSafetyLocation;
328 return fNewSafetyComputed[ navId ];
329}
330
331inline G4double
333 G4double& minSafety )
334{
335 globalCenterPoint = fPreSafetyLocation;
336 minSafety = fPreSafetyMinValue;
337 return fPreSafetyValues[ navId ];
338}
339
340#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double ComputeSafety(const G4ThreeVector &globalPoint)
void EnableParallelNavigation(G4bool enableChoice=true)
G4Navigator * GetNavigator(G4int n) const
void ReLocate(const G4ThreeVector &position)
void ReportMove(const G4ThreeVector &OldV, const G4ThreeVector &NewV, const G4String &Quantity) const
G4double DoNextLinearStep(const G4FieldTrack &FieldTrack, G4double proposedStepLength)
G4bool IsParticleLooping() const
unsigned int GetNumberGeometriesLimitingStep() const
void PushPostSafetyToPreSafety()
G4double ObtainSafety(G4int navId, G4ThreeVector &globalCenterPoint)
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4double DoNextCurvedStep(const G4FieldTrack &FieldTrack, G4double proposedStepLength, G4VPhysicalVolume *pCurrentPhysVolume)
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)
G4bool UseSafetyForOptimization(G4bool)
G4double LastPreSafety(G4int navId, G4ThreeVector &globalCenterPoint, G4double &minSafety)
G4String & LimitedString(ELimited lim)
G4double GetCurrentSafety() const
G4int GetMaxLoopCount() const
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
void SetMaxLoopCount(G4int new_max)
G4double GetMinimumStep() const
static G4PathFinder * GetInstance()
G4int SetVerboseLevel(G4int lev=-1)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
static G4PathFinder * GetInstanceIfExist()
#define G4ThreadLocal
Definition tls.hh:77