Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/navigation/src/G4MultiLevelLocator.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 G4MultiLevelLocator implementation
 27 //
 28 // 27.10.08 - Tatiana Nikitina.
 29 // 04.10.11 - John Apostolakis, revised convergence to use Surface Normal
 30 // ---------------------------------------------------------------------------
 31 
 32 #include <iomanip>
 33 
 34 #include "G4ios.hh"
 35 #include "G4MultiLevelLocator.hh"
 36 #include "G4LocatorChangeRecord.hh"
 37 #include "G4LocatorChangeLogger.hh"
 38 
 39 G4MultiLevelLocator::G4MultiLevelLocator(G4Navigator *theNavigator)
 40   : G4VIntersectionLocator(theNavigator)
 41 {
 42   // In case of too slow progress in finding Intersection Point
 43   // intermediates Points on the Track must be stored.
 44   // Initialise the array of Pointers [max_depth+1] to do this  
 45   
 46   G4ThreeVector zeroV(0.0,0.0,0.0);
 47   for (auto & idepth : ptrInterMedFT)
 48   {
 49     idepth = new G4FieldTrack( zeroV, zeroV, 0., 0., 0., 0.);
 50   }
 51 
 52   if (fCheckMode)
 53   {
 54     //  Trial values       Loose         Medium      Tight
 55     //  To happen:         Infrequent    Occasional  Often
 56     SetMaxSteps(150);   //  300             85          25
 57     SetWarnSteps(80);   //  250             60          15
 58   }  
 59 }      
 60 
 61 G4MultiLevelLocator::~G4MultiLevelLocator()
 62 {
 63   for (auto & idepth : ptrInterMedFT)
 64   {
 65     delete idepth;
 66   }
 67 #ifdef G4DEBUG_FIELD
 68   ReportStatistics(); 
 69 #endif  
 70 }
 71 
 72 
 73 // --------------------------------------------------------------------------
 74 // G4bool G4PropagatorInField::LocateIntersectionPoint( 
 75 //        const G4FieldTrack&       CurveStartPointVelocity,   // A
 76 //        const G4FieldTrack&       CurveEndPointVelocity,     // B
 77 //        const G4ThreeVector&      TrialPoint,                // E
 78 //              G4FieldTrack&       IntersectedOrRecalculated  // Output
 79 //              G4bool&             recalculated )             // Out
 80 // --------------------------------------------------------------------------
 81 //
 82 // Function that returns the intersection of the true path with the surface
 83 // of the current volume (either the external one or the inner one with one
 84 // of the daughters:
 85 //
 86 //     A = Initial point
 87 //     B = another point 
 88 //
 89 // Both A and B are assumed to be on the true path:
 90 //
 91 //     E is the first point of intersection of the chord AB with 
 92 //     a volume other than A  (on the surface of A or of a daughter)
 93 //
 94 // Convention of Use :
 95 //     i) If it returns "true", then IntersectionPointVelocity is set
 96 //       to the approximate intersection point.
 97 //    ii) If it returns "false", no intersection was found.
 98 //        Potential reasons:
 99 //        a) no segment found an intersection
