Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PropagatorInField.hh
이 파일의 문서화 페이지로 가기
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 110728 2018-06-11 06:12:39Z 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 void SetVerboseTrace( G4bool enable ) { fVerbTracePiF = enable; }
103  inline G4bool GetVerboseTrace() { return fVerbTracePiF; }
104  // Tracing key parts of Compute Step
105 
106  inline G4int GetMaxLoopCount() const;
107  inline void SetMaxLoopCount( G4int new_max );
108  // A maximum for the number of steps that a (looping) particle can take.
109 
110  void printStatus( const G4FieldTrack& startFT,
111  const G4FieldTrack& currentFT,
112  G4double requestStep,
113  G4double safety,
114  G4int step,
115  G4VPhysicalVolume* startVolume);
116  // Print Method - useful mostly for debugging.
117 
118  inline G4FieldTrack GetEndState() const;
119 
120  inline G4double GetMinimumEpsilonStep() const; // Min for relative accuracy
121  inline void SetMinimumEpsilonStep( G4double newEpsMin ); // of any step
122  inline G4double GetMaximumEpsilonStep() const;
123  inline void SetMaximumEpsilonStep( G4double newEpsMax );
124  // The 4 above methods are now obsolescent but *for now* will work
125  // They are being replaced by same-name methods in G4FieldManager,
126  // allowing the specialisation in different volumes.
127  // Their new behaviour is to change the values for the global field
128  // manager
129 
130  inline void SetLargestAcceptableStep( G4double newBigDist );
132 
134  // Set the filter that examines & stores 'intermediate'
135  // curved trajectory points. Currently only position is stored.
136 
137  std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
138  // Access the points which have passed by the filter.
139  // Responsibility for deleting the points lies with the client.
140  // This method MUST BE called exactly ONCE per step.
141 
142  void ClearPropagatorState();
143  // Clear all the State of this class and its current associates
144  // --> the current field manager & chord finder will also be called
145 
146  inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
147  // Update this (dangerous) state -- for the time being
148 
149  inline void SetUseSafetyForOptimization( G4bool );
151  // Toggle & view parameter for using safety to discard
152  // unneccesary calls to navigator (thus 'optimising' performance)
153  inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
154  const G4ThreeVector& EndPointB,
155  G4double &NewSafety,
156  G4double &LinearStepLength,
157  G4ThreeVector &IntersectionPoint);
158  // Intersect the chord from StartPointA to EndPointB
159  // and return whether an intersection occurred
160  // NOTE : SAFETY IS CHANGED
161 
165 
167  inline void SetIntersectionLocator(G4VIntersectionLocator *pLocator );
168  // Change or get the object which calculates the exact
169  // intersection point with the next boundary
170 
171  public: // without description
172 
173  inline G4double GetDeltaIntersection() const;
174  inline G4double GetDeltaOneStep() const;
175 
178  // Auxiliary methods - their results can/will change during propagation
179 
180  inline void SetNavigatorForPropagating( G4Navigator *SimpleOrMultiNavigator );
182 
183  inline void SetThresholdNoZeroStep( G4int noAct,
184  G4int noHarsh,
185  G4int noAbandon );
186  inline G4int GetThresholdNoZeroSteps( G4int i );
187 
188  inline G4double GetZeroStepThreshold();
189  inline void SetZeroStepThreshold( G4double newLength );
190 
192  // Update the Locator with parameters from this class
193  // and from current field manager
194 
195  protected: // without description
196 
197  void PrintStepLengthDiagnostic( G4double currentProposedStepLength,
198  G4double decreaseFactor,
199  G4double stepTrial,
200  const G4FieldTrack& aFieldTrack);
201 
202  void ReportLoopingParticle( G4int count, G4double StepTaken,
203  G4double stepRequest, const char* methodName,
204  G4ThreeVector momentumVec,
205  G4VPhysicalVolume* physVol);
206  void ReportStuckParticle( G4int noZeroSteps, G4double proposedStep,
207  G4double lastTriedStep, G4VPhysicalVolume* physVol);
208 
209  private:
210  // ----------------------------------------------------------------------
211  // DATA Members
212  // ----------------------------------------------------------------------
213 
214  // ==================================================================
215  // INVARIANTS - Must not change during tracking
216 
217  // ** PARAMETERS -----------
219  // Limit for the number of sub-steps taken in one call to ComputeStep
221 
222  // Thresholds for identifying "abnormal" cases - which cause looping
223  G4int fActionThreshold_NoZeroSteps; // Threshold # - above it act
224  G4int fSevereActionThreshold_NoZeroSteps; // Threshold # to act harshly
225  G4int fAbandonThreshold_NoZeroSteps; // Threshold # to abandon
227  // Threshold *length* for counting of tiny or 'zero' steps
228 
230  // Maximum size of a step - for optimization (and to avoid problems)
231  // ** End of PARAMETERS -----
232 
234  // Geometrical tolerance defining surface thickness
235 
236  G4bool fAllocatedLocator; // Book-keeping
237 
238  // --------------------------------------------------------
239  // ** Dependent Objects - to which work is delegated
240 
242  // The Field Manager of the whole Detector. (default)
243 
245  // Refines candidate intersection
246 
248  // The filter encapsulates the algorithm which selects which
249  // intermediate points should be stored in a trajectory.
250  // When it is NULL, no intermediate points will be stored.
251  // Else PIF::ComputeStep must submit (all) intermediate
252  // points it calculates, to this filter. (jacek 04/11/2002)
253 
255  // Set externally - only by tracking / run manager
256  //
257  // ** End of Dependent Objects ----------------------------
258 
259  // End of INVARIANTS
260  // ==================================================================
261 
262  // STATE information
263  // -----------------
265  // The Field Manager of the current volume (may be the global)
266  G4bool fSetFieldMgr; // Has it been set for the current step
267 
268  // Parameters of current step
269  G4double fEpsilonStep; // Relative accuracy of current Step
270  G4FieldTrack End_PointAndTangent; // End point storage
272  G4int fNoZeroStep; // Count of zero Steps
273 
274  // State used for Optimisation
277  // Previous step information -- for use in adjust step size
280  // Last safety origin & value: for optimisation
281 
284  // For debugging purposes
285 
289 };
290 
291 // Inline methods.
292 // *******************************
293 
294 #include "G4PropagatorInField.icc"
295 
296 #endif
G4ThreeVector fPreviousSftOrigin
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
void SetMaximumEpsilonStep(G4double newEpsMax)
G4int GetThresholdNoZeroSteps(G4int i)
G4VIntersectionLocator * fIntersectionLocator
void SetZeroStepThreshold(G4double newLength)
void ReportLoopingParticle(G4int count, G4double StepTaken, G4double stepRequest, const char *methodName, G4ThreeVector momentumVec, G4VPhysicalVolume *physVol)
G4int Verbose() const
G4EquationOfMotion * GetCurrentEquationOfMotion()
G4FieldManager * GetCurrentFieldManager()
void SetLargestAcceptableStep(G4double newBigDist)
G4FieldManager * fDetectorFieldMgr
G4ThreeVector EndPosition() const
void SetUseSafetyForOptimization(G4bool)
void SetTrajectoryFilter(G4VCurvedTrajectoryFilter *filter)
void SetVerboseTrace(G4bool enable)
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
G4FieldTrack End_PointAndTangent
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
G4VCurvedTrajectoryFilter * fpTrajectoryFilter
G4FieldTrack GetEndState() const
void SetIntersectionLocator(G4VIntersectionLocator *pLocator)
void SetNavigatorForPropagating(G4Navigator *SimpleOrMultiNavigator)
G4double GetEpsilonStep() const
G4VIntersectionLocator * GetIntersectionLocator()
G4ThreeVector EndMomentumDir() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4double GetDeltaOneStep() const
void SetMaxLoopCount(G4int new_max)
void SetThresholdNoZeroStep(G4int noAct, G4int noHarsh, G4int noAbandon)
void SetDetectorFieldManager(G4FieldManager *newGlobalFieldManager)
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4FieldManager * fCurrentFieldMgr
G4PropagatorInField(G4Navigator *theNavigator, G4FieldManager *detectorFieldMgr, G4VIntersectionLocator *vLocator=0)
G4double GetLargestAcceptableStep()
void SetMinimumEpsilonStep(G4double newEpsMin)
G4Navigator * GetNavigatorForPropagating()
int G4int
Definition: G4Types.hh:78
G4double GetZeroStepThreshold()
void SetEpsilonStep(G4double newEps)
G4bool IsParticleLooping() const
G4int GetMaxLoopCount() const
void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)
G4double GetDeltaIntersection() const
G4ChordFinder * GetChordFinder()
G4bool GetUseSafetyForOptimization()
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
G4int GetVerboseLevel() const
G4int SetVerboseLevel(G4int verbose)
G4double GetMinimumEpsilonStep() const
G4double GetMaximumEpsilonStep() const