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

Last change on this file since 831 was 831, checked in by garnier, 16 years ago

import all except CVS

File size: 11.7 KB
Line 
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//
27// $Id: G4PropagatorInField.hh,v 1.13 2007/06/08 09:49:34 gcosmo Exp $
28// GEANT4 tag $Name:  $
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();
131
132 public:  // with description
133
134   // The following methods are obsolete and will not work --
135   //   as they have been replaced by the same methods in G4FieldManager
136   //   since Geant4 4.0
137   inline G4double  GetDeltaIntersection() const;
138   inline G4double  GetDeltaOneStep() const;
139   inline void    SetAccuraciesWithDeltaOneStep( G4double deltaOneStep ); 
140   inline void    SetDeltaIntersection( G4double deltaIntersection );
141   inline void    SetDeltaOneStep( G4double deltaOneStep ); 
142
143 public:  // without description
144
145   inline G4FieldManager*  GetCurrentFieldManager();
146   inline void             SetNavigatorForPropagating( G4Navigator *SimpleOrMultiNavigator ); 
147   inline G4Navigator*     GetNavigatorForPropagating(); 
148
149 public:  // no description
150
151   inline void SetThresholdNoZeroStep( G4int noAct,
152                                       G4int noHarsh,
153                                       G4int noAbandon );
154   inline G4int GetThresholdNoZeroSteps( G4int i ); 
155
156 public:  // with description
157  //
158  void SetTrajectoryFilter(G4VCurvedTrajectoryFilter* filter);
159  // Set the filter that examines & stores 'intermediate'
160  //  curved trajectory points.  Currently only position is stored.
161
162  std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
163  // Access the points which have passed by the filter.
164  // Responsibility for deleting the points lies with the client.
165  // This method MUST BE called exactly ONCE per step.
166
167  void ClearPropagatorState();
168  // Clear all the State of this class and its current associates
169  //   --> the current field manager & chord finder will also be called
170
171  inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
172      // Update this (dangerous) state -- for the time being
173 
174  inline void   SetUseSafetyForOptimization( G4bool );
175  inline G4bool GetUseSafetyForOptimization();
176      // Toggle & view parameter for using safety to discard
177      //   unneccesary calls to navigator (thus 'optimising' performance)
178
179 protected:  // with description
180
181   G4bool LocateIntersectionPoint( 
182        const  G4FieldTrack&       curveStartPointTangent,  //  A
183        const  G4FieldTrack&       curveEndPointTangent,    //  B
184        const  G4ThreeVector&      trialPoint,              //  E
185               G4FieldTrack&       intersectPointTangent,   // Output
186               G4bool&             recalculatedEndPoint);   // Out:
187
188     // If such an intersection exists, this function
189     // calculate the intersection point of the true path of the particle
190     // with the surface of the current volume (or of one of its daughters).
191     // (Should use lateral displacement as measure of convergence).
192
193   G4bool IntersectChord( G4ThreeVector  StartPointA, 
194                          G4ThreeVector  EndPointB,
195                          G4double      &NewSafety,
196                          G4double      &LinearStepLength,
197                          G4ThreeVector &IntersectionPoint);
198     // Intersect the chord from StartPointA to EndPointB
199     // and return whether an intersection occurred
200
201   G4FieldTrack ReEstimateEndpoint( const G4FieldTrack &CurrentStateA, 
202                                    const G4FieldTrack &EstimtdEndStateB,
203                                          G4double      linearDistSq,
204                                          G4double      curveDist);
205     // Return new estimate for state after curveDist
206     // starting from CurrentStateA,  to replace EstimtdEndStateB,
207     // (and report displacement -- if field is compiled verbose.)
208
209   void PrintStepLengthDiagnostic( G4double      currentProposedStepLength,
210                                   G4double      decreaseFactor,
211                                   G4double      stepTrial,
212                             const G4FieldTrack& aFieldTrack);
213 private:
214
215  // ----------------------------------------------------------------------
216  //  DATA Members
217  // ----------------------------------------------------------------------
218
219   G4FieldManager *fDetectorFieldMgr; 
220     // The  Field Manager of the whole Detector.  (default)
221
222   G4FieldManager *fCurrentFieldMgr;
223     // The  Field Manager of the current volume (may be the one above.)
224
225   G4Navigator   *fNavigator;
226
227   //  STATE information
228   //  -----------------
229
230   G4double    fEpsilonStep;
231     // Relative accuracy for current Step (Calc.)
232
233   G4FieldTrack    End_PointAndTangent;
234     // End point storage
235
236   G4bool      fParticleIsLooping;
237
238   G4int  fVerboseLevel;
239     // For debuging purposes
240
241   G4int  fMax_loop_count;
242     // Limit for the number of sub-steps taken in one call to ComputeStep
243
244   //  Variables to keep track of "abnormal" case - which causes loop
245   //
246   G4int     fNoZeroStep;                        //  Counter of zeroStep
247   G4int     fActionThreshold_NoZeroSteps;       //  Threshold: above this - act
248   G4int     fSevereActionThreshold_NoZeroSteps; //  Threshold to act harshly
249   G4int     fAbandonThreshold_NoZeroSteps;      //  Threshold to abandon
250
251   G4double  fFull_CurveLen_of_LastAttempt; 
252   G4double  fLast_ProposedStepLength; 
253   G4double  fLargestAcceptableStep;
254
255   G4double  fCharge, fInitialMomentumModulus, fMass;
256
257   G4ThreeVector  fPreviousSftOrigin;
258   G4double       fPreviousSafety; 
259   G4bool         fUseSafetyForOptimisation;
260     // Last safety origin & value: for optimisation
261
262   G4bool fSetFieldMgr; 
263     // Flag whether field manager has been set for the current step
264
265   G4double  kCarTolerance;
266     // Geometrical tolerance defining surface thickness
267
268private:
269
270   static const G4int max_depth=4;
271   G4FieldTrack* ptrInterMedFT[max_depth+1];
272     // Used to store intermediate values of tracks in case of
273     // too slow progress
274
275private:
276
277  G4VCurvedTrajectoryFilter* fpTrajectoryFilter;
278    // The filter encapsulates the algorithm which selects which
279    // intermediate points should be stored in a trajectory.
280    // When it is NULL, no intermediate points will be stored.
281    // Else PIF::ComputeStep must submit (all) intermediate
282    // points it calculates, to this filter.  (jacek 04/11/2002)
283};
284
285// ********************************************************************
286// Inline methods.
287// ********************************************************************
288
289#include "G4PropagatorInField.icc"
290
291#endif
Note: See TracBrowser for help on using the repository browser.