source: trunk/source/geometry/navigation/include/G4PropagatorInField.hh@ 883

Last change on this file since 883 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 12.1 KB
RevLine 
[831]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//
[850]27// $Id: G4PropagatorInField.hh,v 1.14 2008/05/28 09:11:59 tnikitin Exp $
28// GEANT4 tag $Name: HEAD $
[831]29//
30// class G4PropagatorInField
31//
32// Class description:
33//
34// This class performs the navigation/propagation of a particle/track
35// in a magnetic field. The field is in general non-uniform.
36// For the calculation of the path, it relies on the class G4ChordFinder.
37//
38// Key Method:
39// ComputeStep(..)
40// History:
41// -------
42// 25.10.96 John Apostolakis, design and implementation
43// 25.03.97 John Apostolakis, adaptation for G4Transportation and cleanup
44// 8.11.02 John Apostolakis, changes to enable use of safety in intersecting
45// ---------------------------------------------------------------------------
46
47#ifndef G4PropagatorInField_hh
48#define G4PropagatorInField_hh 1
49
50#include "G4Types.hh"
51
52#include <vector>
53
54#include "G4FieldTrack.hh"
55#include "G4FieldManager.hh"
56class G4ChordFinder;
57
58class G4Navigator;
59class G4VPhysicalVolume;
60class G4VCurvedTrajectoryFilter;
61
62class G4PropagatorInField
63{
64
65 public: // with description
66
67 G4PropagatorInField( G4Navigator *theNavigator,
68 G4FieldManager *detectorFieldMgr );
69 ~G4PropagatorInField();
70
71 G4double ComputeStep( G4FieldTrack &pFieldTrack,
72 G4double pCurrentProposedStepLength,
73 G4double &pNewSafety,
74 G4VPhysicalVolume *pPhysVol=0 );
75 // Compute the next geometric Step
76
77 inline G4ThreeVector EndPosition() const;
78 inline G4ThreeVector EndMomentumDir() const;
79 inline G4bool IsParticleLooping() const;
80 // Return the state after the Step
81
82 inline G4double GetEpsilonStep() const;
83 // Relative accuracy for current Step (Calc.)
84 inline void SetEpsilonStep(G4double newEps);
85 // The ratio DeltaOneStep()/h_current_step
86
87 inline void SetChargeMomentumMass( G4double charge, // in e+ units
88 G4double momentum, // in Geant4 units
89 G4double pMass );
90 // Inform this and all associated classes of q, p, m
91
92 G4FieldManager* FindAndSetFieldManager(G4VPhysicalVolume* pCurrentPhysVol);
93 // Set (and return) the correct field manager (global or local),
94 // if it exists.
95 // Should be called before ComputeStep is called;
96 // - currently, ComputeStep will call it, if it has not been called.
97
98 inline G4ChordFinder* GetChordFinder();
99
100 G4int SetVerboseLevel( G4int verbose );
101 inline G4int GetVerboseLevel() const;
102 inline G4int Verbose() const;
103
104 inline G4int GetMaxLoopCount() const;
105 inline void SetMaxLoopCount( G4int new_max );
106 // A maximum for the number of steps that a (looping) particle can take.
107
108 void printStatus( const G4FieldTrack& startFT,
109 const G4FieldTrack& currentFT,
110 G4double requestStep,
111 G4double safety,
112 G4int step,
113 G4VPhysicalVolume* startVolume);
114 // Print Method - useful mostly for debugging.
115
116 inline G4FieldTrack GetEndState() const;
117
118 // The following methods are now obsolescent but *for now* will work
119 // They are being replaced by same-name methods in G4FieldManager,
120 // allowing the specialisation in different volumes.
121 // Their new behaviour is to change the values for the global field manager.
122 inline G4double GetMinimumEpsilonStep() const;
123 inline void SetMinimumEpsilonStep( G4double newEpsMin );
124 // Minimum for Relative accuracy of any Step
125
126 inline G4double GetMaximumEpsilonStep() const;
127 inline void SetMaximumEpsilonStep( G4double newEpsMax );
128
129 inline void SetLargestAcceptableStep( G4double newBigDist );
130 inline G4double GetLargestAcceptableStep();
[850]131
132 // Use alternative Locator(based on Brent Method,second order Intersection)
133 inline void SetBrentMethod(G4bool newLocator);
134 inline G4bool GetBrentMethod();
[831]135
136 public: // with description
137
138 // The following methods are obsolete and will not work --
139 // as they have been replaced by the same methods in G4FieldManager
140 // since Geant4 4.0
141 inline G4double GetDeltaIntersection() const;
142 inline G4double GetDeltaOneStep() const;
143 inline void SetAccuraciesWithDeltaOneStep( G4double deltaOneStep );
144 inline void SetDeltaIntersection( G4double deltaIntersection );
145 inline void SetDeltaOneStep( G4double deltaOneStep );
146
147 public: // without description
148
149 inline G4FieldManager* GetCurrentFieldManager();
150 inline void SetNavigatorForPropagating( G4Navigator *SimpleOrMultiNavigator );
151 inline G4Navigator* GetNavigatorForPropagating();
152
153 public: // no description
154
155 inline void SetThresholdNoZeroStep( G4int noAct,
156 G4int noHarsh,
157 G4int noAbandon );
158 inline G4int GetThresholdNoZeroSteps( G4int i );
159
160 public: // with description
161 //
162 void SetTrajectoryFilter(G4VCurvedTrajectoryFilter* filter);
163 // Set the filter that examines & stores 'intermediate'
164 // curved trajectory points. Currently only position is stored.
165
166 std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
167 // Access the points which have passed by the filter.
168 // Responsibility for deleting the points lies with the client.
169 // This method MUST BE called exactly ONCE per step.
170
171 void ClearPropagatorState();
172 // Clear all the State of this class and its current associates
173 // --> the current field manager & chord finder will also be called
174
175 inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
176 // Update this (dangerous) state -- for the time being
177
178 inline void SetUseSafetyForOptimization( G4bool );
179 inline G4bool GetUseSafetyForOptimization();
180 // Toggle & view parameter for using safety to discard
181 // unneccesary calls to navigator (thus 'optimising' performance)
182
183 protected: // with description
184
185 G4bool LocateIntersectionPoint(
186 const G4FieldTrack& curveStartPointTangent, // A
187 const G4FieldTrack& curveEndPointTangent, // B
188 const G4ThreeVector& trialPoint, // E
189 G4FieldTrack& intersectPointTangent, // Output
190 G4bool& recalculatedEndPoint); // Out:
191
192 // If such an intersection exists, this function
193 // calculate the intersection point of the true path of the particle
194 // with the surface of the current volume (or of one of its daughters).
195 // (Should use lateral displacement as measure of convergence).
196
197 G4bool IntersectChord( G4ThreeVector StartPointA,
198 G4ThreeVector EndPointB,
199 G4double &NewSafety,
200 G4double &LinearStepLength,
201 G4ThreeVector &IntersectionPoint);
202 // Intersect the chord from StartPointA to EndPointB
203 // and return whether an intersection occurred
204
205 G4FieldTrack ReEstimateEndpoint( const G4FieldTrack &CurrentStateA,
206 const G4FieldTrack &EstimtdEndStateB,
207 G4double linearDistSq,
208 G4double curveDist);
209 // Return new estimate for state after curveDist
210 // starting from CurrentStateA, to replace EstimtdEndStateB,
211 // (and report displacement -- if field is compiled verbose.)
212
213 void PrintStepLengthDiagnostic( G4double currentProposedStepLength,
214 G4double decreaseFactor,
215 G4double stepTrial,
216 const G4FieldTrack& aFieldTrack);
217 private:
218
219 // ----------------------------------------------------------------------
220 // DATA Members
221 // ----------------------------------------------------------------------
222
223 G4FieldManager *fDetectorFieldMgr;
224 // The Field Manager of the whole Detector. (default)
225
226 G4FieldManager *fCurrentFieldMgr;
227 // The Field Manager of the current volume (may be the one above.)
228
229 G4Navigator *fNavigator;
230
231 // STATE information
232 // -----------------
233
234 G4double fEpsilonStep;
235 // Relative accuracy for current Step (Calc.)
236
237 G4FieldTrack End_PointAndTangent;
238 // End point storage
239
240 G4bool fParticleIsLooping;
241
242 G4int fVerboseLevel;
243 // For debuging purposes
244
245 G4int fMax_loop_count;
246 // Limit for the number of sub-steps taken in one call to ComputeStep
247
248 // Variables to keep track of "abnormal" case - which causes loop
249 //
250 G4int fNoZeroStep; // Counter of zeroStep
251 G4int fActionThreshold_NoZeroSteps; // Threshold: above this - act
252 G4int fSevereActionThreshold_NoZeroSteps; // Threshold to act harshly
253 G4int fAbandonThreshold_NoZeroSteps; // Threshold to abandon
254
255 G4double fFull_CurveLen_of_LastAttempt;
256 G4double fLast_ProposedStepLength;
257 G4double fLargestAcceptableStep;
258
259 G4double fCharge, fInitialMomentumModulus, fMass;
260
261 G4ThreeVector fPreviousSftOrigin;
262 G4double fPreviousSafety;
263 G4bool fUseSafetyForOptimisation;
264 // Last safety origin & value: for optimisation
265
266 G4bool fSetFieldMgr;
267 // Flag whether field manager has been set for the current step
268
269 G4double kCarTolerance;
270 // Geometrical tolerance defining surface thickness
271
[850]272 G4int maxNumberOfStepsForIntersection;
273 G4int maxNumberOfCallsToReIntegration;
274 G4int maxNumberOfCallsToReIntegration_depth;
275 // Counters for Statistics about Location and ReIntegrations
[831]276private:
277
278 static const G4int max_depth=4;
279 G4FieldTrack* ptrInterMedFT[max_depth+1];
280 // Used to store intermediate values of tracks in case of
281 // too slow progress
282private:
[850]283 G4bool fUseBrentLocator;
284
285private:
[831]286
287 G4VCurvedTrajectoryFilter* fpTrajectoryFilter;
288 // The filter encapsulates the algorithm which selects which
289 // intermediate points should be stored in a trajectory.
290 // When it is NULL, no intermediate points will be stored.
291 // Else PIF::ComputeStep must submit (all) intermediate
292 // points it calculates, to this filter. (jacek 04/11/2002)
293};
294
295// ********************************************************************
296// Inline methods.
297// ********************************************************************
298
299#include "G4PropagatorInField.icc"
300
301#endif
Note: See TracBrowser for help on using the repository browser.