Geant4 10.7.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// 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 "G4TouchableHandle.hh"
54#include "G4FieldTrack.hh"
55#include "G4MultiNavigator.hh"
56
58
60{
61
62 public: // with description
63
64 static G4PathFinder* GetInstance();
65 // Retrieve singleton instance and create it if not existing.
66
68 // Retrieve singleton instance pointer.
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 // Compute the next geometric Step -- curved or linear
79 // If it is called with a larger 'stepNo' it will execute a new step;
80 // if 'stepNo' is same as last call, then the results for
81 // the geometry with Id. number 'navigatorId' will be returned.
82
83 void Locate( const G4ThreeVector& position,
84 const G4ThreeVector& direction,
85 G4bool relativeSearch = true);
86 // Make primary relocation of global point in all navigators,
87 // and update them.
88
89 void ReLocate( const G4ThreeVector& position );
90 // Make secondary relocation of global point (within safety only)
91 // in all navigators, and update them.
92
94 const G4ThreeVector& direction,
95 G4VPhysicalVolume* massStartVol = nullptr);
96 // Check and cache set of active navigators.
97
98 void EndTrack();
99 // Signal end of tracking of current track.
100 // Reset internal state
101 // Inform TransportationManager to use 'ordinary' Navigator
102
104 inline G4VPhysicalVolume* GetLocatedVolume( G4int navId ) const;
105
107 const G4ThreeVector& pGlobalPoint,
108 const G4ThreeVector& pDirection,
109 const G4double pCurrentProposedStepLength,
110 G4double* prDistance,
111 G4double* prNewSafety = nullptr) const;
112 // Trial method for checking potential displacement for MS
113
114 // -----------------------------------------------------------------
115
116 inline G4bool IsParticleLooping() const;
117
118 inline G4double GetCurrentSafety() const;
119 // Minimum value of safety after last ComputeStep
120 inline G4double GetMinimumStep() const;
121 // Get the minimum step size from the last ComputeStep call
122 // - in case full step is taken, this is kInfinity
123 inline unsigned int GetNumberGeometriesLimitingStep() const;
124
125 G4double ComputeSafety( const G4ThreeVector& globalPoint);
126 // Recompute safety for the relevant point the endpoint of the last step!!
127 // Maintain vector of individual safety values (for next method)
128
129 G4double ObtainSafety( G4int navId, G4ThreeVector& globalCenterPoint );
130 // Obtain safety for navigator/geometry navId for last point 'computed'
131 // --> last point for which ComputeSafety was called
132 // Returns the point (center) for which this safety is valid
133
134 void EnableParallelNavigation( G4bool enableChoice = true );
135 // Must call it to ensure that PathFinder is prepared,
136 // especially for curved tracks. If true it switches PropagatorInField
137 // to use MultiNavigator. Must call it with false to undo (=PiF use
138 // Navigator for tracking!)
139
140 inline G4int SetVerboseLevel(G4int lev = -1);
141
142 public: // with description
143
144 inline G4int GetMaxLoopCount() const;
145 inline void SetMaxLoopCount( G4int new_max );
146 // A maximum for the number of steps that a (looping) particle can take
147
148 public: // without description
149
150 inline void MovePoint();
151 // Signal that location will be moved -- internal use primarily
152
153 // To provide best compatibility between Coupled and Old Transportation
154 // the next two methods are provided:
155
156 G4double LastPreSafety( G4int navId, G4ThreeVector& globalCenterPoint,
157 G4double& minSafety );
158 // Obtain last safety needed in ComputeStep (for geometry navId)
159 // --> last point at which ComputeStep recalculated safety
160 // Returns the point (center) for which this safety is valid
161 // and also the minimum safety over all navigators (i.e. full)
162
164 // Tell PathFinder to copy PostStep Safety to PreSafety
165 // (for use at next step)
166
168 // Convert ELimited to string
169
171 // Destructor
172
173 protected: // without description
174
175 G4double DoNextLinearStep( const G4FieldTrack& FieldTrack,
176 G4double proposedStepLength);
177
178 G4double DoNextCurvedStep( const G4FieldTrack& FieldTrack,
179 G4double proposedStepLength,
180 G4VPhysicalVolume* pCurrentPhysVolume);
181
182 void WhichLimited();
183 void PrintLimited();
184 // Print key details out for debugging
185
187 // Whether use safety to discard unneccesary calls to navigator
188
189 void ReportMove( const G4ThreeVector& OldV,
190 const G4ThreeVector& NewV,
191 const G4String& Quantity ) const;
192 // Helper method to report movement (likely of initial point)
193
194 protected:
195
196 G4PathFinder(); // Singleton
197
198 inline G4Navigator* GetNavigator(G4int n) const;
199
200 private:
201
202 // ----------------------------------------------------------------------
203 // DATA Members
204 // ----------------------------------------------------------------------
205
206 G4MultiNavigator* fpMultiNavigator;
207 // Object that enables G4PropagatorInField to see many geometries
208
209 G4int fNoActiveNavigators = 0;
210 G4bool fNewTrack = false; // Flag a new track (ensure first step)
211
212 static const G4int fMaxNav = 16;
213
214 // Global state (retained during stepping for one track)
215
216 G4Navigator* fpNavigator[fMaxNav];
217
218 // State changed in a step computation
219
220 ELimited fLimitedStep[fMaxNav];
221 G4bool fLimitTruth[fMaxNav];
222 G4double fCurrentStepSize[fMaxNav];
223 G4int fNoGeometriesLimiting = 0; // How many processes contribute to limit
224
225 G4ThreeVector fPreSafetyLocation;
226 // last initial position for which safety evaluated
227 G4double fPreSafetyMinValue = -1.0;
228 // - corresponding value of full safety
229 G4double fPreSafetyValues[ fMaxNav ];
230 // Safeties for the above point
231
232 // This part of the state can be retained for several calls --> CARE
233
234 G4ThreeVector fPreStepLocation;
235 // point where last ComputeStep called
236 G4double fMinSafety_PreStepPt = -1.0;
237 // - corresponding value of full safety
238 G4double fCurrentPreStepSafety[ fMaxNav ];
239 // Safeties for the above point.
240 // This changes at each step, so it can differ when steps
241 // inside min-safety are made
242
243 G4bool fPreStepCenterRenewed = false;
244 // Whether PreSafety coincides with PreStep point
245
246 G4double fMinStep = -1.0; // As reported by Navigators -- can be kInfinity
247 G4double fTrueMinStep = -1.0; // Corrected in case >= proposed
248
249 // State after calling 'locate'
250 //
251 G4VPhysicalVolume* fLocatedVolume[fMaxNav];
252 G4ThreeVector fLastLocatedPosition;
253
254 // State after calling 'ComputeStep'
255 // (others member variables will be affected)
256 //
257 G4FieldTrack fEndState; // Point, velocity, ... at proposed step end
258 G4bool fFieldExertedForce = false; // In current proposed step
259
260 G4bool fRelocatedPoint = false; // Signals that point was or is being moved
261 // from the position of the last location or
262 // the endpoint resulting from ComputeStep()
263 // -- invalidates fEndState
264
265 // State for 'ComputeSafety' and related methods
266 //
267 G4ThreeVector fSafetyLocation;
268 // point where ComputeSafety is called
269 G4double fMinSafety_atSafLocation = -1.0;
270 // - corresponding value of safety
271 G4double fNewSafetyComputed[ fMaxNav ];
272 // Safeties for last ComputeSafety
273
274 // State for Step numbers
275 //
276 G4int fLastStepNo = -1, fCurrentStepNo = -1;
277
278 G4int fVerboseLevel = 0; // For debugging purposes
279
280 G4TransportationManager* fpTransportManager; // Cache for frequent use
281 G4PropagatorInField* fpFieldPropagator;
282
283 G4double kCarTolerance;
284
285 static G4ThreadLocal G4PathFinder* fpPathFinder;
286};
287
288// ********************************************************************
289// Inline methods.
290// ********************************************************************
291
293{
294 G4VPhysicalVolume* vol = nullptr;
295 if( (navId < fMaxNav) && (navId >= 0) ) { vol= fLocatedVolume[navId]; }
296 return vol;
297}
298
300{
301 G4int old = fVerboseLevel;
302 fVerboseLevel = newLevel;
303 return old;
304}
305
307{
308 return fMinStep;
309}
310
312{
313 unsigned int noGeometries = fNoGeometriesLimiting;
314 return noGeometries;
315}
316
318{
319 return fMinSafety_PreStepPt;
320}
321
323{
324 fRelocatedPoint = true;
325}
326
328{
329 if( (n>fNoActiveNavigators) || (n<0) ) { n=0; }
330 return fpNavigator[n];
331}
332
333inline G4double
335{
336 globalCenterPoint = fSafetyLocation;
337 return fNewSafetyComputed[ navId ];
338}
339
340inline G4double
342 G4double& minSafety )
343{
344 globalCenterPoint = fPreSafetyLocation;
345 minSafety = fPreSafetyMinValue;
346 return fPreSafetyValues[ navId ];
347}
348
349#endif
ELimited
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
void WhichLimited()
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)
void MovePoint()
G4bool UseSafetyForOptimization(G4bool)
G4double LastPreSafety(G4int navId, G4ThreeVector &globalCenterPoint, G4double &minSafety)
G4String & LimitedString(ELimited lim)
G4double GetCurrentSafety() const
void PrintLimited()
G4bool RecheckDistanceToCurrentBoundary(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=nullptr) 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()
Definition: G4PathFinder.cc:52
G4int SetVerboseLevel(G4int lev=-1)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
static G4PathFinder * GetInstanceIfExist()
Definition: G4PathFinder.cc:66
#define G4ThreadLocal
Definition: tls.hh:77