Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SimpleLocator.cc
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// Class G4SimpleLocator implementation
27//
28// Author: Tatiana Nikitina (CERN), 27 October 20008.
29// ---------------------------------------------------------------------------
30
31#include <iomanip>
32
33#include "G4ios.hh"
34#include "G4SimpleLocator.hh"
35
40
42
43// --------------------------------------------------------------------------
44// G4bool G4PropagatorInField::LocateIntersectionPoint(
45// const G4FieldTrack& CurveStartPointVelocity, // A
46// const G4FieldTrack& CurveEndPointVelocity, // B
47// const G4ThreeVector& TrialPoint, // E
48// G4FieldTrack& IntersectedOrRecalculated // Output
49// G4bool& recalculated ) // Out
50// --------------------------------------------------------------------------
51//
52// Function that returns the intersection of the true path with the surface
53// of the current volume (either the external one or the inner one with one
54// of the daughters:
55//
56// A = Initial point
57// B = another point
58//
59// Both A and B are assumed to be on the true path:
60//
61// E is the first point of intersection of the chord AB with
62// a volume other than A (on the surface of A or of a daughter)
63//
64// Convention of Use :
65// i) If it returns "true", then IntersectionPointVelocity is set
66// to the approximate intersection point.
67// ii) If it returns "false", no intersection was found.
68// The validity of IntersectedOrRecalculated depends on 'recalculated'
69// a) if latter is false, then IntersectedOrRecalculated is invalid.
70// b) if latter is true, then IntersectedOrRecalculated is
71// the new endpoint, due to a re-integration.
72// --------------------------------------------------------------------------
73// NOTE: implementation taken from G4PropagatorInField
74//
76 const G4FieldTrack& CurveStartPointVelocity, // A
77 const G4FieldTrack& CurveEndPointVelocity, // B
78 const G4ThreeVector& TrialPoint, // E
79 G4FieldTrack& IntersectedOrRecalculatedFT, // Output
80 G4bool& recalculatedEndPoint, // Out
81 G4double &fPreviousSafety, //In/Out
83{
84 // Find Intersection Point ( A, B, E ) of true path AB - start at E.
85
86 G4bool found_approximate_intersection = false;
87 G4bool there_is_no_intersection = false;
88
89 G4FieldTrack CurrentA_PointVelocity = CurveStartPointVelocity;
90 G4FieldTrack CurrentB_PointVelocity = CurveEndPointVelocity;
91 G4ThreeVector CurrentE_Point = TrialPoint;
92 G4bool validNormalAtE = false;
93 G4ThreeVector NormalAtEntry;
94
95 G4FieldTrack ApproxIntersecPointV(CurveEndPointVelocity); // FT-Def-Construct
96 G4double NewSafety = 0.0;
97 G4bool last_AF_intersection = false;
98 G4bool final_section = true; // Shows whether current section is last
99 // (i.e. B=full end)
100 recalculatedEndPoint = false;
101
102 G4bool restoredFullEndpoint = false;
103
104 G4int substep_no = 0;
105
106 // Limits for substep number
107 //
108 const G4int max_substeps = 100000000; // Test 120 (old value 100 )
109 const G4int warn_substeps = 1000; // 100
110
111 // Statistics for substeps
112 //
113 static G4ThreadLocal G4int max_no_seen= -1;
114
115 NormalAtEntry = GetSurfaceNormal( CurrentE_Point, validNormalAtE);
116
117#ifdef G4DEBUG_FIELD
118 const G4double tolerance = 1.0e-8;
119 G4ThreeVector StartPosition= CurveStartPointVelocity.GetPosition();
120 if( (TrialPoint - StartPosition).mag() < tolerance * CLHEP::mm )
121 {
122 G4Exception("G4SimpleLocator::EstimateIntersectionPoint()",
123 "GeomNav1002", JustWarning,
124 "Intersection point F is exactly at start point A." );
125 }
126#endif
127
128 do
129 {
130 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
131 G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
132
133 // F = a point on true AB path close to point E
134 // (the closest if possible)
135 //
136 ApproxIntersecPointV = GetChordFinderFor()
137 ->ApproxCurvePointV( CurrentA_PointVelocity,
138 CurrentB_PointVelocity,
139 CurrentE_Point,
141 // The above method is the key & most intuitive part ...
142
143#ifdef G4DEBUG_FIELD
144 if( ApproxIntersecPointV.GetCurveLength() >
145 CurrentB_PointVelocity.GetCurveLength() * (1.0 + tolerance) )
146 {
147 G4Exception("G4SimpleLocator::EstimateIntersectionPoint()",
148 "GeomNav0003", FatalException,
149 "Intermediate F point is past end B point!" );
150 }
151#endif
152
153 G4ThreeVector CurrentF_Point = ApproxIntersecPointV.GetPosition();
154
155 // First check whether EF is small - then F is a good approx. point
156 // Calculate the length and direction of the chord AF
157 //
158 G4ThreeVector ChordEF_Vector = CurrentF_Point - CurrentE_Point;
159
160 G4ThreeVector NewMomentumDir = ApproxIntersecPointV.GetMomentumDir();
161 G4double MomDir_dot_Norm = NewMomentumDir.dot( NormalAtEntry ) ;
162
163 G4ThreeVector ChordAB = Point_B - Point_A;
164
165#ifdef G4DEBUG_FIELD
167 ReportTrialStep( substep_no, ChordAB, ChordEF_Vector,
168 NewMomentumDir, NormalAtEntry, validNormalAtE );
169#endif
170 // Check Sign is always exiting !! TODO
171 // Could ( > -epsilon) be used instead?
172 //
173 G4bool adequate_angle = ( MomDir_dot_Norm >= 0.0 )
174 || (! validNormalAtE ); // Invalid
175 G4double EF_dist2= ChordEF_Vector.mag2();
176 if ( ( EF_dist2 <= sqr(fiDeltaIntersection) && ( adequate_angle ) )
177 || ( EF_dist2 <= kCarTolerance*kCarTolerance ) )
178 {
179 found_approximate_intersection = true;
180
181 // Create the "point" return value
182 //
183 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
184 IntersectedOrRecalculatedFT.SetPosition( CurrentE_Point );
185
187 {
188 // Try to Get Correction of IntersectionPoint using SurfaceNormal()
189 //
190 G4ThreeVector IP;
191 G4ThreeVector MomentumDir= ApproxIntersecPointV.GetMomentumDirection();
192 G4bool goodCorrection = AdjustmentOfFoundIntersection( Point_A,
193 CurrentE_Point, CurrentF_Point, MomentumDir,
194 last_AF_intersection, IP, NewSafety,
196
197 if ( goodCorrection )
198 {
199 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
200 IntersectedOrRecalculatedFT.SetPosition(IP);
201 }
202 }
203
204 // Note: in order to return a point on the boundary,
205 // we must return E. But it is F on the curve.
206 // So we must "cheat": we are using the position at point E
207 // and the velocity at point F !!!
208 //
209 // This must limit the length we can allow for displacement!
210 }
211 else // E is NOT close enough to the curve (ie point F)
212 {
213 // Check whether any volumes are encountered by the chord AF
214 // ---------------------------------------------------------
215 // First relocate to restore any Voxel etc information
216 // in the Navigator before calling ComputeStep()
217 //
219
220 G4ThreeVector PointG; // Candidate intersection point
221 G4double stepLengthAF;
222 G4bool usedNavigatorAF = false;
223 G4bool Intersects_AF = IntersectChord( Point_A,
224 CurrentF_Point,
225 NewSafety,
228 stepLengthAF,
229 PointG,
230 &usedNavigatorAF );
231 last_AF_intersection = Intersects_AF;
232 if( Intersects_AF )
233 {
234 // G is our new Candidate for the intersection point.
235 // It replaces "E" and we will repeat the test to see if
236 // it is a good enough approximate point for us.
237 // B <- F
238 // E <- G
239
240 CurrentB_PointVelocity = ApproxIntersecPointV;
241 CurrentE_Point = PointG;
242
243 // Need to recalculate the Exit Normal at the new PointG
244 // Relies on a call to Navigator::ComputeStep in IntersectChord above
245 // If the safety was adequate (for the step) this would NOT be called!
246 // But this must not occur, no intersection can be found in that case,
247 // so this branch, ie if( Intersects_AF ) would not be reached!
248 //
249 G4bool validNormalLast;
250 NormalAtEntry = GetSurfaceNormal( PointG, validNormalLast );
251 validNormalAtE = validNormalLast;
252
253 // By moving point B, must take care if current
254 // AF has no intersection to try current FB!!
255 //
256 final_section = false;
257
258#ifdef G4VERBOSE
259 if( fVerboseLevel > 3 )
260 {
261 G4cout << "G4PiF::LI> Investigating intermediate point"
262 << " at s=" << ApproxIntersecPointV.GetCurveLength()
263 << " on way to full s="
264 << CurveEndPointVelocity.GetCurveLength() << G4endl;
265 }
266#endif
267 }
268 else // not Intersects_AF
269 {
270 // In this case:
271 // There is NO intersection of AF with a volume boundary.
272 // We must continue the search in the segment FB!
273 //
275
276 G4double stepLengthFB;
277 G4ThreeVector PointH;
278 G4bool usedNavigatorFB = false;
279
280 // Check whether any volumes are encountered by the chord FB
281 // ---------------------------------------------------------
282
283 G4bool Intersects_FB = IntersectChord( CurrentF_Point, Point_B,
284 NewSafety,fPreviousSafety,
286 stepLengthFB,
287 PointH, &usedNavigatorFB );
288 if( Intersects_FB )
289 {
290 // There is an intersection of FB with a volume boundary
291 // H <- First Intersection of Chord FB
292
293 // H is our new Candidate for the intersection point.
294 // It replaces "E" and we will repeat the test to see if
295 // it is a good enough approximate point for us.
296
297 // Note that F must be in volume volA (the same as A)
298 // (otherwise AF would meet a volume boundary!)
299 // A <- F
300 // E <- H
301 //
302 CurrentA_PointVelocity = ApproxIntersecPointV;
303 CurrentE_Point = PointH;
304
305 // Need to recalculate the Exit Normal at the new PointG
306 // Relies on call to Navigator::ComputeStep in IntersectChord above
307 // If safety was adequate (for the step) this would NOT be called!
308 // But this must not occur, no intersection found in that case,
309 // so this branch, i.e. if( Intersects_AF ) would not be reached!
310 //
311 G4bool validNormalLast;
312 NormalAtEntry = GetSurfaceNormal( PointH, validNormalLast );
313 validNormalAtE = validNormalLast;
314 }
315 else // not Intersects_FB
316 {
317 // There is NO intersection of FB with a volume boundary
318
319 if( final_section )
320 {
321 // If B is the original endpoint, this means that whatever
322 // volume(s) intersected the original chord, none touch the
323 // smaller chords we have used.
324 // The value of 'IntersectedOrRecalculatedFT' returned is
325 // likely not valid
326
327 there_is_no_intersection = true; // real final_section
328 }
329 else
330 {
331 // We must restore the original endpoint
332
333 CurrentA_PointVelocity = CurrentB_PointVelocity; // Got to B
334 CurrentB_PointVelocity = CurveEndPointVelocity;
335 restoredFullEndpoint = true;
336 }
337 } // Endif (Intersects_FB)
338 } // Endif (Intersects_AF)
339
340 // Ensure that the new endpoints are not further apart in space
341 // than on the curve due to different errors in the integration
342 //
343 G4double linDistSq, curveDist;
344 linDistSq = ( CurrentB_PointVelocity.GetPosition()
345 - CurrentA_PointVelocity.GetPosition() ).mag2();
346 curveDist = CurrentB_PointVelocity.GetCurveLength()
347 - CurrentA_PointVelocity.GetCurveLength();
348
349 // Change this condition for very strict parameters of propagation
350 //
351 if( curveDist*curveDist*(1+2* GetEpsilonStepFor()) < linDistSq )
352 {
353 // Re-integrate to obtain a new B
354 //
355 G4FieldTrack newEndPointFT =
356 ReEstimateEndpoint( CurrentA_PointVelocity,
357 CurrentB_PointVelocity,
358 linDistSq, // to avoid recalculation
359 curveDist );
360 G4FieldTrack oldPointVelB = CurrentB_PointVelocity;
361 CurrentB_PointVelocity = newEndPointFT;
362
363 if( (final_section)) // real final section
364 {
365 recalculatedEndPoint = true;
366 IntersectedOrRecalculatedFT = newEndPointFT;
367 // So that we can return it, if it is the endpoint!
368 }
369 }
370 if( curveDist < 0.0 )
371 {
372 fVerboseLevel = 5; // Print out a maximum of information
373 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
374 -1.0, NewSafety, substep_no );
375 std::ostringstream message;
376 message << "Error in advancing propagation." << G4endl
377 << " Point A (start) is " << CurrentA_PointVelocity
378 << G4endl
379 << " Point B (end) is " << CurrentB_PointVelocity
380 << G4endl
381 << " Curve distance is " << curveDist << G4endl
382 << G4endl
383 << "The final curve point is not further along"
384 << " than the original!" << G4endl;
385
386 if( recalculatedEndPoint )
387 {
388 message << "Recalculation of EndPoint was called with fEpsStep= "
390 }
391 message.precision(20);
392 message << " Point A (Curve start) is " << CurveStartPointVelocity
393 << G4endl
394 << " Point B (Curve end) is " << CurveEndPointVelocity
395 << G4endl
396 << " Point A (Current start) is " << CurrentA_PointVelocity
397 << G4endl
398 << " Point B (Current end) is " << CurrentB_PointVelocity
399 << G4endl
400 << " Point E (Trial Point) is " << CurrentE_Point
401 << G4endl
402 << " Point F (Intersection) is " << ApproxIntersecPointV
403 << G4endl
404 << " LocateIntersection parameters are : Substep no= "
405 << substep_no;
406
407 G4Exception("G4SimpleLocator::EstimateIntersectionPoint()",
408 "GeomNav0003", FatalException, message);
409 }
410
411 if ( restoredFullEndpoint )
412 {
413 final_section = restoredFullEndpoint;
414 restoredFullEndpoint = false;
415 }
416 } // EndIf ( E is close enough to the curve, ie point F. )
417 // tests ChordAF_Vector.mag() <= maximum_lateral_displacement
418
419#ifdef G4DEBUG_LOCATE_INTERSECTION
420 G4int trigger_substepno_print= warn_substeps - 20;
421
422 if( substep_no >= trigger_substepno_print )
423 {
424 G4cout << "Difficulty in converging in "
425 << "G4SimpleLocator::EstimateIntersectionPoint():"
426 << G4endl
427 << " Substep no = " << substep_no << G4endl;
428 if( substep_no == trigger_substepno_print )
429 {
430 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
431 -1.0, NewSafety, 0);
432 }
433 G4cout << " State of point A: ";
434 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
435 -1.0, NewSafety, substep_no-1, 0);
436 G4cout << " State of point B: ";
437 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
438 -1.0, NewSafety, substep_no);
439 }
440#endif
441 ++substep_no;
442
443 } while ( ( ! found_approximate_intersection )
444 && ( ! there_is_no_intersection )
445 && ( substep_no <= max_substeps) ); // UNTIL found or failed
446
447 if( substep_no > max_no_seen )
448 {
449 max_no_seen = substep_no;
450#ifdef G4DEBUG_LOCATE_INTERSECTION
451 if( max_no_seen > warn_substeps )
452 {
453 trigger_substepno_print = max_no_seen-20; // Want to see last 20 steps
454 }
455#endif
456 }
457
458 if( ( substep_no >= max_substeps)
459 && !there_is_no_intersection
460 && !found_approximate_intersection )
461 {
462 G4cout << "ERROR - G4SimpleLocator::EstimateIntersectionPoint()" << G4endl
463 << " Start and Endpoint of Requested Step:" << G4endl;
464 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
465 -1.0, NewSafety, 0);
466 G4cout << G4endl
467 << " Start and end-point of current Sub-Step:" << G4endl;
468 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
469 -1.0, NewSafety, substep_no-1);
470 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
471 -1.0, NewSafety, substep_no);
472
473 std::ostringstream message;
474 message << "Convergence is requiring too many substeps: "
475 << substep_no << G4endl
476 << " Abandoning effort to intersect." << G4endl
477 << " Found intersection = "
478 << found_approximate_intersection << G4endl
479 << " Intersection exists = "
480 << !there_is_no_intersection << G4endl;
481 message.precision(10);
482 G4double done_len = CurrentA_PointVelocity.GetCurveLength();
483 G4double full_len = CurveEndPointVelocity.GetCurveLength();
484 message << " Undertaken only length: " << done_len
485 << " out of " << full_len << " required." << G4endl
486 << " Remaining length = " << full_len-done_len;
487
488 G4Exception("G4SimpleLocator::EstimateIntersectionPoint()",
489 "GeomNav0003", FatalException, message);
490 }
491 else if( substep_no >= warn_substeps )
492 {
493 std::ostringstream message;
494 message.precision(10);
495
496 message << "Many substeps while trying to locate intersection." << G4endl
497 << " Undertaken length: "
498 << CurrentB_PointVelocity.GetCurveLength()
499 << " - Needed: " << substep_no << " substeps." << G4endl
500 << " Warning level = " << warn_substeps
501 << " and maximum substeps = " << max_substeps;
502 G4Exception("G4SimpleLocator::EstimateIntersectionPoint()",
503 "GeomNav1002", JustWarning, message);
504 }
505 return !there_is_no_intersection; // Success or failure
506}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#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
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double mag2() const
double dot(const Hep3Vector &) const
G4FieldTrack ApproxCurvePointV(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4ThreeVector &currentEPoint, G4double epsStep)
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
const G4ThreeVector & GetMomentumDir() const
G4double GetCurveLength() const
G4ThreeVector GetMomentumDirection() const
void SetPosition(const G4ThreeVector &nPos)
G4ThreeVector GetPosition() const
G4Navigator is a class for use by the tracking management, able to obtain/calculate dynamic tracking ...
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
G4SimpleLocator(G4Navigator *aNavigator)
G4bool EstimateIntersectionPoint(const G4FieldTrack &curveStartPointTangent, const G4FieldTrack &curveEndPointTangent, const G4ThreeVector &trialPoint, G4FieldTrack &intersectPointTangent, G4bool &recalculatedEndPoint, G4double &fPreviousSafety, G4ThreeVector &fPreviousSftOrigin) override
~G4SimpleLocator() override
G4Navigator * GetNavigatorFor()
G4ChordFinder * GetChordFinderFor()
G4VIntersectionLocator(G4Navigator *theNavigator)
G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point, G4bool &validNormal)
static void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int stepNum, std::ostream &oss, G4int verboseLevel)
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)
G4double GetEpsilonStepFor()
G4FieldTrack ReEstimateEndpoint(const G4FieldTrack &CurrentStateA, const G4FieldTrack &EstimtdEndStateB, G4double linearDistSq, G4double curveDist)
G4bool GetAdjustementOfFoundIntersection()
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)
T sqr(const T &x)
Definition templates.hh:128
#define G4ThreadLocal
Definition tls.hh:77