Geant4 11.4.0
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// 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// In field, it relies on the class G4PropagatorInField.
36
37// Author: John Apostolakis (CERN), 7 October 2005
38// ---------------------------------------------------------------------
39#ifndef G4PATHFINDER_HH
40#define G4PATHFINDER_HH 1
41
42#include <vector>
43
44#include "G4Types.hh"
45#include "G4FieldTrack.hh"
46#include "G4MultiNavigator.hh"
47#include "G4TouchableHandle.hh"
48
50class G4Navigator;
52
53/**
54 * @brief G4PathFinder directs the lock-stepped propagation of a track in the
55 * 'mass' and other parallel geometries. It ensures that tracking in a magnetic
56 * field sees these parallel geometries at each trial step and that the earliest
57 * boundary limits the step.
58 */
59
60class G4PathFinder
61{
62 public:
63
64 /**
65 * Retrieves singleton instance and creates it if not existing.
66 */
67 static G4PathFinder* GetInstance();
68
69 /**
70 * Retrieve singleton instance pointer.
71 */
72 static G4PathFinder* GetInstanceIfExist();
73
74 /**
75 * Destructor, called only by G4RunManagerKernel.
76 */
78
79 /**
80 * Computes the next geometric Step, curved or linear.
81 * If it is called with a larger 'stepNo' it will execute a new step;
82 * if 'stepNo' is same as last call, then the results for the geometry
83 * with Id number 'navigatorId' will be returned.
84 * @param[in,out] pFieldTrack Field track to be filled.
85 * @param[in] pCurrentProposedStepLength Current proposed step length.
86 * @param[in] navigatorId Identifier of the geometry.
87 * @param[in] stepNo Step number; see next step/check.
88 * @param[in, out] pNewSafety New safety for this geometry.
89 * @param[in, out] limitedStep Step characterisation to be returned.
90 * @param[in, out] EndState Field track end state.
91 * @param[in] currentVolume Pointer to the current volume.
92 * @returns Step length.
93 */
94 G4double ComputeStep( const G4FieldTrack& pFieldTrack,
95 G4double pCurrentProposedStepLength,
96 G4int navigatorId, // Identifies the geometry
97 G4int stepNo, // See next step/check
98 G4double& pNewSafety, // Only for this geometry
99 ELimited& limitedStep,
100 G4FieldTrack& EndState,
101 G4VPhysicalVolume* currentVolume );
102
103 /**
104 * Makes primary relocation of the global point in all navigators,
105 * and updates them.
106 * @param[in] position Point in global coordinates system.
107 * @param[in] direction Global direction vector.
108 * @param[in] relativeSearch If set to true (default), the search begins
109 * is the geometrical hierarchy at the location of the last
110 * located point.
111 */
112 void Locate( const G4ThreeVector& position,
113 const G4ThreeVector& direction,
114 G4bool relativeSearch = true);
115
116 /**
117 * Makes secondary relocation of the global point (within safety only)
118 * in all navigators, and updates them.
119 * @param[in] position Point in global coordinates system.
120 */
121 void ReLocate( const G4ThreeVector& position );
122
123 /**
124 * Checks and caches the set of active navigators.
125 * @param[in] position Point in global coordinates system.
126 * @param[in] direction Global direction vector.
127 * @param[in] massStartVol Pointer to the mass geometry world.
128 */
130 const G4ThreeVector& direction,
131 G4VPhysicalVolume* massStartVol = nullptr);
132
133 /**
134 * Signals the end of tracking of the current track. Resets internal state
135 * and informs G4TransportationManager to use 'ordinary' Navigator.
136 */
137 void EndTrack();
138
139 /**
140 * Creates a touchable handle for the specified navigator.
141 * @param[in] navId The navigator identifier.
142 * @returns A touchable handle of the geometry.
143 */
145
146 /**
147 * Returns the located volume for the specified navigator.
148 * @param[in] navId The navigator identifier.
149 * @returns A pointer to the located volume in the geometry.
150 */
152
153 // -----------------------------------------------------------------
154
155 /**
156 * Returns the minimum value of safety after last ComputeStep().
157 */
159
160 /**
161 * Gets the minimum step size from the last ComputeStep() call.
162 * @note In case full step is taken, this is kInfinity.
163 */
164 inline G4double GetMinimumStep() const;
165
166 /**
167 * Returns the number of all geometries limiting the step.
168 */
169 inline unsigned int GetNumberGeometriesLimitingStep() const;
170
171 /**
172 * Recomputes the safety for the relevant point, i.e. the endpoint of the
173 * last step. Maintains a vector of individual safety values (used by next
174 * method below).
175 * @param[in] globalPoint Point in global coordinates system.
176 * @returns The safety value for the specified point in the geometry.
177 */
178 G4double ComputeSafety( const G4ThreeVector& globalPoint);
179
180 /**
181 * Obtains the safety for the specified navigator/geometry for last point
182 * 'computed' (i.e., the last point for which ComputeSafety() was called).
183 * @param[in] navId The navigator identifier.
184 * @param[in,out] globalCenterPoint The point (center) for which this
185 * safety is valid.
186 * @returns The safety value in the specified geometry.
187 */
188 inline G4double ObtainSafety(G4int navId, G4ThreeVector& globalCenterPoint);
189
190 /**
191 * To enable parallel navigation. Must call it to ensure that G4PathFinder
192 * is prepared, especially for curved tracks.
193 * If true it switches G4PropagatorInField to use G4MultiNavigator.
194 * Must call it with false to undo (i.e. G4PropagatorInField uses classic
195 * G4Navigator for tracking in such case).
196 * @param[in] enableChoice Flag to enable/disable parallel navigation.
197 */
198 void EnableParallelNavigation( G4bool enableChoice = true );
199
200 /**
201 * To control the level of verbosity. Default is no verbosity.
202 */
203 inline G4int SetVerboseLevel(G4int lev = -1);
204
205 /**
206 * To get/set the maximum for the number of steps that a (looping)
207 * particle can take.
208 */
209 inline G4int GetMaxLoopCount() const;
210 inline void SetMaxLoopCount( G4int new_max );
211
212 /**
213 * Signals that the point location will be moved.
214 * @note Internal use primarily.
215 */
216 inline void MovePoint();
217
218 // To provide best compatibility between Coupled and normal Transportation
219 // the next two methods are provided...
220
221 /**
222 * Obtains the last safety needed in ComputeStep() for the specified
223 * geometry 'navId' (i.e. the last point at which ComputeStep() has
224 * recalculated the safety). Returns the point (center) for which this
225 * safety is valid and also the minimum safety over all navigators.
226 * @param[in] navId The navigator identifier.
227 * @param[in,out] globCenterPoint The point (center) for which this
228 * safety is valid.
229 * @param[in,out] minSafety The minimum safety over all navigators.
230 * @returns The safety value in the specified geometry.
231 */
232 inline G4double LastPreSafety( G4int navId, G4ThreeVector& globCenterPoint,
233 G4double& minSafety );
234
235 /**
236 * Tells G4PathFinder to copy PostStep Safety to PreSafety
237 * for use at the next step.
238 */
240
241 /**
242 * Utility to convert ELimited specification to a string.
243 */
245
246 private:
247
248 /**
249 * Private singleton constructor.
250 */
251 G4PathFinder();
252
253 /**
254 * Returns pointer to the specified navigator.
255 */
256 inline G4Navigator* GetNavigator(G4int n) const;
257
258 /**
259 * Performs a linear step.
260 * @param[in,out] FieldTrack Field track to be filled.
261 * @param[in] proposedStepLength Current proposed step length.
262 * @returns The minimum linear step to undertake.
263 */
264 G4double DoNextLinearStep( const G4FieldTrack& FieldTrack,
265 G4double proposedStepLength);
266
267 /**
268 * Performs a curved step.
269 * @param[in,out] FieldTrack Field track to be filled.
270 * @param[in] proposedStepLength Current proposed step length.
271 * @param[in] pCurrentPhysVolume Pointer to the current volume of interest.
272 * @returns The minimum step to undertake.
273 */
274 G4double DoNextCurvedStep( const G4FieldTrack& FieldTrack,
275 G4double proposedStepLength,
276 G4VPhysicalVolume* pCurrentPhysVolume);
277
278 /**
279 * Prints key details out for debugging.
280 */
281 void WhichLimited();
282 void PrintLimited();
283
284 /**
285 * Helper method to report movement (likely of initial point).
286 */
287 void ReportMove( const G4ThreeVector& OldV,
288 const G4ThreeVector& NewV,
289 const G4String& Quantity ) const;
290
291 private:
292
293 // ----------------------------------------------------------------------
294 // DATA Members
295 // ----------------------------------------------------------------------
296
297 /** Object that enables G4PropagatorInField to see many geometries. */
298 G4MultiNavigator* fpMultiNavigator;
299
300 G4int fNoActiveNavigators = 0;
301 G4bool fNewTrack = false; // Flag a new track (ensure first step)
302
303 static const G4int fMaxNav = 16;
304
305 // Global state (retained during stepping for one track)
306
307 G4Navigator* fpNavigator[fMaxNav];
308
309 // ---- State changed in a step computation
310 //
311 ELimited fLimitedStep[fMaxNav];
312 G4bool fLimitTruth[fMaxNav];
313 G4double fCurrentStepSize[fMaxNav];
314 G4int fNoGeometriesLimiting = 0; // How many processes contribute to limit
315
316 /** Last initial position for which safety evaluated. */
317 G4ThreeVector fPreSafetyLocation;
318 /* Corresponding value of full safety. */
319 G4double fPreSafetyMinValue = -1.0;
320
321 /** Safeties for the above point. */
322 G4double fPreSafetyValues[ fMaxNav ];
323
324 // This part of the state can be retained for several calls --> CARE
325
326 /** Point where last ComputeStep() called. */
327 G4ThreeVector fPreStepLocation;
328 /** Corresponding value of full safety. */
329 G4double fMinSafety_PreStepPt = -1.0;
330
331 /** Safeties for the above point.
332 * @note This changes at each step, so it can differ when steps
333 * inside min-safety are made. */
334 G4double fCurrentPreStepSafety[ fMaxNav ];
335
336 /** Whether PreSafety coincides with PreStep point. */
337 G4bool fPreStepCenterRenewed = false;
338
339 G4double fMinStep = -1.0; // As reported by Navigators -- can be kInfinity
340 G4double fTrueMinStep = -1.0; // Corrected in case >= proposed
341
342 // ---- State after calling 'locate'
343 //
344 G4VPhysicalVolume* fLocatedVolume[fMaxNav];
345 G4ThreeVector fLastLocatedPosition;
346
347 // ---- State after calling 'ComputeStep'
348 // (other member variables will be affected)
349 //
350 G4FieldTrack fEndState; // Point, velocity, ... at proposed step end
351 G4bool fFieldExertedForce = false; // In current proposed step
352
353 G4bool fRelocatedPoint = false; // Signals that point was or is being moved
354 // from the position of the last location or
355 // the endpoint resulting from ComputeStep()
356 // -- invalidates fEndState
357
358 // ---- State for 'ComputeSafety' and related methods
359 //
360 /** Point where ComputeSafety() is called. */
361 G4ThreeVector fSafetyLocation;
362 /** Corresponding value of safety. */
363 G4double fMinSafety_atSafLocation = -1.0;
364 /** Safeties for last ComputeSafety(). */
365 G4double fNewSafetyComputed[ fMaxNav ];
366
367 // ---- State for Step numbers
368 //
369 G4int fLastStepNo = -1, fCurrentStepNo = -1;
370
371 G4int fVerboseLevel = 0; // For debugging purposes
372
373 G4TransportationManager* fpTransportManager; // Cache for frequent use
374 G4PropagatorInField* fpFieldPropagator;
375
376 G4double kCarTolerance;
377
378 static G4ThreadLocal G4PathFinder* fpPathFinder;
379};
380
381// ********************************************************************
382// Inline methods.
383// ********************************************************************
384
385#include "G4PathFinder.icc"
386
387#endif
CLHEP::Hep3Vector G4ThreeVector
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
G4TouchableHandle is a type providing reference counting mechanism for any kind of touchable objects....
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
G4MultiNavigator is a utility class for polling the navigators of several geometries to identify the ...
G4Navigator is a class for use by the tracking management, able to obtain/calculate dynamic tracking ...
G4double ComputeSafety(const G4ThreeVector &globalPoint)
void EnableParallelNavigation(G4bool enableChoice=true)
void ReLocate(const G4ThreeVector &position)
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)
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)
void MovePoint()
G4String & LimitedString(ELimited lim)
G4double GetCurrentSafety() const
G4double LastPreSafety(G4int navId, G4ThreeVector &globCenterPoint, G4double &minSafety)
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()
G4PropagatorInField performs the navigation/propagation of a particle/track in a magnetic field....
G4TransportationManager is a singleton class which stores the navigator used by the transportation pr...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
#define G4ThreadLocal
Definition tls.hh:77