Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/magneticfield/src/G4ChordFinder.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 ]

Diff markup

Differences between /geometry/magneticfield/src/G4ChordFinder.cc (Version 11.3.0) and /geometry/magneticfield/src/G4ChordFinder.cc (Version 9.5.p2)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
 26 // G4ChordFinder implementation                << 
 27 //                                                 26 //
 28 // Author: J.Apostolakis - Design and implemen <<  27 // $Id: G4ChordFinder.cc,v 1.53 2009-05-18 14:22:43 gcosmo Exp $
                                                   >>  28 // GEANT4 tag $Name: not supported by cvs2svn $
                                                   >>  29 //
                                                   >>  30 //
                                                   >>  31 // 25.02.97 - John Apostolakis - Design and implementation 
 29 // -------------------------------------------     32 // -------------------------------------------------------------------
 30                                                    33 
 31 #include <iomanip>                                 34 #include <iomanip>
 32                                                    35 
 33 #include "G4ChordFinder.hh"                        36 #include "G4ChordFinder.hh"
 34 #include "G4SystemOfUnits.hh"                  << 
 35 #include "G4MagneticField.hh"                      37 #include "G4MagneticField.hh"
 36 #include "G4Mag_UsualEqRhs.hh"                     38 #include "G4Mag_UsualEqRhs.hh"
 37 #include "G4MagIntegratorDriver.hh"            <<  39 #include "G4ClassicalRK4.hh"
 38 // #include "G4ClassicalRK4.hh"                << 
 39 // #include "G4CashKarpRKF45.hh"               << 
 40 // #include "G4NystromRK4.hh"                  << 
 41 // #include "G4BogackiShampine23.hh"           << 
 42 // #include "G4BogackiShampine45.hh"           << 
 43                                                << 
 44 #include "G4DormandPrince745.hh"               << 
 45                                                << 
 46 // New templated stepper(s) -- avoid virtual c << 
 47 #include "G4TDormandPrince45.hh"               << 
 48                                                << 
 49 // FSAL type driver / steppers -----           << 
 50 #include "G4FSALIntegrationDriver.hh"          << 
 51 #include "G4VFSALIntegrationStepper.hh"        << 
 52 #include "G4RK547FEq1.hh"                      << 
 53 // #include "G4RK547FEq2.hh"                   << 
 54 // #include "G4RK547FEq3.hh"                   << 
 55 // #include "G4FSALBogackiShampine45.hh"       << 
 56 // #include "G4FSALDormandPrince745.hh"        << 
 57                                                << 
 58 // Templated type drivers -----                << 
 59 #include "G4IntegrationDriver.hh"              << 
 60 #include "G4InterpolationDriver.hh"            << 
 61                                                << 
 62 #include "G4HelixHeum.hh"                      << 
 63 #include "G4BFieldIntegrationDriver.hh"        << 
 64                                                << 
 65 #include "G4QSSDriverCreator.hh"               << 
 66                                                    40 
 67 #include "G4CachedMagneticField.hh"            << 
 68                                                    41 
 69 #include <cassert>                             << 
 70 #include <memory>                              << 
 71                                                << 
 72 G4bool G4ChordFinder::gVerboseCtor = false;    << 
 73 // ...........................................     42 // ..........................................................................
 74                                                    43 
 75 G4ChordFinder::G4ChordFinder(G4VIntegrationDri <<  44 G4ChordFinder::G4ChordFinder(G4MagInt_Driver* pIntegrationDriver)
 76   : fDefaultDeltaChord(0.25 * mm), fIntgrDrive <<  45   : fDefaultDeltaChord( 0.25 * mm ),      // Parameters
                                                   >>  46     fDeltaChord( fDefaultDeltaChord ),    //   Internal parameters
                                                   >>  47     fFirstFraction(0.999), fFractionLast(1.00),  fFractionNextEstimate(0.98), 
                                                   >>  48     fMultipleRadius(15.0), 
                                                   >>  49     fStatsVerbose(0),
                                                   >>  50     fDriversStepper(0),                    // Dependent objects 
                                                   >>  51     fAllocatedStepper(false),
                                                   >>  52     fEquation(0),      
                                                   >>  53     fTotalNoTrials_FNC(0), fNoCalls_FNC(0), fmaxTrials_FNC(0)
 77 {                                                  54 {
 78   // Simple constructor -- it does not create      55   // Simple constructor -- it does not create equation
 79   if( gVerboseCtor )                           <<  56   fIntgrDriver= pIntegrationDriver;
 80   {                                            <<  57   fAllocatedStepper= false;
 81     G4cout << "G4ChordFinder: Simple construct <<  58 
 82   }                                            <<  59   fLastStepEstimate_Unconstrained = DBL_MAX;          // Should move q, p to
 83                                                <<  60 
 84   fDeltaChord = fDefaultDeltaChord;       // P <<  61   SetFractions_Last_Next( fFractionLast, fFractionNextEstimate);  
                                                   >>  62     // check the values and set the other parameters
 85 }                                                  63 }
 86                                                    64 
                                                   >>  65 
 87 // ...........................................     66 // ..........................................................................
 88                                                    67 
 89 G4ChordFinder::G4ChordFinder( G4MagneticField*     68 G4ChordFinder::G4ChordFinder( G4MagneticField*        theMagField,
 90                               G4double         <<  69                               G4double                stepMinimum, 
 91                               G4MagIntegratorS <<  70                               G4MagIntegratorStepper* pItsStepper )
 92                               G4int            <<  71   : fDefaultDeltaChord( 0.25 * mm ),     // Constants 
 93   : fDefaultDeltaChord(0.25 * mm)              <<  72     fDeltaChord( fDefaultDeltaChord ),   // Parameters
                                                   >>  73     fFirstFraction(0.999), fFractionLast(1.00),  fFractionNextEstimate(0.98), 
                                                   >>  74     fMultipleRadius(15.0), 
                                                   >>  75     fStatsVerbose(0),
                                                   >>  76     fDriversStepper(0),                  //  Dependent objects
                                                   >>  77     fAllocatedStepper(false),
                                                   >>  78     fEquation(0), 
                                                   >>  79     fTotalNoTrials_FNC(0), fNoCalls_FNC(0), fmaxTrials_FNC(0)  // State - stats
 94 {                                                  80 {
 95   // Construct the Chord Finder                <<  81   //  Construct the Chord Finder
 96   // by creating in inverse order the Driver,  <<  82   //  by creating in inverse order the  Driver, the Stepper and EqRhs ...
 97   constexpr G4int nVar6 = 6;   // Components i << 
 98                                                << 
 99   fDeltaChord = fDefaultDeltaChord;       // P << 
100                                                    83 
101   G4cout << " G4ChordFinder: stepperDriverId:  <<  84   G4Mag_EqRhs *pEquation = new G4Mag_UsualEqRhs(theMagField);
                                                   >>  85   fEquation = pEquation;                            
                                                   >>  86   fLastStepEstimate_Unconstrained = DBL_MAX;          // Should move q, p to
                                                   >>  87                                                      //    G4FieldTrack ??
102                                                    88 
103   G4bool useFSALstepper     = (stepperDriverId <<  89   SetFractions_Last_Next( fFractionLast, fFractionNextEstimate);  
104   G4bool useTemplatedStepper= (stepperDriverId <<  90     // check the values and set the other parameters
105   G4bool useRegularStepper  = (stepperDriverId <<  91 
106   G4bool useBfieldDriver    = (stepperDriverId <<  92   // --->>  Charge    Q = 0 
107   G4bool useG4QSSDriver     = (stepperDriverId <<  93   // --->>  Momentum  P = 1       NOMINAL VALUES !!!!!!!!!!!!!!!!!!
108                                                <<  94 
109   if( stepperDriverId == kQss3DriverType)      <<  95   if( pItsStepper == 0 )
                                                   >>  96   { 
                                                   >>  97      pItsStepper = fDriversStepper = new G4ClassicalRK4(pEquation);
                                                   >>  98      fAllocatedStepper= true;
                                                   >>  99   }
                                                   >> 100   else
110   {                                               101   {
111     stepperDriverId = kQss2DriverType;         << 102      fAllocatedStepper= false; 
112     G4cout << " G4ChordFinder: QSS 3 is curren << 
113   }                                               103   }
                                                   >> 104   fIntgrDriver = new G4MagInt_Driver(stepMinimum, pItsStepper, 
                                                   >> 105                                      pItsStepper->GetNumberOfVariables() );
                                                   >> 106 }
114                                                   107 
115   using EquationType = G4Mag_UsualEqRhs;       << 108 
116                                                << 109 // ......................................................................
117   using TemplatedStepperType =                 << 110 
118          G4TDormandPrince45<EquationType,nVar6 << 111 G4ChordFinder::~G4ChordFinder()
119   const char* TemplatedStepperName =           << 112 {
120       "G4TDormandPrince745 (templated Dormand- << 113   delete   fEquation; // fIntgrDriver->pIntStepper->theEquation_Rhs;
121                                                << 114   if( fAllocatedStepper)
122   using RegularStepperType =                   << 115   { 
123          G4DormandPrince745; // 5th order embe << 116      delete fDriversStepper; 
124          // G4ClassicalRK4;        // The old  << 
125          // G4CashKarpRKF45;       // First em << 
126          // G4BogackiShampine45;   // High eff << 
127          // G4NystromRK4;          // Nystrom  << 
128          // G4RK547FEq1;  // or 2 or 3         << 
129   const char* RegularStepperName =             << 
130       "G4DormandPrince745 (aka DOPRI5): 5th/4t << 
131       // "BogackiShampine 45 (Embedded 5th/4th << 
132       // "Nystrom stepper 4th order";          << 
133                                                << 
134   using NewFsalStepperType = G4DormandPrince74 << 
135                                                << 
136   const char* NewFSALStepperName =             << 
137       "G4RK574FEq1> FSAL 4th/5th order 7-stage << 
138                                                << 
139 #ifdef G4DEBUG_FIELD                           << 
140   static G4bool verboseDebug = true;           << 
141   if( verboseDebug )                           << 
142   {                                            << 
143      G4cout << "G4ChordFinder 2nd Constructor  << 
144      G4cout << " Arguments: " << G4endl        << 
145             << " - min step = " << stepMinimum << 
146             << " - stepper ptr provided : "    << 
147             << ( pItsStepper==nullptr ? " no   << 
148      if( pItsStepper==nullptr )                << 
149         G4cout << " - stepper/driver Id = " << << 
150                << "  useFSAL = " << useFSALste << 
151                << "  , useTemplated = " << use << 
152                << "  , useRegular = " << useRe << 
153                << "  , useFSAL = " << useFSALs << 
154                << G4endl;                      << 
155   }                                               117   }
156 #endif                                         << 118   delete   fIntgrDriver; 
157                                                   119 
158   // useHigherStepper = forceHigherEffiencySte << 120   if( fStatsVerbose ) { PrintStatistics(); }
                                                   >> 121 }
159                                                   122 
160   auto  pEquation = new G4Mag_UsualEqRhs(theMa << 
161   fEquation = pEquation;                       << 
162                                                   123 
163   // G4MagIntegratorStepper* regularStepper =  << 124 // ......................................................................
164   // G4VFSALIntegrationStepper* fsalStepper =  << 
165   // G4MagIntegratorStepper* oldFSALStepper =  << 
166                                                   125 
167   G4bool errorInStepperCreation = false;       << 126 void   
                                                   >> 127 G4ChordFinder::SetFractions_Last_Next( G4double fractLast, G4double fractNext )
                                                   >> 128 { 
                                                   >> 129   // Use -1.0 as request for Default.
                                                   >> 130   if( fractLast == -1.0 )   fractLast = 1.0;   // 0.9;
                                                   >> 131   if( fractNext == -1.0 )   fractNext = 0.98;  // 0.9; 
168                                                   132 
169   std::ostringstream message;  // In case of f << 133   // fFirstFraction  = 0.999; // Orig 0.999 A safe value, range: ~ 0.95 - 0.999
                                                   >> 134   // fMultipleRadius = 15.0;  // For later use, range: ~  2 - 20 
170                                                   135 
171   if( pItsStepper != nullptr )                 << 136   if( fStatsVerbose )
172   {                                            << 137   { 
173      if( gVerboseCtor )                        << 138     G4cout << " ChordFnd> Trying to set fractions: "
174      {                                         << 139            << " first " << fFirstFraction
175        G4cout << " G4ChordFinder: Creating G4I << 140            << " last " <<  fractLast
176               << " stepMinimum = " << stepMini << 141            << " next " <<  fractNext
177               << " numVar= " << pItsStepper->G << 142            << " and multiple " << fMultipleRadius
178      }                                         << 143            << G4endl;
                                                   >> 144   } 
179                                                   145 
180      // Stepper type is not known - so must us << 146   if( (fractLast > 0.0) && (fractLast <=1.0) ) 
181       if(pItsStepper->isQSS())                 << 147   {
182       {                                        << 148     fFractionLast= fractLast;
183          // fIntgrDriver = pItsStepper->build_ << 
184          G4Exception("G4ChordFinder::G4ChordFi << 
185                       "GeomField1001", FatalEx << 
186                       "Cannot provide  QSS ste << 
187       }                                        << 
188       else                                     << 
189       {                                        << 
190          fIntgrDriver = new G4IntegrationDrive << 
191                                   pItsStepper, << 
192          // Stepper type is not known - so mus << 
193          // Non-interpolating driver used by d << 
194          // WAS:  fIntgrDriver = pItsStepper-> << 
195       }                                        << 
196      // -- Older:                              << 
197      // G4cout << " G4ChordFinder: Creating G4 << 
198      // Type is not known - so must use old cl << 
199      // fIntgrDriver = new G4MagInt_Driver( st << 
200      //                                 pItsSt << 
201   }                                               149   }
202   else if ( useTemplatedStepper )              << 150   else
203   {                                               151   {
204      if( gVerboseCtor )                        << 152     G4cerr << "G4ChordFinder::SetFractions_Last_Next: Invalid "
205      {                                         << 153            << " fraction Last = " << fractLast
206         G4cout << " G4ChordFinder: Creating Te << 154            << " must be  0 <  fractionLast <= 1 " << G4endl;
207                << TemplatedStepperName << G4en << 
208      }                                         << 
209      // RegularStepperType* regularStepper = n << 
210      auto templatedStepper = new TemplatedStep << 
211      //                    *** *************** << 
212      //                                        << 
213      // Alternative - for G4NystromRK4:        << 
214      // = new G4NystromRK4(pEquation, 0.1*mm ) << 
215      fRegularStepperOwned = templatedStepper;  << 
216      if( templatedStepper == nullptr )         << 
217      {                                         << 
218         message << "Templated Stepper instanti << 
219         message << "G4ChordFinder: Attempted t << 
220                 << TemplatedStepperName << " t << 
221         errorInStepperCreation = true;         << 
222      }                                         << 
223      else                                      << 
224      {                                         << 
225         fIntgrDriver = new G4IntegrationDriver << 
226            stepMinimum, templatedStepper, nVar << 
227         if( gVerboseCtor )                     << 
228         {                                      << 
229            G4cout << " G4ChordFinder: Using G4 << 
230         }                                      << 
231      }                                         << 
232                                                << 
233   }                                               155   }
234   else if ( useRegularStepper   )  // Plain st << 156   if( (fractNext > 0.0) && (fractNext <1.0) )
235   {                                               157   {
236      auto regularStepper = new RegularStepperT << 158     fFractionNextEstimate = fractNext;
237      //                    *** *************** << 159   }
238      fRegularStepperOwned = regularStepper;    << 160   else
                                                   >> 161   {
                                                   >> 162     G4cerr << "G4ChordFinder:: SetFractions_Last_Next: Invalid "
                                                   >> 163            << " fraction Next = " << fractNext
                                                   >> 164            << " must be  0 <  fractionNext < 1 " << G4endl;
                                                   >> 165   }
                                                   >> 166 }
239                                                   167 
240      if( gVerboseCtor )                        << 168 
241      {                                         << 169 // ......................................................................
242         G4cout << " G4ChordFinder: Creating Dr << 170 
                                                   >> 171 G4double 
                                                   >> 172 G4ChordFinder::AdvanceChordLimited( G4FieldTrack& yCurrent,
                                                   >> 173                                     G4double      stepMax,
                                                   >> 174                                     G4double      epsStep,
                                                   >> 175                                     const G4ThreeVector latestSafetyOrigin,
                                                   >> 176                                     G4double       latestSafetyRadius )
                                                   >> 177 {
                                                   >> 178   G4double stepPossible;
                                                   >> 179   G4double dyErr;
                                                   >> 180   G4FieldTrack yEnd( yCurrent);
                                                   >> 181   G4double  startCurveLen= yCurrent.GetCurveLength();
                                                   >> 182   G4double nextStep;
                                                   >> 183   //            *************
                                                   >> 184   stepPossible= FindNextChord(yCurrent, stepMax, yEnd, dyErr, epsStep,
                                                   >> 185                               &nextStep, latestSafetyOrigin, latestSafetyRadius
                                                   >> 186                              );
                                                   >> 187   //            *************
                                                   >> 188 
                                                   >> 189   G4bool good_advance;
                                                   >> 190 
                                                   >> 191   if ( dyErr < epsStep * stepPossible )
                                                   >> 192   {
                                                   >> 193      // Accept this accuracy.
                                                   >> 194 
                                                   >> 195      yCurrent = yEnd;
                                                   >> 196      good_advance = true; 
                                                   >> 197   }
                                                   >> 198   else
                                                   >> 199   {  
                                                   >> 200      // Advance more accurately to "end of chord"
                                                   >> 201      //                           ***************
                                                   >> 202      good_advance = fIntgrDriver->AccurateAdvance(yCurrent, stepPossible,
                                                   >> 203                                                   epsStep, nextStep);
                                                   >> 204      if ( ! good_advance )
                                                   >> 205      { 
                                                   >> 206        // In this case the driver could not do the full distance
                                                   >> 207        stepPossible= yCurrent.GetCurveLength()-startCurveLen;
243      }                                            208      }
                                                   >> 209   }
                                                   >> 210   return stepPossible;
                                                   >> 211 }
                                                   >> 212 
                                                   >> 213 
                                                   >> 214 // ............................................................................
                                                   >> 215 
                                                   >> 216 G4double
                                                   >> 217 G4ChordFinder::FindNextChord( const  G4FieldTrack& yStart,
                                                   >> 218                                      G4double     stepMax,
                                                   >> 219                                      G4FieldTrack&   yEnd, // Endpoint
                                                   >> 220                                      G4double&   dyErrPos, // Error of endpoint
                                                   >> 221                                      G4double    epsStep,
                                                   >> 222                                      G4double*  pStepForAccuracy, 
                                                   >> 223                               const  G4ThreeVector, //  latestSafetyOrigin,
                                                   >> 224                                      G4double       //  latestSafetyRadius 
                                                   >> 225                                         )
                                                   >> 226 {
                                                   >> 227   // Returns Length of Step taken
                                                   >> 228 
                                                   >> 229   G4FieldTrack yCurrent=  yStart;  
                                                   >> 230   G4double    stepTrial, stepForAccuracy;
                                                   >> 231   G4double    dydx[G4FieldTrack::ncompSVEC]; 
                                                   >> 232 
                                                   >> 233   //  1.)  Try to "leap" to end of interval
                                                   >> 234   //  2.)  Evaluate if resulting chord gives d_chord that is good enough.
                                                   >> 235   // 2a.)  If d_chord is not good enough, find one that is.
                                                   >> 236   
                                                   >> 237   G4bool    validEndPoint= false;
                                                   >> 238   G4double  dChordStep, lastStepLength; //  stepOfLastGoodChord;
                                                   >> 239 
                                                   >> 240   fIntgrDriver-> GetDerivatives( yCurrent, dydx );
                                                   >> 241 
                                                   >> 242   G4int     noTrials=0;
                                                   >> 243   const G4double safetyFactor= fFirstFraction; //  0.975 or 0.99 ? was 0.999
                                                   >> 244 
                                                   >> 245   stepTrial = std::min( stepMax, safetyFactor*fLastStepEstimate_Unconstrained );
                                                   >> 246 
                                                   >> 247   G4double newStepEst_Uncons= 0.0; 
                                                   >> 248   do
                                                   >> 249   { 
                                                   >> 250      G4double stepForChord;  
                                                   >> 251      yCurrent = yStart;    // Always start from initial point
                                                   >> 252     
                                                   >> 253      //            ************
                                                   >> 254      fIntgrDriver->QuickAdvance( yCurrent, dydx, stepTrial, 
                                                   >> 255                                  dChordStep, dyErrPos);
                                                   >> 256      //            ************
244                                                   257      
245      if( regularStepper == nullptr )           << 258      //  We check whether the criterion is met here.
246      {                                         << 259      validEndPoint = AcceptableMissDist(dChordStep);
247         message << "Regular Stepper instantiat << 260 
248         message << "G4ChordFinder: Attempted t << 261      lastStepLength = stepTrial; 
249                 << RegularStepperName << " typ << 262 
250         errorInStepperCreation = true;         << 263      // This method estimates to step size for a good chord.
251      }                                         << 264      stepForChord = NewStep(stepTrial, dChordStep, newStepEst_Uncons );
252      else                                      << 265 
                                                   >> 266      if( ! validEndPoint )
253      {                                            267      {
254         auto dp5= dynamic_cast<G4DormandPrince << 268         if( stepTrial<=0.0 )
255         if( dp5 != nullptr )                   << 269         {
                                                   >> 270           stepTrial = stepForChord;
                                                   >> 271         }
                                                   >> 272         else if (stepForChord <= stepTrial)
256         {                                         273         {
257            fIntgrDriver = new G4InterpolationD << 274           // Reduce by a fraction, possibly up to 20% 
258                                   stepMinimum, << 275           stepTrial = std::min( stepForChord, fFractionLast * stepTrial);
259            if( gVerboseCtor )                  << 
260            {                                   << 
261               G4cout << " Using InterpolationD << 
262            }                                   << 
263         }                                         276         }
264         else                                      277         else
265         {                                         278         {
266            fIntgrDriver = new G4IntegrationDri << 279           stepTrial *= 0.1;
267                                   stepMinimum, << 
268            if( gVerboseCtor )                  << 
269            {                                   << 
270               G4cout << " Using IntegrationDri << 
271            }                                   << 
272         }                                         280         }
273      }                                            281      }
                                                   >> 282      noTrials++; 
274   }                                               283   }
275   else if ( useBfieldDriver )                  << 284   while( ! validEndPoint );   // End of do-while  RKD 
                                                   >> 285 
                                                   >> 286   if( newStepEst_Uncons > 0.0  )
276   {                                               287   {
277      auto regularStepper = new G4DormandPrince << 288      fLastStepEstimate_Unconstrained= newStepEst_Uncons;
278      //                    *** *************** << 
279      //                                        << 
280      fRegularStepperOwned = regularStepper;    << 
281                                                << 
282      {                                         << 
283         using SmallStepDriver = G4Interpolatio << 
284         using LargeStepDriver = G4IntegrationD << 
285                                                << 
286         fLongStepper = std::make_unique<G4Heli << 
287                                                << 
288         fIntgrDriver = new G4BFieldIntegration << 
289           std::make_unique<SmallStepDriver>(st << 
290               regularStepper, regularStepper-> << 
291           std::make_unique<LargeStepDriver>(st << 
292               fLongStepper.get(), regularStepp << 
293                                                << 
294         if( fIntgrDriver == nullptr)           << 
295         {                                      << 
296            message << "Using G4BFieldIntegrati << 
297                    << RegularStepperName << "  << 
298            message << "Driver instantiation FA << 
299            G4Exception("G4ChordFinder::G4Chord << 
300                        "GeomField1001", JustWa << 
301         }                                      << 
302      }                                         << 
303   }                                               289   }
304   else if( useG4QSSDriver )                    << 290 
305   {                                            << 291   AccumulateStatistics( noTrials );
306      if( stepperDriverId == kQss2DriverType )  << 292 
                                                   >> 293   if( pStepForAccuracy )
                                                   >> 294   { 
                                                   >> 295      // Calculate the step size required for accuracy, if it is needed
                                                   >> 296      //
                                                   >> 297      G4double dyErr_relative = dyErrPos/(epsStep*lastStepLength);
                                                   >> 298      if( dyErr_relative > 1.0 )
307      {                                            299      {
308        auto qssStepper2 = G4QSSDriverCreator:: << 300         stepForAccuracy = fIntgrDriver->ComputeNewStepSize( dyErr_relative,
309        if( gVerboseCtor )                      << 301                                                             lastStepLength );
310        {                                       << 
311          G4cout << "-- Created QSS-2 stepper"  << 
312        }                                       << 
313        fIntgrDriver = G4QSSDriverCreator::Crea << 
314      }                                            302      }
315      else                                         303      else
316      {                                            304      {
317        auto qssStepper3 = G4QSSDriverCreator:: << 305         stepForAccuracy = 0.0;   // Convention to show step was ok 
318        if( gVerboseCtor )                      << 
319        {                                       << 
320          G4cout << "-- Created QSS-3 stepper"  << 
321        }                                       << 
322        fIntgrDriver = G4QSSDriverCreator::Crea << 
323      }                                         << 
324      if( gVerboseCtor )                        << 
325      {                                         << 
326        G4cout << "-- G4ChordFinder: Using QSS  << 
327      }                                            306      }
                                                   >> 307      *pStepForAccuracy = stepForAccuracy;
                                                   >> 308   }
                                                   >> 309 
                                                   >> 310 #ifdef  TEST_CHORD_PRINT
                                                   >> 311   static int dbg=0;
                                                   >> 312   if( dbg )
                                                   >> 313   {
                                                   >> 314     G4cout << "ChordF/FindNextChord:  NoTrials= " << noTrials 
                                                   >> 315            << " StepForGoodChord=" << std::setw(10) << stepTrial << G4endl;
                                                   >> 316   }
                                                   >> 317 #endif
                                                   >> 318   yEnd=  yCurrent;  
                                                   >> 319   return stepTrial; 
                                                   >> 320 }
                                                   >> 321 
                                                   >> 322 
                                                   >> 323 // ...........................................................................
                                                   >> 324 
                                                   >> 325 G4double G4ChordFinder::NewStep(G4double  stepTrialOld, 
                                                   >> 326                                 G4double  dChordStep, // Curr. dchord achieved
                                                   >> 327                                 G4double& stepEstimate_Unconstrained )  
                                                   >> 328 {
                                                   >> 329   // Is called to estimate the next step size, even for successful steps,
                                                   >> 330   // in order to predict an accurate 'chord-sensitive' first step
                                                   >> 331   // which is likely to assist in more performant 'stepping'.
                                                   >> 332 
                                                   >> 333   G4double stepTrial;
                                                   >> 334 
                                                   >> 335 #if 1
                                                   >> 336 
                                                   >> 337   if (dChordStep > 0.0)
                                                   >> 338   {
                                                   >> 339     stepEstimate_Unconstrained =
                                                   >> 340                  stepTrialOld*std::sqrt( fDeltaChord / dChordStep );
                                                   >> 341     stepTrial =  fFractionNextEstimate * stepEstimate_Unconstrained;
328   }                                               342   }
329   else                                            343   else
330   {                                               344   {
331      auto fsalStepper=  new NewFsalStepperType << 345     // Should not update the Unconstrained Step estimate: incorrect!
332      //                 *** ****************** << 346     stepTrial =  stepTrialOld * 2.; 
333      fNewFSALStepperOwned = fsalStepper;       << 347   }
334                                                   348 
335      if( fsalStepper == nullptr )              << 349   if( stepTrial <= 0.001 * stepTrialOld)
                                                   >> 350   {
                                                   >> 351      if ( dChordStep > 1000.0 * fDeltaChord )
336      {                                            352      {
337         message << "Stepper instantiation FAIL << 353         stepTrial= stepTrialOld * 0.03;   
338         message << "Attempted to instantiate " << 
339                 << NewFSALStepperName << " typ << 
340         G4Exception("G4ChordFinder::G4ChordFin << 
341                     "GeomField1001", JustWarni << 
342         errorInStepperCreation = true;         << 
343      }                                            354      }
344      else                                         355      else
345      {                                            356      {
346         fIntgrDriver = new                     << 357         if ( dChordStep > 100. * fDeltaChord )
347            G4FSALIntegrationDriver<NewFsalStep << 358         {
348                                           fsal << 359           stepTrial= stepTrialOld * 0.1;   
349            //  ====  Create the driver which k << 360         }
350                                                << 361         else   // Try halving the length until dChordStep OK
351         if( fIntgrDriver == nullptr )          << 
352         {                                         362         {
353            message << "Using G4FSALIntegration << 363           stepTrial= stepTrialOld * 0.5;   
354                    << NewFSALStepperName << G4 << 
355            message << "Integration Driver inst << 
356            G4Exception("G4ChordFinder::G4Chord << 
357                        "GeomField1001", JustWa << 
358         }                                         364         }
359      }                                            365      }
360   }                                               366   }
                                                   >> 367   else if (stepTrial > 1000.0 * stepTrialOld)
                                                   >> 368   {
                                                   >> 369      stepTrial= 1000.0 * stepTrialOld;
                                                   >> 370   }
361                                                   371 
362   // -- Main work is now done                  << 372   if( stepTrial == 0.0 )
363                                                << 
364   //    Now check that no error occured, and r << 
365                                                << 
366   // To test failure to create driver          << 
367   // delete fIntgrDriver;                      << 
368   // fIntgrDriver = nullptr;                   << 
369                                                << 
370   // Detect and report Error conditions        << 
371   //                                           << 
372   if( errorInStepperCreation || (fIntgrDriver  << 
373   {                                               373   {
374      std::ostringstream errmsg;                << 374      stepTrial= 0.000001;
375                                                << 375   }
376      if( errorInStepperCreation )              << 376 
                                                   >> 377 #else
                                                   >> 378 
                                                   >> 379   if ( dChordStep > 1000. * fDeltaChord )
                                                   >> 380   {
                                                   >> 381         stepTrial= stepTrialOld * 0.03;   
                                                   >> 382   }
                                                   >> 383   else
                                                   >> 384   {
                                                   >> 385      if ( dChordStep > 100. * fDeltaChord )
377      {                                            386      {
378         errmsg  << "ERROR> Failure to create S << 387         stepTrial= stepTrialOld * 0.1;   
379                 << "       ------------------- << 
380      }                                            388      }
381      if (fIntgrDriver == nullptr )             << 389      else  // Keep halving the length until dChordStep OK
382      {                                            390      {
383         errmsg  << "ERROR> Failure to create I << 391         stepTrial= stepTrialOld * 0.5;   
384                 << G4endl                      << 
385                 << "       ------------------- << 
386                 << G4endl;                     << 
387      }                                            392      }
388      const std::string BoolName[2]= { "False", << 
389      errmsg << "  Configuration:  (constructor << 
390             << "    provided Stepper = " << pI << 
391             << " stepper/driver Id = " << step << 
392             << "   useTemplated = " << BoolNam << 
393             << "   useRegular = " << BoolName[ << 
394             << "   useFSAL = " << BoolName[use << 
395             << "   using combo BField Driver = << 
396                    BoolName[ ! (useFSALstepper << 
397                                || useRegularSt << 
398             << G4endl;                         << 
399      errmsg << message.str();                  << 
400      errmsg << "Aborting.";                    << 
401      G4Exception("G4ChordFinder::G4ChordFinder << 
402                  "GeomField0003", FatalExcepti << 
403   }                                               393   }
404                                                   394 
405   assert(    ( pItsStepper != nullptr )        << 395 #endif 
406           || ( fRegularStepperOwned != nullptr << 
407           || ( fNewFSALStepperOwned != nullptr << 
408           || useG4QSSDriver                    << 
409      );                                        << 
410   assert( fIntgrDriver != nullptr );           << 
411 }                                              << 
412                                                   396 
413 // ........................................... << 397   // A more sophisticated chord-finder could figure out a better
                                                   >> 398   // stepTrial, from dChordStep and the required d_geometry
                                                   >> 399   //   e.g.
                                                   >> 400   //      Calculate R, r_helix (eg at orig point)
                                                   >> 401   //      if( stepTrial < 2 pi  R )
                                                   >> 402   //          stepTrial = R arc_cos( 1 - fDeltaChord / r_helix )
                                                   >> 403   //      else    
                                                   >> 404   //          ??
414                                                   405 
415 G4ChordFinder::~G4ChordFinder()                << 406   return stepTrial;
416 {                                              << 
417   delete fEquation;                            << 
418   delete fRegularStepperOwned;                 << 
419   delete fNewFSALStepperOwned;                 << 
420   delete fCachedField;                         << 
421   delete fIntgrDriver;                         << 
422 }                                                 407 }
423                                                   408 
                                                   >> 409 
424 // ...........................................    410 // ...........................................................................
425                                                   411 
426 G4FieldTrack                                      412 G4FieldTrack
427 G4ChordFinder::ApproxCurvePointS( const G4Fiel    413 G4ChordFinder::ApproxCurvePointS( const G4FieldTrack&  CurveA_PointVelocity, 
428                                   const G4Fiel    414                                   const G4FieldTrack&  CurveB_PointVelocity, 
429                                   const G4Fiel    415                                   const G4FieldTrack&  ApproxCurveV,
430                                   const G4Thre    416                                   const G4ThreeVector& CurrentE_Point,
431                                   const G4Thre    417                                   const G4ThreeVector& CurrentF_Point,
432                                   const G4Thre    418                                   const G4ThreeVector& PointG,
433                                         G4bool << 419                                        G4bool first, G4double eps_step)
434 {                                                 420 {
435   // ApproxCurvePointS is 2nd implementation o    421   // ApproxCurvePointS is 2nd implementation of ApproxCurvePoint.
436   // Use Brent Algorithm (or InvParabolic) whe    422   // Use Brent Algorithm (or InvParabolic) when possible.
437   // Given a starting curve point A (CurveA_Po    423   // Given a starting curve point A (CurveA_PointVelocity), curve point B
438   // (CurveB_PointVelocity), a point E which i    424   // (CurveB_PointVelocity), a point E which is (generally) not on the curve
439   // and  a point F which is on the curve (fir    425   // and  a point F which is on the curve (first approximation), find new
440   // point S on the curve closer to point E.      426   // point S on the curve closer to point E. 
441   // While advancing towards S utilise 'eps_st    427   // While advancing towards S utilise 'eps_step' as a measure of the
442   // relative accuracy of each Step.              428   // relative accuracy of each Step.
443                                                   429 
444   G4FieldTrack EndPoint(CurveA_PointVelocity);    430   G4FieldTrack EndPoint(CurveA_PointVelocity);
445   if(!first) { EndPoint = ApproxCurveV; }      << 431   if(!first){EndPoint= ApproxCurveV;}
446                                                   432 
447   G4ThreeVector Point_A,Point_B;                  433   G4ThreeVector Point_A,Point_B;
448   Point_A=CurveA_PointVelocity.GetPosition();     434   Point_A=CurveA_PointVelocity.GetPosition();
449   Point_B=CurveB_PointVelocity.GetPosition();     435   Point_B=CurveB_PointVelocity.GetPosition();
450                                                   436 
451   G4double xa,xb,xc,ya,yb,yc;                     437   G4double xa,xb,xc,ya,yb,yc;
452                                                   438  
453   // InverseParabolic. AF Intersects (First Pa    439   // InverseParabolic. AF Intersects (First Part of Curve) 
454                                                   440 
455   if(first)                                       441   if(first)
456   {                                               442   {
457     xa=0.;                                        443     xa=0.;
458     ya=(PointG-Point_A).mag();                    444     ya=(PointG-Point_A).mag();
459     xb=(Point_A-CurrentF_Point).mag();            445     xb=(Point_A-CurrentF_Point).mag();
460     yb=-(PointG-CurrentF_Point).mag();            446     yb=-(PointG-CurrentF_Point).mag();
461     xc=(Point_A-Point_B).mag();                   447     xc=(Point_A-Point_B).mag();
462     yc=-(CurrentE_Point-Point_B).mag();           448     yc=-(CurrentE_Point-Point_B).mag();
463   }                                               449   }    
464   else                                            450   else
465   {                                               451   {
466     xa=0.;                                        452     xa=0.;
467     ya=(Point_A-CurrentE_Point).mag();            453     ya=(Point_A-CurrentE_Point).mag();
468     xb=(Point_A-CurrentF_Point).mag();            454     xb=(Point_A-CurrentF_Point).mag();
469     yb=(PointG-CurrentF_Point).mag();             455     yb=(PointG-CurrentF_Point).mag();
470     xc=(Point_A-Point_B).mag();                   456     xc=(Point_A-Point_B).mag();
471     yc=-(Point_B-PointG).mag();                   457     yc=-(Point_B-PointG).mag();
472     if(xb==0.)                                    458     if(xb==0.)
473     {                                             459     {
474       EndPoint = ApproxCurvePointV(CurveA_Poin << 460       EndPoint=
475                                    CurrentE_Po << 461       ApproxCurvePointV(CurveA_PointVelocity, CurveB_PointVelocity,
                                                   >> 462                         CurrentE_Point, eps_step);
476       return EndPoint;                            463       return EndPoint;
477     }                                             464     }
478   }                                               465   }
479                                                   466 
480   const G4double tolerance = 1.e-12;           << 467   const G4double tolerance= 1.e-12;
481   if(std::abs(ya)<=tolerance||std::abs(yc)<=to    468   if(std::abs(ya)<=tolerance||std::abs(yc)<=tolerance)
482   {                                               469   {
483     ; // What to do for the moment: return the    470     ; // What to do for the moment: return the same point as at start
484       // then PropagatorInField will take care    471       // then PropagatorInField will take care
485   }                                               472   }
486   else                                            473   else
487   {                                               474   {
488     G4double test_step = InvParabolic(xa,ya,xb    475     G4double test_step = InvParabolic(xa,ya,xb,yb,xc,yc);
489     G4double curve;                               476     G4double curve;
490     if(first)                                     477     if(first)
491     {                                             478     {
492       curve=std::abs(EndPoint.GetCurveLength()    479       curve=std::abs(EndPoint.GetCurveLength()
493                     -ApproxCurveV.GetCurveLeng    480                     -ApproxCurveV.GetCurveLength());
494     }                                             481     }
495     else                                          482     else
496     {                                             483     {
497       test_step = test_step - xb;              << 484       test_step=(test_step-xb);
498       curve=std::abs(EndPoint.GetCurveLength()    485       curve=std::abs(EndPoint.GetCurveLength()
499                     -CurveB_PointVelocity.GetC    486                     -CurveB_PointVelocity.GetCurveLength());
500       xb = (CurrentF_Point-Point_B).mag();     << 487       xb=(CurrentF_Point-Point_B).mag();
501     }                                             488     }
502                                                   489       
503     if(test_step<=0)    { test_step=0.1*xb; }     490     if(test_step<=0)    { test_step=0.1*xb; }
504     if(test_step>=xb)   { test_step=0.5*xb; }     491     if(test_step>=xb)   { test_step=0.5*xb; }
505     if(test_step>=curve){ test_step=0.5*curve;    492     if(test_step>=curve){ test_step=0.5*curve; } 
506                                                   493 
507     if(curve*(1.+eps_step)<xb) // Similar to R    494     if(curve*(1.+eps_step)<xb) // Similar to ReEstimate Step from
508     {                          // G4VIntersect    495     {                          // G4VIntersectionLocator
509       test_step=0.5*curve;                        496       test_step=0.5*curve;
510     }                                             497     }
511                                                   498 
512     fIntgrDriver->AccurateAdvance(EndPoint,tes    499     fIntgrDriver->AccurateAdvance(EndPoint,test_step, eps_step);
513                                                   500       
514 #ifdef G4DEBUG_FIELD                              501 #ifdef G4DEBUG_FIELD
515     // Printing Brent and Linear Approximation    502     // Printing Brent and Linear Approximation
516     //                                            503     //
517     G4cout << "G4ChordFinder::ApproxCurvePoint    504     G4cout << "G4ChordFinder::ApproxCurvePointS() - test-step ShF = "
518            << test_step << "  EndPoint = " <<     505            << test_step << "  EndPoint = " << EndPoint << G4endl;
519                                                   506 
520     //  Test Track                                507     //  Test Track
521     //                                            508     //
522     G4FieldTrack TestTrack( CurveA_PointVeloci    509     G4FieldTrack TestTrack( CurveA_PointVelocity);
523     TestTrack = ApproxCurvePointV( CurveA_Poin    510     TestTrack = ApproxCurvePointV( CurveA_PointVelocity, 
524                                    CurveB_Poin    511                                    CurveB_PointVelocity, 
525                                    CurrentE_Po    512                                    CurrentE_Point, eps_step );
526     G4cout.precision(14);                         513     G4cout.precision(14);
527     G4cout << "G4ChordFinder::BrentApprox = "     514     G4cout << "G4ChordFinder::BrentApprox = " << EndPoint  << G4endl;
528     G4cout << "G4ChordFinder::LinearApprox= "     515     G4cout << "G4ChordFinder::LinearApprox= " << TestTrack << G4endl; 
529 #endif                                            516 #endif
530   }                                               517   }
531   return EndPoint;                                518   return EndPoint;
532 }                                                 519 }
533                                                   520 
534                                                   521 
535 // ...........................................    522 // ...........................................................................
536                                                   523 
537 G4FieldTrack G4ChordFinder::                      524 G4FieldTrack G4ChordFinder::
538 ApproxCurvePointV( const G4FieldTrack& CurveA_    525 ApproxCurvePointV( const G4FieldTrack& CurveA_PointVelocity, 
539                    const G4FieldTrack& CurveB_    526                    const G4FieldTrack& CurveB_PointVelocity, 
540                    const G4ThreeVector& Curren    527                    const G4ThreeVector& CurrentE_Point,
541                          G4double eps_step)       528                          G4double eps_step)
542 {                                                 529 {
543   // If r=|AE|/|AB|, and s=true path lenght (A    530   // If r=|AE|/|AB|, and s=true path lenght (AB)
544   // return the point that is r*s along the cu    531   // return the point that is r*s along the curve!
545                                                   532  
546   G4FieldTrack   Current_PointVelocity = Curve    533   G4FieldTrack   Current_PointVelocity = CurveA_PointVelocity; 
547                                                   534 
548   G4ThreeVector  CurveA_Point= CurveA_PointVel    535   G4ThreeVector  CurveA_Point= CurveA_PointVelocity.GetPosition();
549   G4ThreeVector  CurveB_Point= CurveB_PointVel    536   G4ThreeVector  CurveB_Point= CurveB_PointVelocity.GetPosition();
550                                                   537 
551   G4ThreeVector  ChordAB_Vector= CurveB_Point     538   G4ThreeVector  ChordAB_Vector= CurveB_Point   - CurveA_Point;
552   G4ThreeVector  ChordAE_Vector= CurrentE_Poin    539   G4ThreeVector  ChordAE_Vector= CurrentE_Point - CurveA_Point;
553                                                   540 
554   G4double       ABdist= ChordAB_Vector.mag();    541   G4double       ABdist= ChordAB_Vector.mag();
555   G4double  curve_length;  //  A curve length     542   G4double  curve_length;  //  A curve length  of AB
556   G4double  AE_fraction;                          543   G4double  AE_fraction; 
557                                                   544   
558   curve_length= CurveB_PointVelocity.GetCurveL    545   curve_length= CurveB_PointVelocity.GetCurveLength()
559               - CurveA_PointVelocity.GetCurveL    546               - CurveA_PointVelocity.GetCurveLength();  
560                                                   547  
561   G4double integrationInaccuracyLimit= std::ma << 548   G4double  integrationInaccuracyLimit= std::max( perMillion, 0.5*eps_step ); 
562   if( curve_length < ABdist * (1. - integratio    549   if( curve_length < ABdist * (1. - integrationInaccuracyLimit) )
563   {                                               550   { 
564 #ifdef G4DEBUG_FIELD                              551 #ifdef G4DEBUG_FIELD
565     G4cerr << " Warning in G4ChordFinder::Appr    552     G4cerr << " Warning in G4ChordFinder::ApproxCurvePoint: "
566            << G4endl                              553            << G4endl
567            << " The two points are further apa    554            << " The two points are further apart than the curve length "
568            << G4endl                              555            << G4endl
569            << " Dist = "         << ABdist        556            << " Dist = "         << ABdist
570            << " curve length = " << curve_leng    557            << " curve length = " << curve_length 
571            << " relativeDiff = " << (curve_len    558            << " relativeDiff = " << (curve_length-ABdist)/ABdist 
572            << G4endl;                             559            << G4endl;
573     if( curve_length < ABdist * (1. - 10*eps_s    560     if( curve_length < ABdist * (1. - 10*eps_step) )
574     {                                             561     {
575       std::ostringstream message;                 562       std::ostringstream message;
576       message << "Unphysical curve length." <<    563       message << "Unphysical curve length." << G4endl
577               << "The size of the above differ    564               << "The size of the above difference exceeds allowed limits."
578               << G4endl                           565               << G4endl
579               << "Aborting.";                     566               << "Aborting.";
580       G4Exception("G4ChordFinder::ApproxCurveP    567       G4Exception("G4ChordFinder::ApproxCurvePointV()", "GeomField0003",
581                   FatalException, message);       568                   FatalException, message);
582     }                                             569     }
583 #endif                                            570 #endif
584     // Take default corrective action: adjust     571     // Take default corrective action: adjust the maximum curve length. 
585     // NOTE: this case only happens for relati    572     // NOTE: this case only happens for relatively straight paths.
586     // curve_length = ABdist;                     573     // curve_length = ABdist; 
587   }                                               574   }
588                                                   575 
589   G4double new_st_length;                      << 576   G4double  new_st_length; 
590                                                   577 
591   if ( ABdist > 0.0 )                             578   if ( ABdist > 0.0 )
592   {                                               579   {
593      AE_fraction = ChordAE_Vector.mag() / ABdi    580      AE_fraction = ChordAE_Vector.mag() / ABdist;
594   }                                               581   }
595   else                                            582   else
596   {                                               583   {
597      AE_fraction = 0.5;                           584      AE_fraction = 0.5;                         // Guess .. ?; 
598 #ifdef G4DEBUG_FIELD                              585 #ifdef G4DEBUG_FIELD
599      G4cout << "Warning in G4ChordFinder::Appr    586      G4cout << "Warning in G4ChordFinder::ApproxCurvePointV():"
600             << " A and B are the same point!"     587             << " A and B are the same point!" << G4endl
601             << " Chord AB length = " << ChordA    588             << " Chord AB length = " << ChordAE_Vector.mag() << G4endl
602             << G4endl;                            589             << G4endl;
603 #endif                                            590 #endif
604   }                                               591   }
605                                                   592   
606   if( (AE_fraction> 1.0 + perMillion) || (AE_f    593   if( (AE_fraction> 1.0 + perMillion) || (AE_fraction< 0.) )
607   {                                               594   {
608 #ifdef G4DEBUG_FIELD                              595 #ifdef G4DEBUG_FIELD
609     G4cerr << " G4ChordFinder::ApproxCurvePoin    596     G4cerr << " G4ChordFinder::ApproxCurvePointV() - Warning:"
610            << " Anomalous condition:AE > AB or    597            << " Anomalous condition:AE > AB or AE/AB <= 0 " << G4endl
611            << "   AE_fraction = " <<  AE_fract    598            << "   AE_fraction = " <<  AE_fraction << G4endl
612            << "   Chord AE length = " << Chord    599            << "   Chord AE length = " << ChordAE_Vector.mag() << G4endl
613            << "   Chord AB length = " << ABdis    600            << "   Chord AB length = " << ABdist << G4endl << G4endl;
614     G4cerr << " OK if this condition occurs af    601     G4cerr << " OK if this condition occurs after a recalculation of 'B'"
615            << G4endl << " Otherwise it is an e    602            << G4endl << " Otherwise it is an error. " << G4endl ; 
616 #endif                                            603 #endif
617      // This course can now result if B has be    604      // This course can now result if B has been re-evaluated, 
618      // without E being recomputed (1 July 99)    605      // without E being recomputed (1 July 99).
619      // In this case this is not a "real error    606      // In this case this is not a "real error" - but it is undesired
620      // and we cope with it by a default corre    607      // and we cope with it by a default corrective action ...
621      //                                           608      //
622      AE_fraction = 0.5;                           609      AE_fraction = 0.5;                         // Default value
623   }                                               610   }
624                                                   611 
625   new_st_length = AE_fraction * curve_length;  << 612   new_st_length= AE_fraction * curve_length; 
626                                                   613 
627   if ( AE_fraction > 0.0 )                        614   if ( AE_fraction > 0.0 )
628   {                                               615   { 
629      fIntgrDriver->AccurateAdvance(Current_Poi    616      fIntgrDriver->AccurateAdvance(Current_PointVelocity, 
630                                    new_st_leng    617                                    new_st_length, eps_step );
631      //                                           618      //
632      // In this case it does not matter if it     619      // In this case it does not matter if it cannot advance the full distance
633   }                                               620   }
634                                                   621 
635   // If there was a memory of the step_length     622   // If there was a memory of the step_length actually required at the start 
636   // of the integration Step, this could be re    623   // of the integration Step, this could be re-used ...
637                                                   624 
638   G4cout.precision(14);                           625   G4cout.precision(14);
639                                                   626 
640   return Current_PointVelocity;                   627   return Current_PointVelocity;
641 }                                                 628 }
642                                                   629 
643 // ........................................... << 
644                                                   630 
645 std::ostream& operator<<( std::ostream& os, co << 631 // ......................................................................
                                                   >> 632 
                                                   >> 633 void
                                                   >> 634 G4ChordFinder::PrintStatistics()
646 {                                                 635 {
647    // Dumping the state of G4ChordFinder       << 636   // Print Statistics
648    os << "State of G4ChordFinder : " << std::e << 637 
649    os << "   delta_chord   = " <<  cf.fDeltaCh << 638   G4cout << "G4ChordFinder statistics report: " << G4endl;
650    os << "   Default d_c   = " <<  cf.fDefault << 639   G4cout 
                                                   >> 640     << "  No trials: " << fTotalNoTrials_FNC
                                                   >> 641     << "  No Calls: "  << fNoCalls_FNC
                                                   >> 642     << "  Max-trial: " <<  fmaxTrials_FNC
                                                   >> 643     << G4endl; 
                                                   >> 644   G4cout 
                                                   >> 645     << "  Parameters: " 
                                                   >> 646     << "  fFirstFraction "  << fFirstFraction
                                                   >> 647     << "  fFractionLast "   << fFractionLast
                                                   >> 648     << "  fFractionNextEstimate " << fFractionNextEstimate
                                                   >> 649     << G4endl; 
                                                   >> 650 }
                                                   >> 651 
651                                                   652 
652    os << "   stats-verbose = " <<  cf.fStatsVe << 653 // ...........................................................................
653                                                   654 
654    return os;                                  << 655 void G4ChordFinder::TestChordPrint( G4int    noTrials, 
                                                   >> 656                                     G4int    lastStepTrial, 
                                                   >> 657                                     G4double dChordStep, 
                                                   >> 658                                     G4double nextStepTrial )
                                                   >> 659 {
                                                   >> 660      G4int oldprec= G4cout.precision(5);
                                                   >> 661      G4cout << " ChF/fnc: notrial " << std::setw( 3) << noTrials 
                                                   >> 662             << " this_step= "       << std::setw(10) << lastStepTrial;
                                                   >> 663      if( std::fabs( (dChordStep / fDeltaChord) - 1.0 ) < 0.001 )
                                                   >> 664      {
                                                   >> 665        G4cout.precision(8);
                                                   >> 666      }
                                                   >> 667      else
                                                   >> 668      {
                                                   >> 669        G4cout.precision(6);
                                                   >> 670      }
                                                   >> 671      G4cout << " dChordStep=  " << std::setw(12) << dChordStep;
                                                   >> 672      if( dChordStep > fDeltaChord ) { G4cout << " d+"; }
                                                   >> 673      else                           { G4cout << " d-"; }
                                                   >> 674      G4cout.precision(5);
                                                   >> 675      G4cout <<  " new_step= "       << std::setw(10)
                                                   >> 676             << fLastStepEstimate_Unconstrained
                                                   >> 677             << " new_step_constr= " << std::setw(10)
                                                   >> 678             << lastStepTrial << G4endl;
                                                   >> 679      G4cout << " nextStepTrial = " << std::setw(10) << nextStepTrial << G4endl;
                                                   >> 680      G4cout.precision(oldprec);
655 }                                                 681 }
656                                                   682