Geant4  10.02.p01
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: G4VIntersectionLocator.hh 93289 2015-10-15 10:01:15Z gcosmo $
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 stepNum);
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 
113  // Additional inline Set/Get methods for parameters, dependent objects
114 
116  inline G4double GetEpsilonStepFor();
117  inline G4Navigator* GetNavigatorFor();
119 
120  inline void SetSafetyParametersFor(G4bool UseSafety );
121 
122  inline void AddAdjustementOfFoundIntersection(G4bool UseCorrection);
124  // Methods to be made Obsolete - replaced by methods below
125  inline void AdjustIntersections(G4bool UseCorrection);
127  // Change adjustment flag ( New Interface )
128 
129  static void printStatus( const G4FieldTrack& startFT,
130  const G4FieldTrack& currentFT,
131  G4double requestStep,
132  G4double safety,
133  G4int stepNum,
134  std::ostream& oss,
135  G4int verboseLevel );
136  // Print Method for any ostream - e.g. cerr -- and for G4Exception
137 
138  protected: // with description
139 
140  G4FieldTrack ReEstimateEndpoint( const G4FieldTrack &CurrentStateA,
141  const G4FieldTrack &EstimtdEndStateB,
142  G4double linearDistSq, // not used
143  G4double curveDist ); // not used
144  // Return new estimate for state after curveDist starting from
145  // CurrentStateA, to replace EstimtdEndStateB, and report displacement
146  // (if field is compiled verbose)
147 
148  G4bool CheckAndReEstimateEndpoint( const G4FieldTrack& CurrentStartA,
149  const G4FieldTrack& EstimatedEndB,
150  G4FieldTrack& RevisedEndPoint,
151  G4int & errorCode);
152  // Check whether EndB is too far from StartA to be reached
153  // and if, re-estimate new value for EndB (return in RevisedEndPoint)
154  // Report error if EndB is before StartA (in curve length)
155  // In that case return errorCode = 2.
156 
157  G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point,
158  G4bool &validNormal); // const
159  // Position *must* be the intersection point from last call
160  // to G4Navigator's ComputeStep (via IntersectChord )
161  // Will try to use cached (last) value in Navigator for speed,
162  // if it was kept and valid.
163  // Value returned is in global coordinates.
164  // It does NOT guarantee to obtain Normal. This can happen eg if:
165  // - the "Intersection" Point is not on a surface, potentially due to
166  // - inaccuracies in the transformations used, or
167  // - issues with the Solid.
168 
169  G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector &CurrentE_Point,
170  G4bool &validNormal);
171  // Return the SurfaceNormal of Intersecting Solid in global coordinates
172  // Costlier then GetSurfaceNormal
173 
175  const G4ThreeVector &CurrentE_Point,
176  const G4ThreeVector &CurrentF_Point,
177  const G4ThreeVector &MomentumDir,
178  const G4bool IntersectAF,
179  G4ThreeVector &IntersectionPoint,
180  G4double &NewSafety,
181  G4double &fPrevSafety,
182  G4ThreeVector &fPrevSftOrigin );
183  // Optional method for adjustment of located intersection point
184  // using the surface-normal
185 
186  void ReportTrialStep( G4int step_no,
187  const G4ThreeVector& ChordAB_v,
188  const G4ThreeVector& ChordEF_v,
189  const G4ThreeVector& NewMomentumDir,
190  const G4ThreeVector& NormalAtEntry,
191  G4bool validNormal );
192  // Print a three-line report on the current "sub-step", ie trial
193  // intersection
194 
196  // Locate point using navigator - updates state of Navigator.
197  // By default, it assumes that the point is inside the current volume,
198  // and returns true.
199  // In check mode, checks whether the point is *inside* the volume.
200  // If it is inside, it returns true.
201  // If not, issues a warning and returns false.
202 
204  const G4String& CodeLocationInfo,
205  G4int CheckMode );
206  // As above, but report information about code location.
207  // If CheckMode > 1, report extra information.
208 
209  inline void SetCheckMode( G4bool value ) { fCheckMode = value; }
210  inline G4bool GetCheckMode() { return fCheckMode; }
211 
212  protected: // without description
213 
214  // Auxiliary methods -- to report issues
215 
216  void ReportReversedPoints( std::ostringstream& ossMsg,
217  const G4FieldTrack& StartPointVel,
218  const G4FieldTrack& EndPointVel,
219  G4double NewSafety, G4double epsStep,
220  const G4FieldTrack& CurrentA_PointVelocity,
221  const G4FieldTrack& CurrentB_PointVelocity,
222  const G4FieldTrack& SubStart_PointVelocity,
223  const G4ThreeVector& CurrentE_Point,
224  const G4FieldTrack& ApproxIntersecPointV,
225  G4int sbstp_no, G4int sbstp_no_p, G4int depth );
226  // Build error messsage (in ossMsg) to report that point 'B' has
227  // gone past 'A'
228 
229  void ReportProgress( std::ostream& oss,
230  const G4FieldTrack& StartPointVel,
231  const G4FieldTrack& EndPointVel,
232  G4int substep_no,
233  const G4FieldTrack& A_PtVel, // G4double safetyA
234  const G4FieldTrack& B_PtVel,
235  G4double safetyLast,
236  G4int depth= -1 );
237  // Report the current status / progress in finding the first intersection
238 
239  void ReportImmediateHit( const char* MethodName,
240  const G4ThreeVector& StartPosition,
241  const G4ThreeVector& TrialPoint,
242  double tolerance,
243  unsigned long int numCalls );
244  // Report case: trial point is 'close' to start, within tolerance
245 
246  private: // no description
247 
248  G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector &CurrentE_Point,
249  G4bool &validNormal);
250  // Return the SurfaceNormal of Intersecting Solid in local coordinates
251 
252  G4ThreeVector GetLastSurfaceNormal( const G4ThreeVector& intersectPoint,
253  G4bool &validNormal) const;
254  // Position *must* be the intersection point from last call
255  // to G4Navigator's ComputeStep (via IntersectChord )
256  // Temporary - will use the same method in the Navigator
257 
258  protected:
259 
260  G4double kCarTolerance; // Constant
261 
262  G4int fVerboseLevel; // For debugging
263  G4bool fUseNormalCorrection; // Configuration parameter
265 
267 
272  // Parameters set at each physical step by calling method
273  // by G4PropagatorInField
274 
276  // Helper for location
277 
279  // Touchable history hook
280 };
281 
283 
284 #endif
G4FieldTrack ReEstimateEndpoint(const G4FieldTrack &CurrentStateA, const G4FieldTrack &EstimtdEndStateB, G4double linearDistSq, G4double curveDist)
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
static const G4float tolerance
G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point, G4bool &validNormal)
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &PreviousSafety, G4ThreeVector &PreviousSftOrigin, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint, G4bool *calledNavigator=0)
G4ThreeVector GetLastSurfaceNormal(const G4ThreeVector &intersectPoint, G4bool &validNormal) const
#define errorCode
Definition: xmlparse.cc:618
void SetChordFinderFor(G4ChordFinder *fCFinder)
int G4int
Definition: G4Types.hh:78
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int stepNum)
G4double GetEpsilonStepFor()
void ReportImmediateHit(const char *MethodName, const G4ThreeVector &StartPosition, const G4ThreeVector &TrialPoint, double tolerance, unsigned long int numCalls)
double A(double temperature)
bool G4bool
Definition: G4Types.hh:79
void AdjustIntersections(G4bool UseCorrection)
void ReportReversedPoints(std::ostringstream &ossMsg, const G4FieldTrack &StartPointVel, const G4FieldTrack &EndPointVel, G4double NewSafety, G4double epsStep, const G4FieldTrack &CurrentA_PointVelocity, const G4FieldTrack &CurrentB_PointVelocity, const G4FieldTrack &SubStart_PointVelocity, const G4ThreeVector &CurrentE_Point, const G4FieldTrack &ApproxIntersecPointV, G4int sbstp_no, G4int sbstp_no_p, G4int depth)
void ReportTrialStep(G4int step_no, const G4ThreeVector &ChordAB_v, const G4ThreeVector &ChordEF_v, const G4ThreeVector &NewMomentumDir, const G4ThreeVector &NormalAtEntry, G4bool validNormal)
void LocateGlobalPointWithinVolumeCheckAndReport(const G4ThreeVector &pos, const G4String &CodeLocationInfo, G4int CheckMode)
G4Navigator * GetNavigatorFor()
void SetVerboseFor(G4int fVerbose)
G4bool LocateGlobalPointWithinVolumeAndCheck(const G4ThreeVector &pos)
#define fPreviousSftOrigin
G4double GetDeltaIntersectionFor()
G4VIntersectionLocator(G4Navigator *theNavigator)
void AddAdjustementOfFoundIntersection(G4bool UseCorrection)
#define fPreviousSafety
G4bool CheckAndReEstimateEndpoint(const G4FieldTrack &CurrentStartA, const G4FieldTrack &EstimatedEndB, G4FieldTrack &RevisedEndPoint, G4int &errorCode)
G4bool GetAdjustementOfFoundIntersection()
G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
void SetNavigatorFor(G4Navigator *fNavigator)
void SetSafetyParametersFor(G4bool UseSafety)
virtual G4bool EstimateIntersectionPoint(const G4FieldTrack &curveStartPointTangent, const G4FieldTrack &curveEndPointTangent, const G4ThreeVector &trialPoint, G4FieldTrack &intersectPointTangent, G4bool &recalculatedEndPoint, G4double &fPreviousSafety, G4ThreeVector &fPreviousSftOrigin)=0
double G4double
Definition: G4Types.hh:76
G4TouchableHistory * fpTouchable
void ReportProgress(std::ostream &oss, const G4FieldTrack &StartPointVel, const G4FieldTrack &EndPointVel, G4int substep_no, const G4FieldTrack &A_PtVel, const G4FieldTrack &B_PtVel, G4double safetyLast, G4int depth=-1)
void SetEpsilonStepFor(G4double EpsilonStep)
G4ChordFinder * GetChordFinderFor()
static const G4double pos
G4bool AdjustmentOfFoundIntersection(const G4ThreeVector &A, const G4ThreeVector &CurrentE_Point, const G4ThreeVector &CurrentF_Point, const G4ThreeVector &MomentumDir, const G4bool IntersectAF, G4ThreeVector &IntersectionPoint, G4double &NewSafety, G4double &fPrevSafety, G4ThreeVector &fPrevSftOrigin)
void SetDeltaIntersectionFor(G4double deltaIntersection)