Geant4_10
G4PropagatorInField.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 //
27 // $Id: G4PropagatorInField.hh 69711 2013-05-13 09:22:17Z gcosmo $
28 //
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: ComputeStep(..)
39 
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"
57 
58 class G4ChordFinder;
59 
60 class G4Navigator;
61 class G4VPhysicalVolume;
63 
65 {
66 
67  public: // with description
68 
69  G4PropagatorInField( G4Navigator *theNavigator,
70  G4FieldManager *detectorFieldMgr,
71  G4VIntersectionLocator *vLocator=0 );
73 
74  G4double ComputeStep( G4FieldTrack &pFieldTrack,
75  G4double pCurrentProposedStepLength,
76  G4double &pNewSafety,
77  G4VPhysicalVolume *pPhysVol=0 );
78  // Compute the next geometric Step
79 
80  inline G4ThreeVector EndPosition() const;
81  inline G4ThreeVector EndMomentumDir() const;
82  inline G4bool IsParticleLooping() const;
83  // Return the state after the Step
84 
85  inline G4double GetEpsilonStep() const;
86  // Relative accuracy for current Step (Calc.)
87  inline void SetEpsilonStep(G4double newEps);
88  // The ratio DeltaOneStep()/h_current_step
89 
91  // Set (and return) the correct field manager (global or local),
92  // if it exists.
93  // Should be called before ComputeStep is called;
94  // - currently, ComputeStep will call it, if it has not been called.
95 
96  inline G4ChordFinder* GetChordFinder();
97 
98  G4int SetVerboseLevel( G4int verbose );
99  inline G4int GetVerboseLevel() const;
100  inline G4int Verbose() const;
101 
102  inline G4int GetMaxLoopCount() const;
103  inline void SetMaxLoopCount( G4int new_max );
104  // A maximum for the number of steps that a (looping) particle can take.
105 
106  void printStatus( const G4FieldTrack& startFT,
107  const G4FieldTrack& currentFT,
108  G4double requestStep,
109  G4double safety,
110  G4int step,
111  G4VPhysicalVolume* startVolume);
112  // Print Method - useful mostly for debugging.
113 
114  inline G4FieldTrack GetEndState() const;
115 
116  inline G4double GetMinimumEpsilonStep() const; // Min for relative accuracy
117  inline void SetMinimumEpsilonStep( G4double newEpsMin ); // of any step
118  inline G4double GetMaximumEpsilonStep() const;
119  inline void SetMaximumEpsilonStep( G4double newEpsMax );
120  // The 4 above methods are now obsolescent but *for now* will work
121  // They are being replaced by same-name methods in G4FieldManager,
122  // allowing the specialisation in different volumes.
123  // Their new behaviour is to change the values for the global field
124  // manager
125 
126  inline void SetLargestAcceptableStep( G4double newBigDist );
128 
130  // Set the filter that examines & stores 'intermediate'
131  // curved trajectory points. Currently only position is stored.
132 
133  std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
134  // Access the points which have passed by the filter.
135  // Responsibility for deleting the points lies with the client.
136  // This method MUST BE called exactly ONCE per step.
137 
138  void ClearPropagatorState();
139  // Clear all the State of this class and its current associates
140  // --> the current field manager & chord finder will also be called
141 
142  inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
143  // Update this (dangerous) state -- for the time being
144 
145  inline void SetUseSafetyForOptimization( G4bool );
147  // Toggle & view parameter for using safety to discard
148  // unneccesary calls to navigator (thus 'optimising' performance)
149  inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
150  const G4ThreeVector& EndPointB,
151  G4double &NewSafety,
152  G4double &LinearStepLength,
153  G4ThreeVector &IntersectionPoint);
154  // Intersect the chord from StartPointA to EndPointB
155  // and return whether an intersection occurred
156  // NOTE : SAFETY IS CHANGED
157 
159  inline void SetIntersectionLocator(G4VIntersectionLocator *pLocator );
160  // Change or get the object which calculates the exact
161  // intersection point with the next boundary
162 
163  public: // without description
164 
165  inline G4double GetDeltaIntersection() const;
166  inline G4double GetDeltaOneStep() const;
167 
169  inline void SetNavigatorForPropagating( G4Navigator *SimpleOrMultiNavigator );
171 
172  inline void SetThresholdNoZeroStep( G4int noAct,
173  G4int noHarsh,
174  G4int noAbandon );
175  inline G4int GetThresholdNoZeroSteps( G4int i );
176 
177  inline G4double GetZeroStepThreshold();
178  inline void SetZeroStepThreshold( G4double newLength );
179 
181  // Update the Locator with parameters from this class
182  // and from current field manager
183 
184  protected: // with description
185 
186  void PrintStepLengthDiagnostic( G4double currentProposedStepLength,
187  G4double decreaseFactor,
188  G4double stepTrial,
189  const G4FieldTrack& aFieldTrack);
190  private:
191  // ----------------------------------------------------------------------
192  // DATA Members
193  // ----------------------------------------------------------------------
194 
195  // ==================================================================
196  // INVARIANTS - Must not change during tracking
197 
198  // ** PARAMETERS -----------
199  G4int fMax_loop_count;
200  // Limit for the number of sub-steps taken in one call to ComputeStep
201  G4bool fUseSafetyForOptimisation;
202 
203  // Thresholds for identifying "abnormal" cases - which cause looping
204  G4int fActionThreshold_NoZeroSteps; // Threshold # - above it act
205  G4int fSevereActionThreshold_NoZeroSteps; // Threshold # to act harshly
206  G4int fAbandonThreshold_NoZeroSteps; // Threshold # to abandon
207  G4double fZeroStepThreshold;
208  // Threshold *length* for counting of tiny or 'zero' steps
209 
210  G4double fLargestAcceptableStep;
211  // Maximum size of a step - for optimization (and to avoid problems)
212  // ** End of PARAMETERS -----
213 
214  G4double kCarTolerance;
215  // Geometrical tolerance defining surface thickness
216 
217  G4bool fAllocatedLocator; // Book-keeping
218 
219  // --------------------------------------------------------
220  // ** Dependent Objects - to which work is delegated
221 
222  G4FieldManager *fDetectorFieldMgr;
223  // The Field Manager of the whole Detector. (default)
224 
225  G4VIntersectionLocator *fIntersectionLocator;
226  // Refines candidate intersection
227 
228  G4VCurvedTrajectoryFilter* fpTrajectoryFilter;
229  // The filter encapsulates the algorithm which selects which
230  // intermediate points should be stored in a trajectory.
231  // When it is NULL, no intermediate points will be stored.
232  // Else PIF::ComputeStep must submit (all) intermediate
233  // points it calculates, to this filter. (jacek 04/11/2002)
234 
235  G4Navigator *fNavigator;
236  // Set externally - only by tracking / run manager
237  //
238  // ** End of Dependent Objects ----------------------------
239 
240  // End of INVARIANTS
241  // ==================================================================
242 
243  // STATE information
244  // -----------------
245  G4FieldManager *fCurrentFieldMgr;
246  // The Field Manager of the current volume (may be the global)
247  G4bool fSetFieldMgr; // Has it been set for the current step
248 
249  // Parameters of current step
250  G4double fEpsilonStep; // Relative accuracy of current Step
251  G4FieldTrack End_PointAndTangent; // End point storage
252  G4bool fParticleIsLooping;
253  G4int fNoZeroStep; // Count of zero Steps
254 
255  // State used for Optimisation
256  G4double fFull_CurveLen_of_LastAttempt;
257  G4double fLast_ProposedStepLength;
258  // Previous step information -- for use in adjust step size
259  G4ThreeVector fPreviousSftOrigin;
260  G4double fPreviousSafety;
261  // Last safety origin & value: for optimisation
262 
263  G4int fVerboseLevel;
264  // For debuging purposes
265 
266  private:
267 };
268 
269 // Inline methods.
270 // *******************************
271 
272 #include "G4PropagatorInField.icc"
273 
274 #endif
G4double GetZeroStepThreshold()
pid_t filter
Definition: tracer.cxx:30
void SetEpsilonStep(G4double newEps)
void SetNavigatorForPropagating(G4Navigator *SimpleOrMultiNavigator)
G4int GetMaxLoopCount() const
G4Navigator * GetNavigatorForPropagating()
void SetIntersectionLocator(G4VIntersectionLocator *pLocator)
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4int Verbose() const
void SetThresholdNoZeroStep(G4int noAct, G4int noHarsh, G4int noAbandon)
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
G4ThreeVector EndPosition() const
void SetUseSafetyForOptimization(G4bool)
G4double GetDeltaIntersection() const
G4ThreeVector EndMomentumDir() const
G4PropagatorInField(G4Navigator *theNavigator, G4FieldManager *detectorFieldMgr, G4VIntersectionLocator *vLocator=0)
int G4int
Definition: G4Types.hh:78
G4double GetMaximumEpsilonStep() const
G4int SetVerboseLevel(G4int verbose)
G4FieldTrack GetEndState() const
bool G4bool
Definition: G4Types.hh:79
G4int GetVerboseLevel() const
G4bool GetUseSafetyForOptimization()
void SetZeroStepThreshold(G4double newLength)
G4double GetDeltaOneStep() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
void SetMaxLoopCount(G4int new_max)
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
G4int GetThresholdNoZeroSteps(G4int i)
G4bool IsParticleLooping() const
G4ChordFinder * GetChordFinder()
G4FieldManager * GetCurrentFieldManager()
void SetMinimumEpsilonStep(G4double newEpsMin)
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
void SetTrajectoryFilter(G4VCurvedTrajectoryFilter *filter)
G4double GetMinimumEpsilonStep() const
void SetLargestAcceptableStep(G4double newBigDist)
double G4double
Definition: G4Types.hh:76
void SetDetectorFieldManager(G4FieldManager *newGlobalFieldManager)
void SetMaximumEpsilonStep(G4double newEpsMax)
G4double GetLargestAcceptableStep()
G4VIntersectionLocator * GetIntersectionLocator()
G4double GetEpsilonStep() const
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)