100 //        b) too many steps were required - after that it abandoned the effort
101 //           and is returning how far it could go.  (New - 29 Oct 2015)
102 //           (If so, it must set 'recalculated' to true.) 
103 //        TODO/idea: add a new flag: 'unfinished' to identify these cases.
104 //
105 //        IntersectedOrRecalculated means different things:
106 //        a) if it is the same curve lenght along, it is a revision of the
107 //            original enpdoint due to the need for re-integration.
108 //        b) if it is at a shorter curve length, it is 'end of what it could do'
109 //           i.e. as far as it could go, because it took too many steps!
110 //        Note: IntersectedOrRecalculated is valid only if 'recalculated' is
111 //       'true'.
112 // --------------------------------------------------------------------------
113 // NOTE: implementation taken from G4PropagatorInField
114 //
115 G4bool G4MultiLevelLocator::EstimateIntersectionPoint( 
116          const  G4FieldTrack&       CurveStartPointVelocity,       // A
117          const  G4FieldTrack&       CurveEndPointVelocity,         // B
118          const  G4ThreeVector&      TrialPoint,                    // E
119                 G4FieldTrack&       IntersectedOrRecalculatedFT,   // Output
120                 G4bool&             recalculatedEndPoint,          // Out
121                 G4double&           previousSafety,                // In/Out
122                 G4ThreeVector&      previousSftOrigin)             // In/Out
123 {
124   // Find Intersection Point ( A, B, E )  of true path AB - start at E.
125   const char* MethodName= "G4MultiLevelLocator::EstimateIntersectionPoint()";
126   
127   G4bool found_approximate_intersection = false;
128   G4bool there_is_no_intersection       = false;
129 
130   G4FieldTrack  CurrentA_PointVelocity = CurveStartPointVelocity; 
131   G4FieldTrack  CurrentB_PointVelocity = CurveEndPointVelocity;
132   G4ThreeVector CurrentE_Point = TrialPoint;
133   G4bool        validNormalAtE = false;
134   G4ThreeVector NormalAtEntry;
135 
136   G4FieldTrack  ApproxIntersecPointV(CurveEndPointVelocity); // FT-Def-Construct
137   G4bool        validIntersectP= true;   // Is it current ?
138   G4double      NewSafety = 0.0;
139   G4bool        last_AF_intersection = false;   
140 
141   auto    integrDriver = GetChordFinderFor()->GetIntegrationDriver();
142   G4bool  driverReIntegrates = integrDriver->DoesReIntegrate();
143   
144   G4bool first_section = true;
145   recalculatedEndPoint = false; 
146   G4bool restoredFullEndpoint = false;
147 
148   unsigned int substep_no = 0;
149 
150   // Statistics for substeps
151   static G4ThreadLocal unsigned int max_no_seen= 0;
152 
153 #ifdef G4DEBUG_FIELD
154   unsigned int trigger_substepno_print = 0;
155   const G4double tolerance = 1.0e-8 * CLHEP::mm;
156   unsigned int biggest_depth = 0;
157     // using kInitialisingCL = G4LocatorChangeRecord::kInitialisingCL;
158 #endif
159 
160   // Log the location, iteration of changes in A,B 
161   //----------------------------------------------
162   static G4ThreadLocal G4LocatorChangeLogger endChangeA("StartPointA"),
163     endChangeB("EndPointB"), recApproxPoint("ApproxPoint"),
164     pointH_logger("Trial points 'E': position, normal");
165   unsigned int eventCount = 0;
166 
167   if (fCheckMode)
168   {
169     // Clear previous call's data
170     endChangeA.clear();
171     endChangeB.clear();
172     recApproxPoint.clear();
173     pointH_logger.clear();
174 
175     // Record the initialisation
176     ++eventCount;
177     endChangeA.AddRecord( G4LocatorChangeRecord::kInitialisingCL, substep_no,
178                           eventCount, CurrentA_PointVelocity );
179     endChangeB.AddRecord( G4LocatorChangeRecord::kInitialisingCL, substep_no,
180                           eventCount, CurrentB_PointVelocity );
181   }
182 
183   //--------------------------------------------------------------------------  
184   //  Algorithm for the case if progress in founding intersection is too slow.
185   //  Process is defined too slow if after N=param_substeps advances on the
186   //  path, it will be only 'fraction_done' of the total length.
187   //  In this case the remaining length is divided in two half and 
188   //  the loop is restarted for each half.  
189   //  If progress is still too slow, the division in two halfs continue
190   //  until 'max_depth'.
191   //--------------------------------------------------------------------------
192 
193   const G4int param_substeps = 5;  // Test value for the maximum number
194                                    // of substeps
195   const G4double fraction_done = 0.3;
196 
197   G4bool Second_half = false;    // First half or second half of divided step
198 
199   // We need to know this for the 'final_section':
200   // real 'final_section' or first half 'final_section'
201   // In algorithm it is considered that the 'Second_half' is true
202   // and it becomes false only if we are in the first-half of level
203   // depthness or if we are in the first section
204 
205   unsigned int depth = 0; // Depth counts subdivisions of initial step made
206   ++fNumCalls;
207   
208   NormalAtEntry = GetSurfaceNormal(CurrentE_Point, validNormalAtE);
209 
210   if (fCheckMode)
211   {
212     pointH_logger.AddRecord( G4LocatorChangeRecord::kInitialisingCL,
213                              substep_no, eventCount,
214                              G4FieldTrack( CurrentE_Point,0.,NormalAtEntry,0.,
215                                            0., 1.,G4ThreeVector(0.),0.,0.) );
216   #if (G4DEBUG_FIELD>1) 
217     G4ThreeVector  StartPosition = CurveStartPointVelocity.GetPosition(); 
218     if( (TrialPoint - StartPosition).mag2() < sqr(tolerance)) 
219     {
220        ReportImmediateHit( MethodName, StartPosition, TrialPoint,
221                            tolerance, fNumCalls); 
222     }
223   #endif
224   }
225 
226   // Intermediates Points on the Track = Subdivided Points must be stored.
227   // Give the initial values to 'InterMedFt'
228   // Important is 'ptrInterMedFT[0]', it saves the 'EndCurvePoint'
229   //
230   *ptrInterMedFT[0] = CurveEndPointVelocity;
231   for ( auto idepth=1; idepth<max_depth+1; ++idepth )
232   {
233     *ptrInterMedFT[idepth] = CurveStartPointVelocity;
234   }
235 
236   // Final_section boolean store
237   //
238   G4bool fin_section_depth[max_depth];
239   for (bool & idepth : fin_section_depth)
240   {
241     idepth = true;
242   }
243   // 'SubStartPoint' is needed to calculate the length of the divided step
244   //
245   G4FieldTrack SubStart_PointVelocity = CurveStartPointVelocity;
246    
247   do  // Loop checking, 07.10.2016, J.Apostolakis
248   {
249     unsigned int substep_no_p = 0;
250     G4bool sub_final_section = false; // the same as final_section,
251                                       // but for 'sub_section'
252     SubStart_PointVelocity = CurrentA_PointVelocity;
253  
254     do // Loop checking, 07.10.2016, J.Apostolakis
255     { // REPEAT param
256       G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();  
257       G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
258        
259 #ifdef G4DEBUG_FIELD
260       const G4double lenA = CurrentA_PointVelocity.GetCurveLength() ;
261       const G4double lenB = CurrentB_PointVelocity.GetCurveLength() ;            
262       G4double curv_lenAB = lenB - lenA;
263       G4double     distAB = (Point_B - Point_A).mag();
264       if( curv_lenAB < distAB * ( 1. - 10.*fiEpsilonStep ) )
265       {
266         G4cerr << "ERROR> (Start) Point A coincides with or has gone past (end) point B"           
267                << "MLL: iters = " << substep_no << G4endl;
268         G4long op=G4cerr.precision(6);
269         G4cerr << "       Difference = " << distAB - curv_lenAB
270                << " exceeds limit of relative dist (10*epsilon)= " << 10*fiEpsilonStep
271                << "  i.e. limit = " << 10 * fiEpsilonStep * distAB << G4endl;
272         G4cerr.precision(9);        
273         G4cerr << "        Len A, B = " << lenA << " " << lenB << G4endl
274                << "        Position A: " << Point_A << G4endl
275                << "        Position B: " << Point_B << G4endl;
276         G4cerr.precision(op);
277         // G4LocatorChangeRecord::ReportVector(G4cerr, "endPointB", endChangeB );
278         // G4cerr<<"EndPoints A(start) and B(end): combined changes " << G4endl;
279         if (fCheckMode)
280         {
281           G4LocatorChangeLogger::ReportEndChanges(G4cerr, endChangeA, endChangeB);
282         }
283       }
284 
285       if( !validIntersectP )
286       {
287         G4ExceptionDescription errmsg;
288         errmsg << "Assertion FAILURE - invalid (stale) Interection point. Substep: "
289                << substep_no << " call: " << fNumCalls << G4endl;
290         if (fCheckMode)
291         {
292           G4LocatorChangeRecord::ReportEndChanges(errmsg, endChangeA, endChangeB );
293         }
294         G4Exception("G4MultiLevelLocator::EstimateIntersectionPoint", "GeomNav0004",
295                     JustWarning, errmsg);
296       }
297 #endif    
298 
299       // F = a point on true AB path close to point E 
300       // (the closest if possible)
301       //
302       ApproxIntersecPointV = GetChordFinderFor()
303                              ->ApproxCurvePointV( CurrentA_PointVelocity, 
304                                                   CurrentB_PointVelocity, 
305                                                   CurrentE_Point,
306                                                   GetEpsilonStepFor() );
307         // The above method is the key & most intuitive part ...
308 
309 #ifdef G4DEBUG_FIELD      
310       recApproxPoint.push_back(G4LocatorChangeRecord(G4LocatorChangeRecord::kInvalidCL,
311                                                   substep_no, eventCount, ApproxIntersecPointV ) );
312       G4double lenIntsc= ApproxIntersecPointV.GetCurveLength();
313       G4double checkVsEnd= lenB - lenIntsc;
314 
315       if( lenIntsc > lenB ) 
316       {
317         std::ostringstream errmsg;
318         errmsg.precision(17);
319         G4double ratio    = checkVsEnd / lenB;
320         G4double ratioTol = std::fabs(ratio) / tolerance;
321         errmsg << "Intermediate F point is past end B point" << G4endl
322         << "   l( intersection ) = " << lenIntsc << G4endl
323         << "   l( endpoint     ) = " << lenB     << G4endl;
324         errmsg.precision(8);
325         errmsg << "   l_end - l_inters  = " << checkVsEnd << G4endl
326         << "          / l_end      = " << ratio << G4endl
327         << "   ratio  / tolerance  = " << ratioTol  << G4endl;
328         if( ratioTol < 1.0 )
329           G4Exception(MethodName, "GeomNav0003", JustWarning, errmsg );
330         else
331           G4Exception(MethodName, "GeomNav0003", FatalException, errmsg );
332       }
333 #endif
334 
335       G4ThreeVector CurrentF_Point= ApproxIntersecPointV.GetPosition();
336 
337       // First check whether EF is small - then F is a good approx. point 
338       // Calculate the length and direction of the chord AF
339       //
340       G4ThreeVector  ChordEF_Vector = CurrentF_Point - CurrentE_Point;
341 
342       G4ThreeVector  NewMomentumDir = ApproxIntersecPointV.GetMomentumDir(); 
343       G4double       MomDir_dot_Norm = NewMomentumDir.dot( NormalAtEntry );
344 
345 #ifdef G4DEBUG_FIELD
346       if( fVerboseLevel > 3 )
347       { 
348          G4ThreeVector  ChordAB           = Point_B - Point_A;
349          G4double       ABchord_length    = ChordAB.mag(); 
350          G4double       MomDir_dot_ABchord;
351          MomDir_dot_ABchord = (1.0 / ABchord_length)
352                             * NewMomentumDir.dot( ChordAB );
353          G4VIntersectionLocator::ReportTrialStep( substep_no, ChordAB,
354               ChordEF_Vector, NewMomentumDir, NormalAtEntry, validNormalAtE ); 
355          G4cout << " dot( MomentumDir, ABchord_unit ) = "
356                 << MomDir_dot_ABchord << G4endl;
357       }
358 #endif
359       G4bool adequate_angle =
360              ( MomDir_dot_Norm >= 0.0 ) // Can use ( > -epsilon) instead
361           || (! validNormalAtE );       // Invalid, cannot use this criterion
362       G4double EF_dist2 = ChordEF_Vector.mag2();
363       if ( ( EF_dist2 <= sqr(fiDeltaIntersection) && ( adequate_angle ) )
364         || ( EF_dist2 <= kCarTolerance*kCarTolerance ) )
365       { 
366         found_approximate_intersection = true;
367 
368         // Create the "point" return value
369         //
370         IntersectedOrRecalculatedFT = ApproxIntersecPointV;
371         IntersectedOrRecalculatedFT.SetPosition( CurrentE_Point );
372 
373         if ( GetAdjustementOfFoundIntersection() )
374         {
375           // Try to Get Correction of IntersectionPoint using SurfaceNormal()
376           //  
377           G4ThreeVector IP;
378           G4ThreeVector MomentumDir=ApproxIntersecPointV.GetMomentumDirection();
379           G4bool goodCorrection = AdjustmentOfFoundIntersection(Point_A,
380                                     CurrentE_Point, CurrentF_Point, MomentumDir,
381                                     last_AF_intersection, IP, NewSafety,
382                                     previousSafety, previousSftOrigin );
383           if ( goodCorrection )
384           {
385             IntersectedOrRecalculatedFT = ApproxIntersecPointV;
386             IntersectedOrRecalculatedFT.SetPosition(IP);
387           }
388         }
389         // Note: in order to return a point on the boundary, 
390         //       we must return E. But it is F on the curve.
391         //       So we must "cheat": we are using the position at point E
392         //       and the velocity at point F !!!
393         //
394         // This must limit the length we can allow for displacement!
395       }
396       else  // E is NOT close enough to the curve (ie point F)
397       {
398         // Check whether any volumes are encountered by the chord AF
399         // ---------------------------------------------------------
400         // First relocate to restore any Voxel etc information
401         // in the Navigator before calling ComputeStep()
402         //
403         GetNavigatorFor()->LocateGlobalPointWithinVolume( Point_A );
404 
405         G4ThreeVector PointG;   // Candidate intersection point
406         G4double stepLengthAF; 
407         G4bool Intersects_FB = false;
408         G4bool Intersects_AF = IntersectChord( Point_A,   CurrentF_Point,
409                                                NewSafety, previousSafety,
410                                                previousSftOrigin,
411                                                stepLengthAF,
412                                                PointG );
413         last_AF_intersection = Intersects_AF;
414         if( Intersects_AF )
415         {
416           // G is our new Candidate for the intersection point.
417           // It replaces  "E" and we will see if it's good enough.
418           CurrentB_PointVelocity = ApproxIntersecPointV;  // B  <- F
419           CurrentE_Point = PointG;                        // E  <- G
420 
421           validIntersectP = true;    // 'E' has been updated.
422 
423           G4bool validNormalLast; 
424           NormalAtEntry  = GetSurfaceNormal( PointG, validNormalLast ); 
425           validNormalAtE = validNormalLast; 
426 
427           // As we move point B, must take care in case the current
428           // AF has no intersection to try current FB!!
429           fin_section_depth[depth] = false;
430 
431           if (fCheckMode)
432           {
433             ++eventCount;
434             endChangeB.push_back(
435               G4LocatorChangeRecord(G4LocatorChangeRecord::kIntersectsAF,
436                              substep_no, eventCount, CurrentB_PointVelocity) );
437           }
438 #ifdef G4VERBOSE
439           if( fVerboseLevel > 3 )
440           {
441             G4cout << "G4PiF::LI> Investigating intermediate point"
442                    << " at s=" << ApproxIntersecPointV.GetCurveLength()
443                    << " on way to full s="
444                    << CurveEndPointVelocity.GetCurveLength() << G4endl;
445           }
446 #endif
447         }
448         else  // not Intersects_AF
449         {  
450           // In this case:
451           // There is NO intersection of AF with a volume boundary.
452           // We must continue the search in the segment FB!
453           //
454           GetNavigatorFor()->LocateGlobalPointWithinVolume( CurrentF_Point );
455 
456           G4double stepLengthFB;
457           G4ThreeVector PointH;
458 
459           // Check whether any volumes are encountered by the chord FB
460           // ---------------------------------------------------------
461 
462           Intersects_FB = IntersectChord( CurrentF_Point, Point_B, 
463                                           NewSafety, previousSafety,
464                                           previousSftOrigin,
465                                           stepLengthFB,
466                                           PointH );
467           if( Intersects_FB )
468           {
469             // There is an intersection of FB with a volume boundary
470             // H <- First Intersection of Chord FB 
471 
472             // H is our new Candidate for the intersection point.
473             // It replaces  "E" and we will repeat the test to see if
474             // it is a good enough approximate point for us.
475 
476             // Note that F must be in volume volA  (the same as A)
477             // (otherwise AF would meet a volume boundary!)
478             //   A    <- F 
479             //   E    <- H
480             //
481             CurrentA_PointVelocity = ApproxIntersecPointV;
482             CurrentE_Point = PointH;
483 
484             validIntersectP = true;    // 'E' has been updated.
485 
486             G4bool validNormalH;
487             NormalAtEntry  = GetSurfaceNormal( PointH, validNormalH ); 
488             validNormalAtE = validNormalH;
489 
490             if (fCheckMode)
491             {
492               ++eventCount;
493               endChangeA.push_back( 
494                  G4LocatorChangeRecord(G4LocatorChangeRecord::kIntersectsFB,
495                              substep_no, eventCount, CurrentA_PointVelocity) );
496               G4FieldTrack intersectH_pn('0');  // Point and normal
497                                                 // nothing else will be valid
498               intersectH_pn.SetPosition( PointH );
499               intersectH_pn.SetMomentum( NormalAtEntry ); 
500               pointH_logger.AddRecord(G4LocatorChangeRecord::kIntersectsFB,
501                                     substep_no, eventCount, intersectH_pn );
502             }
503           }
504           else  // not Intersects_FB
505           {
506             validIntersectP = false;    // Intersections are now stale
507             if( fin_section_depth[depth] )
508             {
509               // If B is the original endpoint, this means that whatever
510               // volume(s) intersected the original chord, none touch the
511               // smaller chords we have used.
512               // The value of 'IntersectedOrRecalculatedFT' returned is
513               // likely not valid 
514 
515               // Check on real final_section or SubEndSection
516               //
517               if( ((Second_half)&&(depth==0)) || (first_section) )
518               {
519                 there_is_no_intersection = true;   // real final_section
520               }
521               else
522               {
523                 // end of subsection, not real final section 
524                 // exit from the and go to the depth-1 level 
525                 substep_no_p = param_substeps+2;  // exit from the loop
526 
527                 // but 'Second_half' is still true because we need to find
528                 // the 'CurrentE_point' for the next loop
529                 Second_half = true; 
530                 sub_final_section = true;
531               }
532             }
533             else
534             {
535               CurrentA_PointVelocity = CurrentB_PointVelocity;  // Got to B
536               CurrentB_PointVelocity = (depth==0) ? CurveEndPointVelocity 
537                                                   : *ptrInterMedFT[depth] ;
538               SubStart_PointVelocity = CurrentA_PointVelocity;
539               restoredFullEndpoint = true;
540 
541               validIntersectP = false;  // 'E' has NOT been updated.
542 
543               if (fCheckMode)
544               {
545                 ++eventCount;
546                 endChangeA.push_back(
547                   G4LocatorChangeRecord(
548                     G4LocatorChangeRecord::kNoIntersectAForFB,
549                              substep_no, eventCount, CurrentA_PointVelocity) ); 
550                 endChangeB.push_back(
551                   G4LocatorChangeRecord (
552                     G4LocatorChangeRecord::kNoIntersectAForFB,
553                              substep_no, eventCount, CurrentB_PointVelocity) );
554               }
555             }
556           } // Endif (Intersects_FB)
557         } // Endif (Intersects_AF)
558 
559         G4int errorEndPt = 0; // Default: no error (if not calling CheckAnd...
560 
561         G4bool recalculatedB= false;
562         if( driverReIntegrates )
563         {
564           G4FieldTrack RevisedB_FT = CurrentB_PointVelocity;
565           recalculatedB= CheckAndReEstimateEndpoint(CurrentA_PointVelocity,
566                                                     CurrentB_PointVelocity,
567                                                     RevisedB_FT,
568                                                     errorEndPt );
569           if( recalculatedB )
570           {
571             CurrentB_PointVelocity = RevisedB_FT;  // Use it !
572             // Do not invalidate intersection F -- it is still roughly OK.
573             //
574             // The best course would be to invalidate (reset validIntersectP)
575             // BUT if we invalidate it, we must re-estimate it somewhere! E.g.
576             //   validIntersectP = false;  // 'E' has NOT been updated.
577  
578             if ( (fin_section_depth[depth])           // real final section
579                &&( first_section || ((Second_half)&&(depth==0)) ) )
580             {
581               recalculatedEndPoint = true;
582               IntersectedOrRecalculatedFT = RevisedB_FT;
583               // So that we can return it, if it is the endpoint!
584             }
585             // else
586             //  Move forward the other points 
587             //   - or better flag it, so that they are re-computed when next used
588             //     [ Implementation: a counter for # of recomputations
589             //       => avoids extra work]
590           }
591           if (fCheckMode)
592           {
593             ++eventCount;
594             endChangeB.push_back(
595               G4LocatorChangeRecord( G4LocatorChangeRecord::kRecalculatedB,
596                                      substep_no, eventCount, RevisedB_FT ) );
597           }
598         }
599         else
600         {
601           if( CurrentB_PointVelocity.GetCurveLength()
602             < CurrentA_PointVelocity.GetCurveLength() )
603           {
604             errorEndPt = 2;
605           }
606         }
607 
608         if( errorEndPt > 1 )  // errorEndPt = 1 is milder, just: len(B)=len(A)
609         {
610           std::ostringstream errmsg;
611 
612           ReportReversedPoints(errmsg,
613                                CurveStartPointVelocity, CurveEndPointVelocity,
614                                NewSafety, fiEpsilonStep,
615                                CurrentA_PointVelocity, CurrentB_PointVelocity,
616                                SubStart_PointVelocity, CurrentE_Point,
617                                ApproxIntersecPointV, substep_no, substep_no_p, depth);
618           
619           if (fCheckMode)
620           {
621             G4LocatorChangeRecord::ReportEndChanges(errmsg, endChangeA, endChangeB );
622           }
623 
624           errmsg << G4endl << " * Location: " << MethodName
625                  << "- After EndIf(Intersects_AF)" << G4endl;
626           errmsg << " * Bool flags:  Recalculated = " << recalculatedB
627                  << "   Intersects_AF = " << Intersects_AF
628                  << "   Intersects_FB = " << Intersects_FB << G4endl;
629           errmsg << " * Number of calls to MLL:EIP= " << fNumCalls << G4endl;
630           G4Exception(MethodName, "GeomNav0003", FatalException, errmsg);
631         }
632         if( restoredFullEndpoint )
633         {
634           fin_section_depth[depth] = restoredFullEndpoint;
635           restoredFullEndpoint = false;
636         }
637       } // EndIf ( E is close enough to the curve, ie point F. )
638         // tests ChordAF_Vector.mag() <= maximum_lateral_displacement 
639 
640 #ifdef G4DEBUG_FIELD
641       if( trigger_substepno_print == 0)
642       {
643         trigger_substepno_print= fWarnSteps - 20;
644       }
645 
646       if( substep_no >= trigger_substepno_print )
647       {
648         G4cout << "Difficulty in converging in " << MethodName
649                << "    Substep no = " << substep_no << G4endl;
650         if( substep_no == trigger_substepno_print )
651         {
652           G4cout << " Start:   ";
653           printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
654                        -1.0, NewSafety, 0 );
655           if( fCheckMode ) {
656             G4LocatorChangeRecord::ReportEndChanges(G4cout, endChangeA, endChangeB );
657           } else {
658             G4cout << " ** For more information enable 'check mode' in G4MultiLevelLocator "
659                    << "-- (it saves and can output change records) " << G4endl;
660           }
661         }
662         G4cout << " Point A: "; 
663         printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
664                      -1.0, NewSafety, substep_no-1 );
665         G4cout << " Point B: ";
666         printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
667                      -1.0, NewSafety, substep_no );
668       }
669 #endif
670       ++substep_no; 
671       ++substep_no_p;
672 
673     } while (  ( ! found_approximate_intersection )
674             && ( ! there_is_no_intersection )     
675             && validIntersectP        // New condition:  must refresh intersection !!
676             && ( substep_no_p <= param_substeps) );  // UNTIL found or
677                                                      // failed param substep
678 
679     if( (!found_approximate_intersection) && (!there_is_no_intersection) )
680     {
681       G4double did_len = std::abs( CurrentA_PointVelocity.GetCurveLength()
682                        - SubStart_PointVelocity.GetCurveLength()); 
683       G4double all_len = std::abs( CurrentB_PointVelocity.GetCurveLength()
684                        - SubStart_PointVelocity.GetCurveLength());
685    
686       G4double distAB = -1;
687       //
688       // Is progress is too slow, and is it possible to go deeper?
689       // If so, then *halve the step*
690       //              ==============
691       if(   (did_len < fraction_done*all_len)
692          && (depth<max_depth) && (!sub_final_section) )
693       {
694 #ifdef G4DEBUG_FIELD         
695         static G4ThreadLocal unsigned int numSplits=0;   // For debugging only 
696         biggest_depth = std::max(depth, biggest_depth);
697         ++numSplits;
698 #endif
699         Second_half = false;
700         ++depth;
701         first_section = false;
702 
703         G4double Sub_len = (all_len-did_len)/(2.);
704         G4FieldTrack midPoint = CurrentA_PointVelocity;
705         G4bool fullAdvance=              
706            integrDriver->AccurateAdvance(midPoint, Sub_len, fiEpsilonStep);
707                          
708         ++fNumAdvanceTrials;
709         if( fullAdvance )  { ++fNumAdvanceFull; }
710 
711         G4double lenAchieved=
712            midPoint.GetCurveLength()-CurrentA_PointVelocity.GetCurveLength();
713         
714         const G4double adequateFraction = (1.0-CLHEP::perThousand);
715         G4bool goodAdvance = (lenAchieved >= adequateFraction * Sub_len);
716         if ( goodAdvance )  { ++fNumAdvanceGood; }
717 
718 #ifdef G4DEBUG_FIELD
719         else  //  !goodAdvance
720         {
721            G4cout << "MLL> AdvanceChordLimited not full at depth=" << depth 
722                   << "  total length achieved = " << lenAchieved << " of "
723                   << Sub_len << "  fraction= ";
724            if (Sub_len != 0.0 ) { G4cout << lenAchieved / Sub_len; }
725            else                 { G4cout << "DivByZero"; }
726            G4cout << "  Good-enough fraction = " << adequateFraction;
727            G4cout << "  Number of call to mll = " << fNumCalls
728                   << "   iteration " << substep_no
729                   << "  inner =  " << substep_no_p << G4endl;
730            G4cout << "  Epsilon of integration = " << fiEpsilonStep << G4endl;
731            G4cout << "  State at start is      = " << CurrentA_PointVelocity
732                   << G4endl
733                   << "        at end (midpoint)= " << midPoint << G4endl;
734            G4cout << "  Particle mass = " << midPoint.GetRestMass() << G4endl;
735 
736            G4EquationOfMotion *equation = integrDriver->GetEquationOfMotion();
737            ReportFieldValue( CurrentA_PointVelocity, "start", equation );
738            ReportFieldValue( midPoint, "midPoint", equation );            
739            G4cout << "  Original Start = "
740                   << CurveStartPointVelocity << G4endl;
741            G4cout << "  Original   End = "
742                   << CurveEndPointVelocity << G4endl;
743            G4cout << "  Original TrialPoint = "
744                   << TrialPoint << G4endl;
745            G4cout << "  (this is STRICT mode) "
746                   << "  num Splits= " << numSplits;           
747            G4cout << G4endl;
748         }
749 #endif
750 
751         *ptrInterMedFT[depth] =  midPoint;
752         CurrentB_PointVelocity = midPoint; 
753 
754         if (fCheckMode)
755         {
756           ++eventCount;
757           endChangeB.push_back(
758             G4LocatorChangeRecord( G4LocatorChangeRecord::kInsertingMidPoint,
759                                    substep_no, eventCount, midPoint) );
760         }        
761         
762         // Adjust 'SubStartPoint' to calculate the 'did_length' in next loop
763         //
764         SubStart_PointVelocity = CurrentA_PointVelocity;
765 
766         // Find new trial intersection point needed at start of the loop
767         //
768         G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
769         G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();   
770 
771         G4ThreeVector PointGe;
772         GetNavigatorFor()->LocateGlobalPointWithinVolume(Point_A);
773         G4bool Intersects_AB = IntersectChord(Point_A, Point_B,
774                                               NewSafety, previousSafety,
775                                               previousSftOrigin, distAB,
776                                               PointGe);
777         if( Intersects_AB )
778         {
779           last_AF_intersection = Intersects_AB;
780           CurrentE_Point = PointGe;
781           fin_section_depth[depth] = true;
782 
783           validIntersectP = true;  // 'E' has been updated.
784 
785           G4bool validNormalAB; 
786           NormalAtEntry  = GetSurfaceNormal( PointGe, validNormalAB ); 
787           validNormalAtE = validNormalAB;
788         }
789         else
790         {
791           // No intersection found for first part of curve
792           // (CurrentA,InterMedPoint[depth]). Go to the second part
793           //
794           Second_half = true;
795 
796           validIntersectP=  false;  // No new 'E' chord intersection found
797         }
798       } // if did_len
799 
800       G4bool unfinished = Second_half;
801       while ( unfinished && (depth>0) )  // Loop checking, 07.10.2016, JA
802       {
803         // Second part of curve (InterMed[depth],Intermed[depth-1])) 
804         // On the depth-1 level normally we are on the 'second_half'
805 
806         //  Find new trial intersection point needed at start of the loop
807         //
808         SubStart_PointVelocity = *ptrInterMedFT[depth];
809         CurrentA_PointVelocity = *ptrInterMedFT[depth];
810         CurrentB_PointVelocity = *ptrInterMedFT[depth-1];
811 
812         if (fCheckMode)
813         {
814           ++eventCount;
815           G4LocatorChangeRecord chngPop_a( G4LocatorChangeRecord::kLevelPop,
816                               substep_no, eventCount, CurrentA_PointVelocity);
817           endChangeA.push_back( chngPop_a );
818           G4LocatorChangeRecord chngPop_b( G4LocatorChangeRecord::kLevelPop,
819                               substep_no, eventCount, CurrentB_PointVelocity);
820           endChangeB.push_back( chngPop_b );
821         }
822         
823         // Ensure that the new endpoints are not further apart in space
824         // than on the curve due to different errors in the integration
825         //
826         G4int  errorEndPt = -1;
827         G4bool recalculatedB= false;
828         if( driverReIntegrates )
829         {
830           // Ensure that the new endpoints are not further apart in space
831           // than on the curve due to different errors in the integration
832           //
833           G4FieldTrack RevisedEndPointFT = CurrentB_PointVelocity;
834           recalculatedB = 
835              CheckAndReEstimateEndpoint( CurrentA_PointVelocity,  
836                                          CurrentB_PointVelocity,
837                                          RevisedEndPointFT,
838                                          errorEndPt );
839           if( recalculatedB )
840           {
841             CurrentB_PointVelocity = RevisedEndPointFT;  // Use it !
842             
843             if ( depth == 1 )
844             {
845                recalculatedEndPoint = true;
846                IntersectedOrRecalculatedFT = RevisedEndPointFT;
847                // So that we can return it, if it is the endpoint!
848             }
849           }
850           else
851           {          
852             if( CurrentB_PointVelocity.GetCurveLength()
853               < CurrentA_PointVelocity.GetCurveLength() )
854             {
855               errorEndPt = 2;
856             }
857           }
858           
859           if (fCheckMode)
860           {
861             ++eventCount;
862             endChangeB.push_back(
863               G4LocatorChangeRecord(G4LocatorChangeRecord::kRecalculatedBagn,
864                                     substep_no, eventCount, RevisedEndPointFT));
865           }
866         }
867         if( errorEndPt > 1 )  // errorEndPt = 1 is milder, just: len(B)=len(A)
868         {
869           std::ostringstream errmsg; 
870           ReportReversedPoints(errmsg, 
871                     CurveStartPointVelocity, CurveEndPointVelocity,
872                     NewSafety, fiEpsilonStep, 
873                     CurrentA_PointVelocity, CurrentB_PointVelocity,
874                     SubStart_PointVelocity, CurrentE_Point,
875                     ApproxIntersecPointV, substep_no, substep_no_p, depth);
876           errmsg << " * Location:  " << MethodName << "- Second-Half" << G4endl;
877           errmsg << " * Recalculated = " << recalculatedB << G4endl; // false
878           G4Exception(MethodName, "GeomNav0003", FatalException, errmsg);
879         }
880         
881         G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
882         G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
883         G4ThreeVector PointGi;     
884         GetNavigatorFor()->LocateGlobalPointWithinVolume(Point_A);
885         G4bool Intersects_AB = IntersectChord(Point_A, Point_B, NewSafety,
886                                               previousSafety,
887                                               previousSftOrigin, distAB,
888                                               PointGi);
889         if( Intersects_AB )
890         {
891           last_AF_intersection = Intersects_AB;
892           CurrentE_Point = PointGi;
893 
894           validIntersectP = true;  // 'E' has been updated.
895           NormalAtEntry  = GetSurfaceNormal( PointGi, validNormalAtE ); 
896         }
897         else
898         {
899           validIntersectP = false;  // No new 'E' chord intersection found
900           if( depth == 1)
901           {
902             there_is_no_intersection = true;
903           }
904         }
905         depth--;
906         fin_section_depth[depth] = true;
907         unfinished = !validIntersectP;
908       }
909 #ifdef G4DEBUG_FIELD
910       if( ! ( validIntersectP || there_is_no_intersection ) )
911       {
912          // What happened ??
913          G4cout << "MLL - WARNING Potential FAILURE: Conditions not met!"
914                 << G4endl
915                 << " Depth = " << depth << G4endl
916                 << " Num Substeps= " << substep_no << G4endl;
917          G4cout << " Found intersection= " << found_approximate_intersection
918                 << G4endl;
919          G4cout << " Progress report: -- " << G4endl;
920          ReportProgress(G4cout, 
921                         CurveStartPointVelocity, CurveEndPointVelocity,
922                         substep_no, CurrentA_PointVelocity,
923                         CurrentB_PointVelocity,
924                         NewSafety, depth ); 
925          G4cout << G4endl; 
926       }
927 #endif      
928     }  // if(!found_aproximate_intersection)
929 
930     assert( validIntersectP || there_is_no_intersection
931                             || found_approximate_intersection);
932 
933   } while ( ( ! found_approximate_intersection )
934             && ( ! there_is_no_intersection )     
935             && ( substep_no <= fMaxSteps) ); // UNTIL found or failed
936 
937   if( substep_no > max_no_seen )
938   {
939     max_no_seen = substep_no; 
940 #ifdef G4DEBUG_FIELD
941     if( max_no_seen > fWarnSteps )
942     {
943       trigger_substepno_print = max_no_seen-20; // Want to see last 20 steps 
944     } 
945 #endif
946   }
947 
948   if( !there_is_no_intersection && !found_approximate_intersection )
949   {
950      if( substep_no >= fMaxSteps)
951      {
952         // Since we cannot go further (yet), we return as far as we have gone
953 
954         recalculatedEndPoint = true;
955         IntersectedOrRecalculatedFT = CurrentA_PointVelocity;
956         found_approximate_intersection = false; 
957      
958         std::ostringstream message;
959         message << G4endl;
960         message << "Convergence is requiring too many substeps: "
961                 << substep_no << "  ( limit = "<<  fMaxSteps << ")"
962                 << G4endl 
963                 << "  Abandoning effort to intersect. " << G4endl << G4endl;
964         message << "    Number of calls to MLL: " <<   fNumCalls;
965         message << "  iteration = " << substep_no <<G4endl << G4endl;
966         
967         message.precision( 12 ); 
968         G4double done_len = CurrentA_PointVelocity.GetCurveLength(); 
969         G4double full_len = CurveEndPointVelocity.GetCurveLength();
970         message << "        Undertaken only length: " << done_len
971                 << " out of " << full_len << " required." << G4endl
972                 << "        Remaining length = " << full_len - done_len; 
973         
974         message << "     Start and end-point of requested Step:" << G4endl;
975         printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
976                      -1.0, NewSafety, 0,     message, -1 );
977         message << "        Start and end-point of current Sub-Step:" << G4endl;
978         printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
979                      -1.0, NewSafety, substep_no-1, message, -1 );
980         printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
981                      -1.0, NewSafety, substep_no, message, -1 );
982         
983         G4Exception(MethodName, "GeomNav0003", JustWarning, message);
984      }
985      else if( substep_no >= fWarnSteps)
986      {  
987         std::ostringstream message;
988         message << "Many substeps while trying to locate intersection."
989                 << G4endl
990                 << "          Undertaken length: "  
991                 << CurrentB_PointVelocity.GetCurveLength()
992                 << " - Needed: "  << substep_no << " substeps." << G4endl
993                 << "          Warning number = " << fWarnSteps
994                 << " and maximum substeps = " << fMaxSteps;
995         G4Exception(MethodName, "GeomNav1002", JustWarning, message);
996      }
997   }
998   
999   return  (!there_is_no_intersection) && found_approximate_intersection;
1000     //  Success or failure
1001 }
1002 
1003 void G4MultiLevelLocator::ReportStatistics()
1004 {
1005    G4cout << " Number of calls = " << fNumCalls << G4endl;
1006    G4cout << " Number of split level ('advances'):  "
1007           << fNumAdvanceTrials << G4endl;
1008    G4cout << " Number of full advances:             "
1009           << fNumAdvanceGood << G4endl;
1010    G4cout << " Number of good advances:             "
1011           << fNumAdvanceFull << G4endl;
1012 }
1013 
1014 void G4MultiLevelLocator::ReportFieldValue( const G4FieldTrack& locationPV,
1015                                             const char* nameLoc,
1016                                             const G4EquationOfMotion* equation )
1017 {
1018    enum { maxNumFieldComp = 24 };
1019 
1020    G4ThreeVector position = locationPV.GetPosition();   
1021    G4double startPoint[4] = { position.x(), position.y(), position.z(),
1022                               locationPV.GetLabTimeOfFlight() };
1023    G4double FieldVec[maxNumFieldComp]; // 24 ;
1024    for (double & i : FieldVec)
1025    {
1026      i = 0.0;
1027    }
1028    equation->GetFieldValue( startPoint, FieldVec);
1029    G4cout << "  B-field value (" << nameLoc << ")=   "
1030           << FieldVec[0] << " " << FieldVec[1] << " " << FieldVec[2];
1031    G4double Emag2= G4ThreeVector( FieldVec[3],
1032                                   FieldVec[4],
1033                                   FieldVec[5] ).mag2(); 
1034    if( Emag2 > 0.0 )
1035    {
1036       G4cout << " Electric = " << FieldVec[3] << " "
1037                                << FieldVec[4] << " "
1038                                << FieldVec[5]<< G4endl;
1039    }
1040    return;
1041 }
1042