Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VIntersectionLocator.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 G4VIntersectionLocator
31 //
32 // class description:
33 //
34 // Base class for the calculation of the intersection point with a boundary
35 // when PropagationInField is used.
36 // Gives possibility to choose the method of intersection; concrete locators
37 // implemented are: G4SimpleLocator, G4MultiLevelLocator, G4BrentLocator.
38 //
39 // Key Method: EstimateIntersectionPoint()
40 
41 // History:
42 // -------
43 // 27.10.08 - John Apostolakis, Tatiana Nikitina: Design and implementation
44 // ---------------------------------------------------------------------------
45 
46 #ifndef G4VINTERSECTIONLOCATOR_HH
47 #define G4VINTERSECTIONLOCATOR_HH
48 
49 #include "G4Types.hh"
50 #include "G4ThreeVector.hh"
51 #include "G4FieldTrack.hh"
52 
53 #include "G4Navigator.hh"
54 #include "G4ChordFinder.hh"
55 
57  {
58  public: // with description
59 
60  G4VIntersectionLocator(G4Navigator *theNavigator);
61  // Constructor
62  virtual ~G4VIntersectionLocator();
63  // Default destructor
64 
66  const G4FieldTrack& curveStartPointTangent, // A
67  const G4FieldTrack& curveEndPointTangent, // B
68  const G4ThreeVector& trialPoint, // E
69  G4FieldTrack& intersectPointTangent, // Output
70  G4bool& recalculatedEndPoint, // Out
71  G4double& fPreviousSafety, // In/Out
72  G4ThreeVector& fPreviousSftOrigin) = 0; // In/Out
73  // If such an intersection exists, this function calculates the
74  // intersection point of the true path of the particle with the surface
75  // of the current volume (or of one of its daughters).
76  // Should use lateral displacement as measure of convergence
77  // NOTE: changes the safety!
78 
79  void printStatus( const G4FieldTrack& startFT,
80  const G4FieldTrack& currentFT,
81  G4double requestStep,
82  G4double safety,
83  G4int step);
84  // Print Method, useful mostly for debugging
85 
86  inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
87  const G4ThreeVector& EndPointB,
88  G4double &NewSafety,
89  G4double &PreviousSafety, // In/Out
90  G4ThreeVector &PreviousSftOrigin, // In/Out
91  G4double &LinearStepLength,
92  G4ThreeVector &IntersectionPoint,
93  G4bool *calledNavigator=0 );
94  // Intersect the chord from StartPointA to EndPointB and return
95  // whether an intersection occurred. NOTE: changes the Safety!
96 
97  inline void SetEpsilonStepFor( G4double EpsilonStep );
98  inline void SetDeltaIntersectionFor( G4double deltaIntersection );
99  inline void SetNavigatorFor( G4Navigator *fNavigator );
100  inline void SetChordFinderFor(G4ChordFinder *fCFinder );
101  // These parameters must be set at each step, in case they were changed
102 
103  // Note: This simple approach ensures that all scenarios are considered.
104  // [ Future refinement may identify which are invariant during a
105  // track, run or event ]
106 
107  inline void SetVerboseFor(G4int fVerbose);
108  inline G4int GetVerboseFor();
109  // Controling verbosity enables checking of the locating of intersections
110 
111  public: // without description
112  // Additional inline Set/Get methods for parameters, dependent objects
113 
115  inline G4double GetEpsilonStepFor();
116  inline G4Navigator* GetNavigatorFor();
118 
119  inline void SetSafetyParametersFor(G4bool UseSafety );
120 
121  inline void AddAdjustementOfFoundIntersection(G4bool UseCorrection);
123  // Methods to be made Obsolete - replaced by methods below
124  inline void AdjustIntersections(G4bool UseCorrection);
126  // Change adjustment flag ( New Interface )
127 
128  protected: // with description
129 
130  G4FieldTrack ReEstimateEndpoint( const G4FieldTrack &CurrentStateA,
131  const G4FieldTrack &EstimtdEndStateB,
132  G4double linearDistSq,
133  G4double curveDist);
134  // Return new estimate for state after curveDist starting from
135  // CurrentStateA, to replace EstimtdEndStateB, and report displacement
136  // (if field is compiled verbose)
137 
138  G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point,
139  G4bool &validNormal); // const
140  // Position *must* be the intersection point from last call
141  // to G4Navigator's ComputeStep (via IntersectChord )
142  // Will try to use cached (last) value in Navigator for speed,
143  // if it was kept and valid.
144  // Value returned is in global coordinates.
145  // It does NOT guarantee to obtain Normal. This can happen eg if:
146  // - the "Intersection" Point is not on a surface, potentially due to
147  // - inaccuracies in the transformations used, or
148  // - issues with the Solid.
149 
150  G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector &CurrentE_Point,
151  G4bool &validNormal);
152  // Return the SurfaceNormal of Intersecting Solid in global coordinates
153  // Costlier then GetSurfaceNormal
154 
156  const G4ThreeVector &CurrentE_Point,
157  const G4ThreeVector &CurrentF_Point,
158  const G4ThreeVector &MomentumDir,
159  const G4bool IntersectAF,
160  G4ThreeVector &IntersectionPoint,
161  G4double &NewSafety,
162  G4double &fPrevSafety,
163  G4ThreeVector &fPrevSftOrigin );
164  // Optional method for adjustment of located intersection point
165  // using the surface-normal
166 
167  void ReportTrialStep( G4int step_no,
168  const G4ThreeVector& ChordAB_v,
169  const G4ThreeVector& ChordEF_v,
170  const G4ThreeVector& NewMomentumDir,
171  const G4ThreeVector& NormalAtEntry,
172  G4bool validNormal );
173  // Print a three-line report on the current "sub-step", ie trial intersection
174 
175  private: // no description
176 
177  G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector &CurrentE_Point,
178  G4bool &validNormal);
179  // Return the SurfaceNormal of Intersecting Solid in local coordinates
180 
181  G4ThreeVector GetLastSurfaceNormal( const G4ThreeVector& intersectPoint,
182  G4bool &validNormal) const;
183  // Position *must* be the intersection point from last call
184  // to G4Navigator's ComputeStep (via IntersectChord )
185  // Temporary - will use the same method in the Navigator
186 
187  protected:
188 
189  G4double kCarTolerance; // Constant
190 
191  G4int fVerboseLevel; // For debugging
192  G4bool fUseNormalCorrection; // Configuration parameter
193 
195 
200  // Parameters set at each physical step by calling method
201  // by G4PropagatorInField
202 
204  // Helper for location
205 
207  // Touchable history hook
208 };
209 
210 #include "G4VIntersectionLocator.icc"
211 
212 #endif