Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BrentLocator.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 G4BrentLocator implementation
27//
28// Author: Tatiana Nikitina (CERN), 27 October 2008
29// --------------------------------------------------------------------
30
31#include <iomanip>
32
33#include "G4BrentLocator.hh"
34#include "G4ios.hh"
35
37 : G4VIntersectionLocator(theNavigator)
38{
39 // In case of too slow progress in finding Intersection Point
40 // intermediates Points on the Track must be stored.
41 // Initialise the array of Pointers [max_depth+1] to do this
42
43 G4ThreeVector zeroV(0.0,0.0,0.0);
44 for (auto & idepth : ptrInterMedFT)
45 {
46 idepth = new G4FieldTrack( zeroV, zeroV, 0., 0., 0., 0.);
47 }
48}
49
51{
52 for (auto & idepth : ptrInterMedFT)
53 {
54 delete idepth;
55 }
56}
57
58// --------------------------------------------------------------------------
59// G4bool G4PropagatorInField::LocateIntersectionPoint(
60// const G4FieldTrack& CurveStartPointVelocity, // A
61// const G4FieldTrack& CurveEndPointVelocity, // B
62// const G4ThreeVector& TrialPoint, // E
63// G4FieldTrack& IntersectedOrRecalculated // Output
64// G4bool& recalculated) // Out
65// --------------------------------------------------------------------------
66//
67// Function that returns the intersection of the true path with the surface
68// of the current volume (either the external one or the inner one with one
69// of the daughters:
70//
71// A = Initial point
72// B = another point
73//
74// Both A and B are assumed to be on the true path:
75//
76// E is the first point of intersection of the chord AB with
77// a volume other than A (on the surface of A or of a daughter)
78//
79// Convention of Use :
80// i) If it returns "true", then IntersectionPointVelocity is set
81// to the approximate intersection point.
82// ii) If it returns "false", no intersection was found.
83// The validity of IntersectedOrRecalculated depends on 'recalculated'
84// a) if latter is false, then IntersectedOrRecalculated is invalid.
85// b) if latter is true, then IntersectedOrRecalculated is
86// the new endpoint, due to a re-integration.
87// --------------------------------------------------------------------------
88// NOTE: implementation taken from G4PropagatorInField
89// New second order locator is added
90//
92 const G4FieldTrack& CurveStartPointVelocity, // A
93 const G4FieldTrack& CurveEndPointVelocity, // B
94 const G4ThreeVector& TrialPoint, // E
95 G4FieldTrack& IntersectedOrRecalculatedFT, // Output
96 G4bool& recalculatedEndPoint, // Out
97 G4double& fPreviousSafety, // In/Out
99
100{
101 // Find Intersection Point ( A, B, E ) of true path AB - start at E.
102
103 G4bool found_approximate_intersection = false;
104 G4bool there_is_no_intersection = false;
105
106 G4FieldTrack CurrentA_PointVelocity = CurveStartPointVelocity;
107 G4FieldTrack CurrentB_PointVelocity = CurveEndPointVelocity;
108 G4ThreeVector CurrentE_Point = TrialPoint;
109 G4bool validNormalAtE = false;
110 G4ThreeVector NormalAtEntry;
111
112 G4FieldTrack ApproxIntersecPointV(CurveEndPointVelocity); // FT-Def-Construct
113 G4double NewSafety = 0.0;
114 G4bool last_AF_intersection = false;
115
116 // G4bool final_section= true; // Shows whether current section is last
117 // (i.e. B=full end)
118 G4bool first_section = true;
119 recalculatedEndPoint = false;
120
121 G4bool restoredFullEndpoint = false;
122
123 G4long oldprc; // cout, cerr precision
124 G4int substep_no = 0;
125
126 // Limits for substep number
127 //
128 const G4int max_substeps= 10000; // Test 120 (old value 100 )
129 const G4int warn_substeps= 1000; // 100
130
131 // Statistics for substeps
132 //
133 static G4ThreadLocal G4int max_no_seen= -1;
134
135 // Counter for restarting Bintermed
136 //
137 G4int restartB = 0;
138
139 //--------------------------------------------------------------------------
140 // Algorithm for the case if progress in founding intersection is too slow.
141 // Process is defined too slow if after N=param_substeps advances on the
142 // path, it will be only 'fraction_done' of the total length.
143 // In this case the remaining length is divided in two half and
144 // the loop is restarted for each half.
145 // If progress is still too slow, the division in two halfs continue
146 // until 'max_depth'.
147 //--------------------------------------------------------------------------
148
149 const G4int param_substeps = 50; // Test value for the maximum number
150 // of substeps
151 const G4double fraction_done = 0.3;
152
153 G4bool Second_half = false; // First half or second half of divided step
154
155 NormalAtEntry = GetSurfaceNormal(CurrentE_Point, validNormalAtE);
156
157 // We need to know this for the 'final_section':
158 // real 'final_section' or first half 'final_section'
159 // In algorithm it is considered that the 'Second_half' is true
160 // and it becomes false only if we are in the first-half of level
161 // depthness or if we are in the first section
162
163 G4int depth = 0; // Depth counts how many subdivisions of initial step made
164
165#ifdef G4DEBUG_FIELD
166 const G4double tolerance = 1.0e-8;
167 G4ThreeVector StartPosition = CurveStartPointVelocity.GetPosition();
168 if( (TrialPoint - StartPosition).mag() < tolerance * CLHEP::mm )
169 {
170 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
171 "GeomNav1002", JustWarning,
172 "Intersection point F is exactly at start point A." );
173 }
174#endif
175
176 // Intermediates Points on the Track = Subdivided Points must be stored.
177 // Give the initial values to 'InterMedFt'
178 // Important is 'ptrInterMedFT[0]', it saves the 'EndCurvePoint'
179 //
180 *ptrInterMedFT[0] = CurveEndPointVelocity;
181 for (auto idepth=1; idepth<max_depth+1; ++idepth )
182 {
183 *ptrInterMedFT[idepth] = CurveStartPointVelocity;
184 }
185
186 //Final_section boolean store
187 G4bool fin_section_depth[max_depth];
188 for (bool & idepth : fin_section_depth)
189 {
190 idepth = true;
191 }
192
193 // 'SubStartPoint' is needed to calculate the length of the divided step
194 //
195 G4FieldTrack SubStart_PointVelocity = CurveStartPointVelocity;
196
197 do // Loop checking, 07.10.2016, J.Apostolakis
198 {
199 G4int substep_no_p = 0;
200 G4bool sub_final_section = false; // the same as final_section,
201 // but for 'sub_section'
202 SubStart_PointVelocity = CurrentA_PointVelocity;
203
204 do // Loop checking, 07.10.2016, J.Apostolakis
205 { // REPEAT param
206 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
207 G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
208
209 // F = a point on true AB path close to point E
210 // (the closest if possible)
211 //
212 if(substep_no_p==0)
213 {
214 ApproxIntersecPointV = GetChordFinderFor()
215 ->ApproxCurvePointV( CurrentA_PointVelocity,
216 CurrentB_PointVelocity,
217 CurrentE_Point,
219 // The above method is the key & most intuitive part ...
220 }
221#ifdef G4DEBUG_FIELD
222 if( ApproxIntersecPointV.GetCurveLength() >
223 CurrentB_PointVelocity.GetCurveLength() * (1.0 + tolerance) )
224 {
225 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
226 "GeomNav0003", FatalException,
227 "Intermediate F point is past end B point" );
228 }
229#endif
230
231 G4ThreeVector CurrentF_Point = ApproxIntersecPointV.GetPosition();
232
233 // First check whether EF is small - then F is a good approx. point
234 // Calculate the length and direction of the chord AF
235 //
236 G4ThreeVector ChordEF_Vector = CurrentF_Point - CurrentE_Point;
237 G4ThreeVector NewMomentumDir = ApproxIntersecPointV.GetMomentumDir();
238 G4double MomDir_dot_Norm = NewMomentumDir.dot( NormalAtEntry ) ;
239
240#ifdef G4DEBUG_FIELD
241 G4ThreeVector ChordAB = Point_B - Point_A;
242
243 G4VIntersectionLocator::ReportTrialStep( substep_no, ChordAB,
244 ChordEF_Vector, NewMomentumDir, NormalAtEntry, validNormalAtE );
245#endif
246
247 G4bool adequate_angle;
248 adequate_angle = ( MomDir_dot_Norm >= 0.0 ) // Can use -epsilon instead.
249 || (! validNormalAtE ); // Makes criterion invalid
250 G4double EF_dist2 = ChordEF_Vector.mag2();
251 if ( ( EF_dist2 <= sqr(fiDeltaIntersection) && ( adequate_angle ) )
252 || ( EF_dist2 <= kCarTolerance*kCarTolerance ) )
253 {
254 found_approximate_intersection = true;
255
256 // Create the "point" return value
257 //
258 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
259 IntersectedOrRecalculatedFT.SetPosition( CurrentE_Point );
260
262 {
263 // Try to Get Correction of IntersectionPoint using SurfaceNormal()
264 //
265 G4ThreeVector IP;
266 G4ThreeVector MomentumDir=ApproxIntersecPointV.GetMomentumDirection();
267 G4bool goodCorrection = AdjustmentOfFoundIntersection( Point_A,
268 CurrentE_Point, CurrentF_Point, MomentumDir,
269 last_AF_intersection, IP, NewSafety,
271 if ( goodCorrection )
272 {
273 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
274 IntersectedOrRecalculatedFT.SetPosition(IP);
275 }
276 }
277
278 // Note: in order to return a point on the boundary,
279 // we must return E. But it is F on the curve.
280 // So we must "cheat": we are using the position at point E
281 // and the velocity at point F !!!
282 //
283 // This must limit the length we can allow for displacement!
284 }
285 else // E is NOT close enough to the curve (ie point F)
286 {
287 // Check whether any volumes are encountered by the chord AF
288 // ---------------------------------------------------------
289 // First relocate to restore any Voxel etc information
290 // in the Navigator before calling ComputeStep()
291 //
293
294 G4ThreeVector PointG; // Candidate intersection point
295 G4double stepLengthAF;
296 G4bool usedNavigatorAF = false;
297 G4bool Intersects_AF = IntersectChord( Point_A, CurrentF_Point,
298 NewSafety,fPreviousSafety,
300 stepLengthAF,
301 PointG,
302 &usedNavigatorAF);
303 last_AF_intersection = Intersects_AF;
304 if( Intersects_AF )
305 {
306 // G is our new Candidate for the intersection point.
307 // It replaces "E" and we will repeat the test to see if
308 // it is a good enough approximate point for us.
309 // B <- F
310 // E <- G
311 //
312 G4FieldTrack EndPoint = ApproxIntersecPointV;
313 ApproxIntersecPointV = GetChordFinderFor()->ApproxCurvePointS(
314 CurrentA_PointVelocity, CurrentB_PointVelocity,
315 EndPoint,CurrentE_Point, CurrentF_Point,PointG,
316 true, GetEpsilonStepFor() );
317 CurrentB_PointVelocity = EndPoint;
318 CurrentE_Point = PointG;
319
320 // Need to recalculate the Exit Normal at the new PointG
321 // Know that a call was made to Navigator::ComputeStep in
322 // IntersectChord above.
323 //
324 G4bool validNormalLast;
325 NormalAtEntry = GetSurfaceNormal( PointG, validNormalLast );
326 validNormalAtE = validNormalLast;
327
328 // By moving point B, must take care if current
329 // AF has no intersection to try current FB!!
330 //
331 fin_section_depth[depth] = false;
332#ifdef G4VERBOSE
333 if( fVerboseLevel > 3 )
334 {
335 G4cout << "G4PiF::LI> Investigating intermediate point"
336 << " at s=" << ApproxIntersecPointV.GetCurveLength()
337 << " on way to full s="
338 << CurveEndPointVelocity.GetCurveLength() << G4endl;
339 }
340#endif
341 }
342 else // not Intersects_AF
343 {
344 // In this case:
345 // There is NO intersection of AF with a volume boundary.
346 // We must continue the search in the segment FB!
347 //
349
350 G4double stepLengthFB;
351 G4ThreeVector PointH;
352 G4bool usedNavigatorFB = false;
353
354 // Check whether any volumes are encountered by the chord FB
355 // ---------------------------------------------------------
356
357 G4bool Intersects_FB = IntersectChord( CurrentF_Point, Point_B,
358 NewSafety,fPreviousSafety,
360 stepLengthFB,
361 PointH,
362 &usedNavigatorFB);
363 if( Intersects_FB )
364 {
365 // There is an intersection of FB with a volume boundary
366 // H <- First Intersection of Chord FB
367
368 // H is our new Candidate for the intersection point.
369 // It replaces "E" and we will repeat the test to see if
370 // it is a good enough approximate point for us.
371
372 // Note that F must be in volume volA (the same as A)
373 // (otherwise AF would meet a volume boundary!)
374 // A <- F
375 // E <- H
376 //
377 G4FieldTrack InterMed = ApproxIntersecPointV;
378 ApproxIntersecPointV = GetChordFinderFor()->ApproxCurvePointS(
379 CurrentA_PointVelocity,CurrentB_PointVelocity,
380 InterMed,CurrentE_Point,CurrentF_Point,PointH,
381 false,GetEpsilonStepFor());
382 CurrentA_PointVelocity = InterMed;
383 CurrentE_Point = PointH;
384
385 // Need to recalculate the Exit Normal at the new PointG
386 //
387 G4bool validNormalLast;
388 NormalAtEntry = GetSurfaceNormal( PointH, validNormalLast );
389 validNormalAtE = validNormalLast;
390 }
391 else // not Intersects_FB
392 {
393 // There is NO intersection of FB with a volume boundary
394
395 if( fin_section_depth[depth] )
396 {
397 // If B is the original endpoint, this means that whatever
398 // volume(s) intersected the original chord, none touch the
399 // smaller chords we have used.
400 // The value of 'IntersectedOrRecalculatedFT' returned is
401 // likely not valid
402
403 // Check on real final_section or SubEndSection
404 //
405 if( ((Second_half)&&(depth==0)) || (first_section) )
406 {
407 there_is_no_intersection = true; // real final_section
408 }
409 else
410 {
411 // end of subsection, not real final section
412 // exit from the and go to the depth-1 level
413
414 substep_no_p = param_substeps+2; // exit from the loop
415
416 // but 'Second_half' is still true because we need to find
417 // the 'CurrentE_point' for the next loop
418 //
419 Second_half = true;
420 sub_final_section = true;
421 }
422 }
423 else
424 {
425 if( depth==0 )
426 {
427 // We must restore the original endpoint
428 //
429 CurrentA_PointVelocity = CurrentB_PointVelocity; // Got to B
430 CurrentB_PointVelocity = CurveEndPointVelocity;
431 SubStart_PointVelocity = CurrentA_PointVelocity;
432 ApproxIntersecPointV = GetChordFinderFor()
433 ->ApproxCurvePointV( CurrentA_PointVelocity,
434 CurrentB_PointVelocity,
435 CurrentE_Point,
437
438 restoredFullEndpoint = true;
439 ++restartB; // counter
440 }
441 else
442 {
443 // We must restore the depth endpoint
444 //
445 CurrentA_PointVelocity = CurrentB_PointVelocity; // Got to B
446 CurrentB_PointVelocity = *ptrInterMedFT[depth];
447 SubStart_PointVelocity = CurrentA_PointVelocity;
448 ApproxIntersecPointV = GetChordFinderFor()
449 ->ApproxCurvePointV( CurrentA_PointVelocity,
450 CurrentB_PointVelocity,
451 CurrentE_Point,
453 restoredFullEndpoint = true;
454 ++restartB; // counter
455 }
456 }
457 } // Endif (Intersects_FB)
458 } // Endif (Intersects_AF)
459
460 // Ensure that the new endpoints are not further apart in space
461 // than on the curve due to different errors in the integration
462 //
463 G4double linDistSq, curveDist;
464 linDistSq = ( CurrentB_PointVelocity.GetPosition()
465 - CurrentA_PointVelocity.GetPosition() ).mag2();
466 curveDist = CurrentB_PointVelocity.GetCurveLength()
467 - CurrentA_PointVelocity.GetCurveLength();
468
469 // Change this condition for very strict parameters of propagation
470 //
471 if( curveDist*curveDist*(1+2* GetEpsilonStepFor()) < linDistSq )
472 {
473 // Re-integrate to obtain a new B
474 //
475 G4FieldTrack newEndPointFT=
476 ReEstimateEndpoint( CurrentA_PointVelocity,
477 CurrentB_PointVelocity,
478 linDistSq, // to avoid recalculation
479 curveDist );
480 G4FieldTrack oldPointVelB = CurrentB_PointVelocity;
481 CurrentB_PointVelocity = newEndPointFT;
482
483 if ( (fin_section_depth[depth]) // real final section
484 &&( first_section || ((Second_half)&&(depth==0)) ) )
485 {
486 recalculatedEndPoint = true;
487 IntersectedOrRecalculatedFT = newEndPointFT;
488 // So that we can return it, if it is the endpoint!
489 }
490 }
491 if( curveDist < 0.0 )
492 {
493 fVerboseLevel = 5; // Print out a maximum of information
494 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
495 -1.0, NewSafety, substep_no );
496 std::ostringstream message;
497 message << "Error in advancing propagation." << G4endl
498 << " Error in advancing propagation." << G4endl
499 << " Point A (start) is " << CurrentA_PointVelocity
500 << G4endl
501 << " Point B (end) is " << CurrentB_PointVelocity
502 << G4endl
503 << " Curve distance is " << curveDist << G4endl
504 << G4endl
505 << "The final curve point is not further along"
506 << " than the original!" << G4endl;
507
508 if( recalculatedEndPoint )
509 {
510 message << "Recalculation of EndPoint was called with fEpsStep= "
512 }
513 oldprc = G4cerr.precision(20);
514 message << " Point A (Curve start) is " << CurveStartPointVelocity
515 << G4endl
516 << " Point B (Curve end) is " << CurveEndPointVelocity
517 << G4endl
518 << " Point A (Current start) is " << CurrentA_PointVelocity
519 << G4endl
520 << " Point B (Current end) is " << CurrentB_PointVelocity
521 << G4endl
522 << " Point S (Sub start) is " << SubStart_PointVelocity
523 << G4endl
524 << " Point E (Trial Point) is " << CurrentE_Point
525 << G4endl
526 << " Old Point F(Intersection) is " << CurrentF_Point
527 << G4endl
528 << " New Point F(Intersection) is " << ApproxIntersecPointV
529 << G4endl
530 << " LocateIntersection parameters are : Substep no= "
531 << substep_no << G4endl
532 << " Substep depth no= "<< substep_no_p << " Depth= "
533 << depth << G4endl
534 << " Restarted no= "<< restartB << " Epsilon= "
535 << GetEpsilonStepFor() <<" DeltaInters= "
537 G4cerr.precision( oldprc );
538
539 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
540 "GeomNav0003", FatalException, message);
541 }
542
543 if( restoredFullEndpoint )
544 {
545 fin_section_depth[depth] = restoredFullEndpoint;
546 restoredFullEndpoint = false;
547 }
548 } // EndIf ( E is close enough to the curve, ie point F. )
549 // tests ChordAF_Vector.mag() <= maximum_lateral_displacement
550
551#ifdef G4DEBUG_LOCATE_INTERSECTION
552 G4int trigger_substepno_print= warn_substeps - 20 ;
553
554 if( substep_no >= trigger_substepno_print )
555 {
556 G4cout << "Difficulty in converging in "
557 << "G4BrentLocator::EstimateIntersectionPoint()"
558 << G4endl
559 << " Substep no = " << substep_no << G4endl;
560 if( substep_no == trigger_substepno_print )
561 {
562 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
563 -1.0, NewSafety, 0);
564 }
565 G4cout << " State of point A: ";
566 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
567 -1.0, NewSafety, substep_no-1, 0);
568 G4cout << " State of point B: ";
569 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
570 -1.0, NewSafety, substep_no);
571 }
572#endif
573 ++substep_no;
574 ++substep_no_p;
575
576 } while ( ( ! found_approximate_intersection )
577 && ( ! there_is_no_intersection )
578 && ( substep_no_p <= param_substeps) ); // UNTIL found or
579 // failed param substep
580 first_section = false;
581
582 if( (!found_approximate_intersection) && (!there_is_no_intersection) )
583 {
584 G4double did_len = std::abs( CurrentA_PointVelocity.GetCurveLength()
585 - SubStart_PointVelocity.GetCurveLength());
586 G4double all_len = std::abs( CurrentB_PointVelocity.GetCurveLength()
587 - SubStart_PointVelocity.GetCurveLength());
588
589 G4double stepLengthAB;
590 G4ThreeVector PointGe;
591
592 // Check if progress is too slow and if it possible to go deeper,
593 // then halve the step if so
594 //
595 if ( ( did_len < fraction_done*all_len )
596 && (depth < max_depth) && (!sub_final_section) )
597 {
598 Second_half=false;
599 ++depth;
600
601 G4double Sub_len = (all_len-did_len)/(2.);
602 G4FieldTrack start = CurrentA_PointVelocity;
603 auto integrDriver =
605 integrDriver->AccurateAdvance(start, Sub_len, GetEpsilonStepFor());
606 *ptrInterMedFT[depth] = start;
607 CurrentB_PointVelocity = *ptrInterMedFT[depth];
608
609 // Adjust 'SubStartPoint' to calculate the 'did_length' in next loop
610 //
611 SubStart_PointVelocity = CurrentA_PointVelocity;
612
613 // Find new trial intersection point needed at start of the loop
614 //
615 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
616 G4ThreeVector SubE_point = CurrentB_PointVelocity.GetPosition();
617
619 G4bool Intersects_AB = IntersectChord(Point_A, SubE_point,
620 NewSafety, fPreviousSafety,
621 fPreviousSftOrigin,stepLengthAB,
622 PointGe);
623 if( Intersects_AB )
624 {
625 last_AF_intersection = Intersects_AB;
626 CurrentE_Point = PointGe;
627 fin_section_depth[depth] = true;
628
629 // Need to recalculate the Exit Normal at the new PointG
630 //
631 G4bool validNormalAB;
632 NormalAtEntry = GetSurfaceNormal( PointGe, validNormalAB );
633 validNormalAtE = validNormalAB;
634 }
635 else
636 {
637 // No intersection found for first part of curve
638 // (CurrentA,InterMedPoint[depth]). Go to the second part
639 //
640 Second_half = true;
641 }
642 } // if did_len
643
644 if( (Second_half)&&(depth!=0) )
645 {
646 // Second part of curve (InterMed[depth],Intermed[depth-1]) )
647 // On the depth-1 level normally we are on the 'second_half'
648
649 Second_half = true;
650
651 // Find new trial intersection point needed at start of the loop
652 //
653 SubStart_PointVelocity = *ptrInterMedFT[depth];
654 CurrentA_PointVelocity = *ptrInterMedFT[depth];
655 CurrentB_PointVelocity = *ptrInterMedFT[depth-1];
656 // Ensure that the new endpoints are not further apart in space
657 // than on the curve due to different errors in the integration
658 //
659 G4double linDistSq, curveDist;
660 linDistSq = ( CurrentB_PointVelocity.GetPosition()
661 - CurrentA_PointVelocity.GetPosition() ).mag2();
662 curveDist = CurrentB_PointVelocity.GetCurveLength()
663 - CurrentA_PointVelocity.GetCurveLength();
664 if( curveDist*curveDist*(1+2*GetEpsilonStepFor() ) < linDistSq )
665 {
666 // Re-integrate to obtain a new B
667 //
668 G4FieldTrack newEndPointFT =
669 ReEstimateEndpoint( CurrentA_PointVelocity,
670 CurrentB_PointVelocity,
671 linDistSq, // to avoid recalculation
672 curveDist );
673 G4FieldTrack oldPointVelB = CurrentB_PointVelocity;
674 CurrentB_PointVelocity = newEndPointFT;
675 if ( depth==1 )
676 {
677 recalculatedEndPoint = true;
678 IntersectedOrRecalculatedFT = newEndPointFT;
679 // So that we can return it, if it is the endpoint!
680 }
681 }
682
683
684 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
685 G4ThreeVector SubE_point = CurrentB_PointVelocity.GetPosition();
687 G4bool Intersects_AB = IntersectChord(Point_A, SubE_point, NewSafety,
689 fPreviousSftOrigin,stepLengthAB, PointGe);
690 if( Intersects_AB )
691 {
692 last_AF_intersection = Intersects_AB;
693 CurrentE_Point = PointGe;
694
695 G4bool validNormalAB;
696 NormalAtEntry = GetSurfaceNormal( PointGe, validNormalAB );
697 validNormalAtE = validNormalAB;
698 }
699
700 depth--;
701 fin_section_depth[depth]=true;
702 }
703 } // if(!found_aproximate_intersection)
704
705 } while ( ( ! found_approximate_intersection )
706 && ( ! there_is_no_intersection )
707 && ( substep_no <= max_substeps) ); // UNTIL found or failed
708
709 if( substep_no > max_no_seen )
710 {
711 max_no_seen = substep_no;
712#ifdef G4DEBUG_LOCATE_INTERSECTION
713 if( max_no_seen > warn_substeps )
714 {
715 trigger_substepno_print = max_no_seen-20; // Want to see last 20 steps
716 }
717#endif
718 }
719
720 if( ( substep_no >= max_substeps)
721 && !there_is_no_intersection
722 && !found_approximate_intersection )
723 {
724 G4cout << "ERROR - G4BrentLocator::EstimateIntersectionPoint()" << G4endl
725 << " Start and end-point of requested Step:" << G4endl;
726 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
727 -1.0, NewSafety, 0);
728 G4cout << " Start and end-point of current Sub-Step:" << G4endl;
729 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
730 -1.0, NewSafety, substep_no-1);
731 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
732 -1.0, NewSafety, substep_no);
733 std::ostringstream message;
734 message << "Too many substeps!" << G4endl
735 << " Convergence is requiring too many substeps: "
736 << substep_no << G4endl
737 << " Abandoning effort to intersect. " << G4endl
738 << " Found intersection = "
739 << found_approximate_intersection << G4endl
740 << " Intersection exists = "
741 << !there_is_no_intersection << G4endl;
742 oldprc = G4cout.precision( 10 );
743 G4double done_len = CurrentA_PointVelocity.GetCurveLength();
744 G4double full_len = CurveEndPointVelocity.GetCurveLength();
745 message << " Undertaken only length: " << done_len
746 << " out of " << full_len << " required." << G4endl
747 << " Remaining length = " << full_len - done_len;
748 G4cout.precision( oldprc );
749
750 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
751 "GeomNav0003", FatalException, message);
752 }
753 else if( substep_no >= warn_substeps )
754 {
755 oldprc = G4cout.precision( 10 );
756 std::ostringstream message;
757 message << "Many substeps while trying to locate intersection."
758 << G4endl
759 << " Undertaken length: "
760 << CurrentB_PointVelocity.GetCurveLength()
761 << " - Needed: " << substep_no << " substeps." << G4endl
762 << " Warning level = " << warn_substeps
763 << " and maximum substeps = " << max_substeps;
764 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
765 "GeomNav1002", JustWarning, message);
766 G4cout.precision( oldprc );
767 }
768 return !there_is_no_intersection; // Success or failure
769}
@ 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
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double mag2() const
double dot(const Hep3Vector &) const
G4bool EstimateIntersectionPoint(const G4FieldTrack &curveStartPointTangent, const G4FieldTrack &curveEndPointTangent, const G4ThreeVector &trialPoint, G4FieldTrack &intersectPointTangent, G4bool &recalculatedEndPoint, G4double &fPreviousSafety, G4ThreeVector &fPreviousSftOrigin) override
G4BrentLocator(G4Navigator *theNavigator)
~G4BrentLocator() override
G4FieldTrack ApproxCurvePointV(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4ThreeVector &currentEPoint, G4double epsStep)
G4FieldTrack ApproxCurvePointS(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4FieldTrack &ApproxCurveV, const G4ThreeVector &currentEPoint, const G4ThreeVector &currentFPoint, const G4ThreeVector &PointG, G4bool first, G4double epsStep)
G4VIntegrationDriver * GetIntegrationDriver()
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)
virtual G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0)=0
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()
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)
T sqr(const T &x)
Definition templates.hh:128
#define G4ThreadLocal
Definition tls.hh:77