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

Last change on this file since 1202 was 921, checked in by garnier, 15 years ago

en test de gl2ps. Problemes de libraries

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