Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// G4VIntersectionLocator
27//
28// Class description:
29//
30// Base class for the calculation of the intersection point with a boundary
31// when PropagationInField is used.
32// Gives possibility to choose the method of intersection; concrete locators
33// implemented are: G4SimpleLocator, G4MultiLevelLocator, G4BrentLocator.
34// Key Method: EstimateIntersectionPoint().
35
36// Authors: John Apostolakis, Tatiana Nikitina (CERN), 27 October 2008
37// --------------------------------------------------------------------
38#ifndef G4VINTERSECTIONLOCATOR_HH
39#define G4VINTERSECTIONLOCATOR_HH 1
40
41#include "G4Types.hh"
42#include "G4ThreeVector.hh"
43#include "G4FieldTrack.hh"
44
45#include "G4Navigator.hh"
46#include "G4ChordFinder.hh"
47
48/**
49 * @brief G4VIntersectionLocator is a base class for the calculation of the
50 * intersection point with a boundary when PropagationInField is used.
51 * It gives the possibility to choose the method of intersection.
52 */
53
55{
56 public:
57
58 /**
59 * Constructor and virtual Destructor.
60 */
63
64 /**
65 * If such an intersection exists, this method calculates the intersection
66 * point of the true path of the particle with the surface of the current
67 * volume (or of one of its daughters).
68 * Should use lateral displacement as measure of convergence.
69 * @note Changes the safety!
70 * @param[in] curveStartPointTangent Start point tangent track.
71 * @param[in] curveEndPointTangent End point tangent track.
72 * @param[in] trialPoint Trial point.
73 * @param[out] intersectPointTangent Intersection point tangent track.
74 * @param[out] recalculatedEndPoint Flagging if end point was recomputed.
75 * @param[in,out] fPreviousSafety Previous safety distance.
76 * @param[in,out] fPreviousSftOrigin Previous safety point origin.
77 * @returns Whether intersection exists or not.
78 */
80 const G4FieldTrack& curveStartPointTangent, // A
81 const G4FieldTrack& curveEndPointTangent, // B
82 const G4ThreeVector& trialPoint, // E
83 G4FieldTrack& intersectPointTangent, // Output
84 G4bool& recalculatedEndPoint, // Out
85 G4double& fPreviousSafety, // In/Out
86 G4ThreeVector& fPreviousSftOrigin) = 0; // In/Out
87
88 /**
89 * Intersects the chord from StartPointA to EndPointB and returns
90 * whether an intersection occurred.
91 * @note Changes the safety!
92 * @param[in] StartPointA Chord starting point.
93 * @param[in] EndPointB Chord end point.
94 * @param[out] NewSafety New calculated safety distance.
95 * @param[in,out] PreviousSafety Previous safety distance.
96 * @param[in,out] PreviousSftOrigin Previous safety point origin.
97 * @param[out] LinearStepLength Linear chord length.
98 * @param[out] IntersectionPoint Intersection point.
99 * @param[in,out] calledNavigator Pointer to flag indicating if the
100 * navigator has been called or not.
101 * @returns Whether intersection exists or not.
102 */
103 inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
104 const G4ThreeVector& EndPointB,
105 G4double& NewSafety,
106 G4double& PreviousSafety, // In/Out
107 G4ThreeVector& PreviousSftOrigin, // In/Out
108 G4double& LinearStepLength,
109 G4ThreeVector& IntersectionPoint,
110 G4bool* calledNavigator = nullptr );
111
112 /**
113 * Setters for parameters which must be set at each step, in case they are
114 * changed.
115 * @note This simple approach ensures that all scenarios are considered.
116 * Future refinement may identify which are invariant during a
117 * track, run or event.
118 */
119 inline void SetEpsilonStepFor( G4double EpsilonStep );
120 inline void SetDeltaIntersectionFor( G4double deltaIntersection );
121 inline void SetNavigatorFor( G4Navigator* fNavigator );
122 inline void SetChordFinderFor(G4ChordFinder* fCFinder );
123
124 /**
125 * Verbosity control.
126 * Controlling verbosity enables checking of the locating of intersections.
127 */
128 inline void SetVerboseFor(G4int fVerbose);
130
131 /**
132 * Additional inline Get/Set methods for parameters, dependent objects.
133 */
138 inline void SetSafetyParametersFor(G4bool UseSafety );
139
140 /**
141 * Adjustment flag accessor/modifier.
142 */
143 inline void AdjustIntersections(G4bool UseCorrection);
145
146 /**
147 * Adjustment flag accessor/modifier.
148 * @deprecated Replaced by methods above. To be removed in future releases.
149 */
150 inline void AddAdjustementOfFoundIntersection(G4bool UseCorrection);
152
153 /**
154 * Dumps status of propagator to any ostream.
155 */
156 static void printStatus( const G4FieldTrack& startFT,
157 const G4FieldTrack& currentFT,
158 G4double requestStep,
159 G4double safety,
160 G4int stepNum,
161 std::ostream& oss,
162 G4int verboseLevel );
163
164 /**
165 * Dumps status of propagator to cout, useful mostly for debugging.
166 */
167 void printStatus( const G4FieldTrack& startFT,
168 const G4FieldTrack& currentFT,
169 G4double requestStep,
170 G4double safety,
171 G4int stepNum);
172
173 /**
174 * Sets/gets check mode.
175 * When enabled, uses additional verifications and stricter condictions
176 * for ensuring correctness. Effective only when G4VERBOSE is enabled.
177 */
178 inline void SetCheckMode( G4bool value ) { fCheckMode = value; }
179 inline G4bool GetCheckMode() { return fCheckMode; }
180
181 protected:
182
183 /**
184 * Returns new estimate for state after curveDist starting from
185 * CurrentStateA, to replace EstimtdEndStateB, and reports displacement
186 * (if field is compiled verbose).
187 * @param[in] CurrentStateA Start point tangent track.
188 * @param[in] EstimtdEndStateB Estimated end point tangent track.
189 * @param[in] linearDistSq Not used.
190 * @param[in] curveDist Not used.
191 * @returns New estimate for state.
192 */
193 G4FieldTrack ReEstimateEndpoint( const G4FieldTrack& CurrentStateA,
194 const G4FieldTrack& EstimtdEndStateB,
195 G4double linearDistSq, // not used
196 G4double curveDist ); // not used
197
198 /**
199 * Checks whether EndB is too far from StartA to be reached and if,
200 * re-estimates new value for EndB (return in RevisedEndPoint).
201 * Reports error if EndB is before StartA (in curve length)
202 * In that case return errorCode = 2.
203 * @param[in] CurrentStartA Start point tangent track.
204 * @param[in] EstimatedEndB Estimated end point tangent track.
205 * @param[in] RevisedEndPoint Revised end point tangent track.
206 * @param[in] errorCode Error code (0=OK, 1=coincident points, 2=error).
207 * @returns true if end point has been revised.
208 */
209 G4bool CheckAndReEstimateEndpoint( const G4FieldTrack& CurrentStartA,
210 const G4FieldTrack& EstimatedEndB,
211 G4FieldTrack& RevisedEndPoint,
212 G4int& errorCode);
213
214 /**
215 * Returns the surface normal. Position *must* be the intersection point
216 * from last call to G4Navigator's ComputeStep() (via IntersectChord).
217 * It tries to use cached (last) value in Navigator for speed, if it was
218 * kept and valid. The value returned is in global coordinates.
219 * @note It does NOT guarantee to obtain Normal. This can happen e.g.
220 * if the "Intersection" Point is not on a surface, potentially
221 * due to either inaccuracies in the transformations used, or
222 * issues with the Solid.
223 * @param[in] CurrentInt_Point Current point.
224 * @param[in,out] validNormal Flagging if normal is a valid vector.
225 * @returns The surface normal vector in local coordinates.
226 */
227 G4ThreeVector GetSurfaceNormal(const G4ThreeVector& CurrentInt_Point,
228 G4bool& validNormal);
229
230 /**
231 * Returns the surface normal of the Intersecting Solid in global
232 * coordinates.
233 * @note This method is costlier then GetSurfaceNormal().
234 * @param[in] CurrentInt_Point Current point.
235 * @param[in,out] validNormal Flagging if normal is a valid vector.
236 * @returns The surface normal vector in global coordinates.
237 */
239 G4bool& validNormal);
240 /**
241 * Optional method for adjustment of located intersection point using
242 * the surface-normal.
243 * @param[in] A Chord starting point.
244 * @param[in] CurrentE_Point E Chord point.
245 * @param[in] CurrentF_Point F Chord point.
246 * @param[in] MomentumDir Momentum direction.
247 * @param[in] IntersectAF First part intersecting?
248 * @param[in,out] IntersectionPoint Intersection point tangent track.
249 * @param[in,out] NewSafety New safety distance.
250 * @param[in,out] fPrevSafety Previous safety distance.
251 * @param[in,out] fPrevSftOrigin Previous safety point origin.
252 * @returns Whether intersection exists or not.
253 */
255 const G4ThreeVector& CurrentE_Point,
256 const G4ThreeVector& CurrentF_Point,
257 const G4ThreeVector& MomentumDir,
258 const G4bool IntersectAF,
259 G4ThreeVector& IntersectionPoint,
260 G4double& NewSafety,
261 G4double& fPrevSafety,
262 G4ThreeVector& fPrevSftOrigin );
263
264 /**
265 * Prints a three-line report on the current "sub-step",
266 * i.e. trial intersection.
267 * @param[in] step_no Step number.
268 * @param[in] ChordAB_v AB chord.
269 * @param[in] ChordEF_v EF chord.
270 * @param[in] NewMomentumDir Momentum direction.
271 * @param[in] NormalAtEntry Normal vector at entry.
272 * @param[in] validNormal Validity flag for normal vector at E.
273 */
274 void ReportTrialStep( G4int step_no,
275 const G4ThreeVector& ChordAB_v,
276 const G4ThreeVector& ChordEF_v,
277 const G4ThreeVector& NewMomentumDir,
278 const G4ThreeVector& NormalAtEntry,
279 G4bool validNormal );
280
281 /**
282 * Locates a point using the navigator and updates the state of Navigator.
283 * By default, it assumes that the point is inside the current volume,
284 * and returns true.
285 * In check mode, it checks whether the point is *inside* the volume.
286 * If it is inside, it returns true.
287 * If not, issues a warning and returns false.
288 * @param[in] pos The point to locate.
289 * @returns If a point is inside the volume or not.
290 */
292
293 /**
294 * Locates a point using the navigator and updates the state of Navigator,
295 * but report information about code location.
296 * If CheckMode > 1, report extra information.
297 * @param[in] pos The point to locate.
298 * @param[in] CodeLocationInfo String for code location info.
299 * @param[in] CheckMode Not used.
300 */
302 const G4String& CodeLocationInfo,
303 G4int CheckMode );
304
305 // Auxiliary methods -- to report issues
306
307 /**
308 * Builds error message (in ossMsg) to report that point 'B' has
309 * gone past 'A'.
310 */
311 void ReportReversedPoints( std::ostringstream& ossMsg,
312 const G4FieldTrack& StartPointVel,
313 const G4FieldTrack& EndPointVel,
314 G4double NewSafety, G4double epsStep,
315 const G4FieldTrack& CurrentA_PointVelocity,
316 const G4FieldTrack& CurrentB_PointVelocity,
317 const G4FieldTrack& SubStart_PointVelocity,
318 const G4ThreeVector& CurrentE_Point,
319 const G4FieldTrack& ApproxIntersecPointV,
320 G4int sbstp_no, G4int sbstp_no_p, G4int depth );
321
322 /**
323 * Reports the current status / progress in finding the first intersection.
324 */
325 void ReportProgress( std::ostream& oss,
326 const G4FieldTrack& StartPointVel,
327 const G4FieldTrack& EndPointVel,
328 G4int substep_no,
329 const G4FieldTrack& A_PtVel, // G4double safetyA
330 const G4FieldTrack& B_PtVel,
331 G4double safetyLast,
332 G4int depth= -1 );
333
334 /**
335 * Report case: trial point is 'close' to start, within tolerance.
336 */
337 void ReportImmediateHit( const char* MethodName,
338 const G4ThreeVector& StartPosition,
339 const G4ThreeVector& TrialPoint,
340 G4double tolerance,
341 unsigned long int numCalls );
342
343 private:
344
345 /**
346 * Returns the SurfaceNormal of the Intersecting Solid in local coordinates.
347 */
348 G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector& CurrentE_Point,
349 G4bool& validNormal);
350
351 /**
352 * Position *must* be the intersection point from last call
353 * to G4Navigator's ComputeStep (via IntersectChord).
354 */
355 G4ThreeVector GetLastSurfaceNormal( const G4ThreeVector& intersectPoint,
356 G4bool& validNormal) const;
357
358 protected:
359
361
362 G4int fVerboseLevel = 0; // For debugging
363 G4bool fUseNormalCorrection = false; // Configuration parameter
365 G4bool fiUseSafety = false; // Whether to use safety for 'fast steps'
366
368
369 /**
370 * Parameters set at each physical step by G4PropagatorInField.
371 */
372 G4ChordFinder* fiChordFinder = nullptr; // Overridden at each step
373 G4double fiEpsilonStep = -1.0; // Overridden at each step
374 G4double fiDeltaIntersection = -1.0; // Overridden at each step
375
376 /** Helper for location. */
378
379 /** Touchable history hook. */
381};
382
383#include "G4VIntersectionLocator.icc"
384
385#endif
#define fPreviousSftOrigin
#define fPreviousSafety
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4ChordFinder is a class that provides Runge-Kutta integration of motion ODE and also has a method th...
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
G4Navigator is a class for use by the tracking management, able to obtain/calculate dynamic tracking ...
G4TouchableHistory is an object representing a touchable detector element, and its history in the geo...
G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
G4Navigator * GetNavigatorFor()
void AddAdjustementOfFoundIntersection(G4bool UseCorrection)
G4ChordFinder * GetChordFinderFor()
G4VIntersectionLocator(G4Navigator *theNavigator)
G4TouchableHistory * fpTouchable
void AdjustIntersections(G4bool UseCorrection)
void SetDeltaIntersectionFor(G4double deltaIntersection)
G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point, G4bool &validNormal)
void SetVerboseFor(G4int fVerbose)
void SetNavigatorFor(G4Navigator *fNavigator)
static void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int stepNum, std::ostream &oss, G4int verboseLevel)
void LocateGlobalPointWithinVolumeCheckAndReport(const G4ThreeVector &pos, const G4String &CodeLocationInfo, G4int CheckMode)
void ReportTrialStep(G4int step_no, const G4ThreeVector &ChordAB_v, const G4ThreeVector &ChordEF_v, const G4ThreeVector &NewMomentumDir, const G4ThreeVector &NormalAtEntry, G4bool validNormal)
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &PreviousSafety, G4ThreeVector &PreviousSftOrigin, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint, G4bool *calledNavigator=nullptr)
G4bool CheckAndReEstimateEndpoint(const G4FieldTrack &CurrentStartA, const G4FieldTrack &EstimatedEndB, G4FieldTrack &RevisedEndPoint, G4int &errorCode)
virtual G4bool EstimateIntersectionPoint(const G4FieldTrack &curveStartPointTangent, const G4FieldTrack &curveEndPointTangent, const G4ThreeVector &trialPoint, G4FieldTrack &intersectPointTangent, G4bool &recalculatedEndPoint, G4double &fPreviousSafety, G4ThreeVector &fPreviousSftOrigin)=0
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)
G4double GetEpsilonStepFor()
void ReportImmediateHit(const char *MethodName, const G4ThreeVector &StartPosition, const G4ThreeVector &TrialPoint, G4double tolerance, unsigned long int numCalls)
G4FieldTrack ReEstimateEndpoint(const G4FieldTrack &CurrentStateA, const G4FieldTrack &EstimtdEndStateB, G4double linearDistSq, G4double curveDist)
G4bool GetAdjustementOfFoundIntersection()
G4bool LocateGlobalPointWithinVolumeAndCheck(const G4ThreeVector &pos)
void SetSafetyParametersFor(G4bool UseSafety)
G4double GetDeltaIntersectionFor()
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 SetChordFinderFor(G4ChordFinder *fCFinder)
void SetEpsilonStepFor(G4double EpsilonStep)
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)