Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
58class G4ChordFinder;
59
60class G4Navigator;
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;
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
90 inline void SetChargeMomentumMass( G4double charge, // in e+ units
91 G4double momentum, // in Geant4 units
92 G4double pMass );
93 // Inform this and all associated classes of q, p, m
94
96 // Set (and return) the correct field manager (global or local),
97 // if it exists.
98 // Should be called before ComputeStep is called;
99 // - currently, ComputeStep will call it, if it has not been called.
100
102
103 G4int SetVerboseLevel( G4int verbose );
104 inline G4int GetVerboseLevel() const;
105 inline G4int Verbose() const;
106
107 inline G4int GetMaxLoopCount() const;
108 inline void SetMaxLoopCount( G4int new_max );
109 // A maximum for the number of steps that a (looping) particle can take.
110
111 void printStatus( const G4FieldTrack& startFT,
112 const G4FieldTrack& currentFT,
113 G4double requestStep,
114 G4double safety,
115 G4int step,
116 G4VPhysicalVolume* startVolume);
117 // Print Method - useful mostly for debugging.
118
119 inline G4FieldTrack GetEndState() const;
120
121 inline G4double GetMinimumEpsilonStep() const; // Min for relative accuracy
122 inline void SetMinimumEpsilonStep( G4double newEpsMin ); // of any step
124 inline void SetMaximumEpsilonStep( G4double newEpsMax );
125 // The 4 above methods are now obsolescent but *for now* will work
126 // They are being replaced by same-name methods in G4FieldManager,
127 // allowing the specialisation in different volumes.
128 // Their new behaviour is to change the values for the global field
129 // manager
130
131 inline void SetLargestAcceptableStep( G4double newBigDist );
133
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
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
152 // Toggle & view parameter for using safety to discard
153 // unneccesary calls to navigator (thus 'optimising' performance)
154 inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
155 const 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
165 // Change or get the object which calculates the exact
166 // intersection point with the next boundary
167
168 public: // without description
169
171 inline G4double GetDeltaOneStep() const;
172
174 inline void SetNavigatorForPropagating( G4Navigator *SimpleOrMultiNavigator );
176
177 inline void SetThresholdNoZeroStep( G4int noAct,
178 G4int noHarsh,
179 G4int noAbandon );
181
183 inline void SetZeroStepThreshold( G4double newLength );
184
186 // Update the Locator with parameters from this class
187 // and from current field manager
188
189 protected: // with description
190
191 void PrintStepLengthDiagnostic( G4double currentProposedStepLength,
192 G4double decreaseFactor,
193 G4double stepTrial,
194 const G4FieldTrack& aFieldTrack);
195 private:
196 // ----------------------------------------------------------------------
197 // DATA Members
198 // ----------------------------------------------------------------------
199
200 // ==================================================================
201 // INVARIANTS - Must not change during tracking
202
203 // ** PARAMETERS -----------
204 G4int fMax_loop_count;
205 // Limit for the number of sub-steps taken in one call to ComputeStep
206 G4bool fUseSafetyForOptimisation;
207
208 // Thresholds for identifying "abnormal" cases - which cause looping
209 G4int fActionThreshold_NoZeroSteps; // Threshold # - above it act
210 G4int fSevereActionThreshold_NoZeroSteps; // Threshold # to act harshly
211 G4int fAbandonThreshold_NoZeroSteps; // Threshold # to abandon
212 G4double fZeroStepThreshold;
213 // Threshold *length* for counting of tiny or 'zero' steps
214
215 G4double fLargestAcceptableStep;
216 // Maximum size of a step - for optimization (and to avoid problems)
217 // ** End of PARAMETERS -----
218
219 G4double kCarTolerance;
220 // Geometrical tolerance defining surface thickness
221
222 G4bool fAllocatedLocator; // Book-keeping
223
224 // --------------------------------------------------------
225 // ** Dependent Objects - to which work is delegated
226
227 G4FieldManager *fDetectorFieldMgr;
228 // The Field Manager of the whole Detector. (default)
229
230 G4VIntersectionLocator *fIntersectionLocator;
231 // Refines candidate intersection
232
233 G4VCurvedTrajectoryFilter* fpTrajectoryFilter;
234 // The filter encapsulates the algorithm which selects which
235 // intermediate points should be stored in a trajectory.
236 // When it is NULL, no intermediate points will be stored.
237 // Else PIF::ComputeStep must submit (all) intermediate
238 // points it calculates, to this filter. (jacek 04/11/2002)
239
240 G4Navigator *fNavigator;
241 // Set externally - only by tracking / run manager
242 //
243 // ** End of Dependent Objects ----------------------------
244
245 // End of INVARIANTS
246 // ==================================================================
247
248 // STATE information
249 // -----------------
250 G4FieldManager *fCurrentFieldMgr;
251 // The Field Manager of the current volume (may be the global)
252 G4bool fSetFieldMgr; // Has it been set for the current step
253
254 // Parameters of current step
255 G4double fCharge, fInitialMomentumModulus, fMass;
256 G4double fEpsilonStep; // Relative accuracy of current Step
257 G4FieldTrack End_PointAndTangent; // End point storage
258 G4bool fParticleIsLooping;
259 G4int fNoZeroStep; // Count of zero Steps
260
261 // State used for Optimisation
262 G4double fFull_CurveLen_of_LastAttempt;
263 G4double fLast_ProposedStepLength;
264 // Previous step information -- for use in adjust step size
265 G4ThreeVector fPreviousSftOrigin;
266 G4double fPreviousSafety;
267 // Last safety origin & value: for optimisation
268
269 G4int fVerboseLevel;
270 // For debuging purposes
271
272 private:
273};
274
275// Inline methods.
276// *******************************
277
278#include "G4PropagatorInField.icc"
279
280#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4ThreeVector EndPosition() const
void SetThresholdNoZeroStep(G4int noAct, G4int noHarsh, G4int noAbandon)
G4int GetThresholdNoZeroSteps(G4int i)
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4double GetDeltaOneStep() const
G4FieldManager * GetCurrentFieldManager()
G4int SetVerboseLevel(G4int verbose)
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
G4double GetZeroStepThreshold()
void SetMaxLoopCount(G4int new_max)
void SetNavigatorForPropagating(G4Navigator *SimpleOrMultiNavigator)
G4int Verbose() const
G4bool GetUseSafetyForOptimization()
G4bool IsParticleLooping() const
void SetZeroStepThreshold(G4double newLength)
void SetDetectorFieldManager(G4FieldManager *newGlobalFieldManager)
G4double GetMaximumEpsilonStep() const
G4ChordFinder * GetChordFinder()
void SetUseSafetyForOptimization(G4bool)
G4double GetDeltaIntersection() const
void SetChargeMomentumMass(G4double charge, G4double momentum, G4double pMass)
G4double GetEpsilonStep() const
void SetMaximumEpsilonStep(G4double newEpsMax)
G4Navigator * GetNavigatorForPropagating()
G4FieldTrack GetEndState() const
G4ThreeVector EndMomentumDir() const
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
G4VIntersectionLocator * GetIntersectionLocator()
G4int GetVerboseLevel() const
void SetMinimumEpsilonStep(G4double newEpsMin)
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4int GetMaxLoopCount() const
void SetIntersectionLocator(G4VIntersectionLocator *pLocator)
G4double GetLargestAcceptableStep()
void SetTrajectoryFilter(G4VCurvedTrajectoryFilter *filter)
G4double GetMinimumEpsilonStep() const
void SetLargestAcceptableStep(G4double newBigDist)
void SetEpsilonStep(G4double newEps)