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