Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/magneticfield/src/G4MagIntegratorDriver.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/G4MagIntegratorDriver.cc (Version 11.3.0) and /geometry/magneticfield/src/G4MagIntegratorDriver.cc (Version 10.6)


  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 // G4MagInt_Driver implementation                  26 // G4MagInt_Driver implementation
 27 //                                                 27 //
 28 // V.Grichine, 07.10.1996 - Created                28 // V.Grichine, 07.10.1996 - Created
 29 // W.Wander, 28.01.1998 - Added ability for lo     29 // W.Wander, 28.01.1998 - Added ability for low order integrators
 30 // J.Apostolakis, 08.11.2001 - Respect minimum     30 // J.Apostolakis, 08.11.2001 - Respect minimum step in AccurateAdvance
 31 // -------------------------------------------     31 // --------------------------------------------------------------------
 32                                                    32 
 33 #include <iomanip>                                 33 #include <iomanip>
 34                                                    34 
 35 #include "globals.hh"                              35 #include "globals.hh"
 36 #include "G4SystemOfUnits.hh"                      36 #include "G4SystemOfUnits.hh"
 37 #include "G4GeometryTolerance.hh"                  37 #include "G4GeometryTolerance.hh"
 38 #include "G4MagIntegratorDriver.hh"                38 #include "G4MagIntegratorDriver.hh"
 39 #include "G4FieldTrack.hh"                         39 #include "G4FieldTrack.hh"
 40                                                    40 
 41 #ifdef   G4DEBUG_FIELD                         << 
 42 #include "G4DriverReporter.hh"                 << 
 43 #endif                                         << 
 44                                                << 
 45 // -------------------------------------------     41 // ---------------------------------------------------------
 46                                                    42 
 47 //  Constructor                                    43 //  Constructor
 48 //                                                 44 //
 49 G4MagInt_Driver::G4MagInt_Driver( G4double         45 G4MagInt_Driver::G4MagInt_Driver( G4double                hminimum, 
 50                                   G4MagIntegra     46                                   G4MagIntegratorStepper* pStepper,
 51                                   G4int            47                                   G4int                   numComponents,
 52                                   G4int            48                                   G4int                   statisticsVerbose)
 53   : fNoIntegrationVariables(numComponents),        49   : fNoIntegrationVariables(numComponents), 
 54     fNoVars( std::max( fNoIntegrationVariables     50     fNoVars( std::max( fNoIntegrationVariables, fMinNoVars )),
 55     fStatisticsVerboseLevel(statisticsVerbose)     51     fStatisticsVerboseLevel(statisticsVerbose)
 56 {                                                  52 {  
 57   // In order to accomodate "Laboratory Time",     53   // In order to accomodate "Laboratory Time", which is [7], fMinNoVars=8
 58   // is required. For proper time of flight an     54   // is required. For proper time of flight and spin,  fMinNoVars must be 12
 59                                                    55 
 60   RenewStepperAndAdjust( pStepper );               56   RenewStepperAndAdjust( pStepper );
 61   fMinimumStep = hminimum;                         57   fMinimumStep = hminimum;
 62                                                    58 
 63   fMaxNoSteps = fMaxStepBase / pIntStepper->In     59   fMaxNoSteps = fMaxStepBase / pIntStepper->IntegratorOrder();
 64 #ifdef G4DEBUG_FIELD                               60 #ifdef G4DEBUG_FIELD
 65   fVerboseLevel=2;                                 61   fVerboseLevel=2;
 66 #endif                                             62 #endif
 67                                                    63 
 68   if( (fVerboseLevel > 0) || (fStatisticsVerbo     64   if( (fVerboseLevel > 0) || (fStatisticsVerboseLevel > 1) )
 69   {                                                65   {
 70     G4cout << "MagIntDriver version: Accur-Adv     66     G4cout << "MagIntDriver version: Accur-Adv: "
 71            << "invE_nS, QuickAdv-2sqrt with St     67            << "invE_nS, QuickAdv-2sqrt with Statistics "
 72 #ifdef G4FLD_STATS                                 68 #ifdef G4FLD_STATS
 73            << " enabled "                          69            << " enabled "
 74 #else                                              70 #else
 75            << " disabled "                         71            << " disabled "
 76 #endif                                             72 #endif
 77            << G4endl;                              73            << G4endl;
 78   }                                                74   }
 79 }                                                  75 }
 80                                                    76 
 81 // -------------------------------------------     77 // ---------------------------------------------------------
 82                                                    78 
 83 //  Destructor                                     79 //  Destructor
 84 //                                                 80 //
 85 G4MagInt_Driver::~G4MagInt_Driver()                81 G4MagInt_Driver::~G4MagInt_Driver()
 86 {                                                  82 { 
 87   if( fStatisticsVerboseLevel > 1 )                83   if( fStatisticsVerboseLevel > 1 )
 88   {                                                84   {
 89     PrintStatisticsReport();                       85     PrintStatisticsReport();
 90   }                                                86   }
 91 }                                                  87 }
 92                                                    88 
 93 // -------------------------------------------     89 // ---------------------------------------------------------
 94                                                    90 
 95 G4bool                                             91 G4bool
 96 G4MagInt_Driver::AccurateAdvance(G4FieldTrack&     92 G4MagInt_Driver::AccurateAdvance(G4FieldTrack& y_current,
 97                                  G4double          93                                  G4double      hstep,
 98                                  G4double          94                                  G4double      eps,
 99                                  G4double          95                                  G4double      hinitial )
100 {                                                  96 {
101   // Runge-Kutta driver with adaptive stepsize     97   // Runge-Kutta driver with adaptive stepsize control. Integrate starting
102   // values at y_current over hstep x2 with ac     98   // values at y_current over hstep x2 with accuracy eps. 
103   // On output ystart is replaced by values at     99   // On output ystart is replaced by values at the end of the integration 
104   // interval. RightHandSide is the right-hand    100   // interval. RightHandSide is the right-hand side of ODE system. 
105   // The source is similar to odeint routine f    101   // The source is similar to odeint routine from NRC p.721-722 .
106                                                   102 
107   G4int nstp, i;                               << 103   G4int nstp, i, no_warnings = 0;
108   G4double x, hnext, hdid, h;                     104   G4double x, hnext, hdid, h;
109                                                   105 
110 #ifdef G4DEBUG_FIELD                              106 #ifdef G4DEBUG_FIELD
111   G4int no_warnings = 0;                       << 
112   static G4int dbg = 1;                           107   static G4int dbg = 1;
113   static G4int nStpPr = 50;   // For debug pri    108   static G4int nStpPr = 50;   // For debug printing of long integrations
114   G4double ySubStepStart[G4FieldTrack::ncompSV    109   G4double ySubStepStart[G4FieldTrack::ncompSVEC];
115   G4FieldTrack  yFldTrkStart(y_current);          110   G4FieldTrack  yFldTrkStart(y_current);
116 #endif                                            111 #endif
117                                                   112 
118   G4double y[G4FieldTrack::ncompSVEC] = {0., 0 << 113   G4double y[G4FieldTrack::ncompSVEC], dydx[G4FieldTrack::ncompSVEC];
119   G4double dydx[G4FieldTrack::ncompSVEC] = {0. << 114   G4double ystart[G4FieldTrack::ncompSVEC], yEnd[G4FieldTrack::ncompSVEC]; 
120   G4double ystart[G4FieldTrack::ncompSVEC] = { << 
121   G4double yEnd[G4FieldTrack::ncompSVEC] = {0. << 
122   G4double  x1, x2;                               115   G4double  x1, x2;
123   G4bool succeeded = true;                     << 116   G4bool succeeded = true, lastStepSucceeded;
124                                                   117 
125   G4double startCurveLength;                      118   G4double startCurveLength;
126                                                   119 
                                                   >> 120   G4int noFullIntegr = 0, noSmallIntegr = 0;
                                                   >> 121   static G4ThreadLocal G4int  noGoodSteps = 0;  // Bad = chord > curve-len 
127   const G4int nvar = fNoVars;                     122   const G4int nvar = fNoVars;
128                                                   123 
129   G4FieldTrack yStartFT(y_current);               124   G4FieldTrack yStartFT(y_current);
130                                                   125 
131   //  Ensure that hstep > 0                       126   //  Ensure that hstep > 0
132   //                                              127   //
133   if( hstep <= 0.0 )                              128   if( hstep <= 0.0 )
134   {                                               129   { 
135     if( hstep == 0.0 )                            130     if( hstep == 0.0 )
136     {                                             131     {
137       std::ostringstream message;                 132       std::ostringstream message;
138       message << "Proposed step is zero; hstep    133       message << "Proposed step is zero; hstep = " << hstep << " !";
139       G4Exception("G4MagInt_Driver::AccurateAd    134       G4Exception("G4MagInt_Driver::AccurateAdvance()", 
140                   "GeomField1001", JustWarning    135                   "GeomField1001", JustWarning, message);
141       return succeeded;                           136       return succeeded; 
142     }                                             137     }
143                                                << 138     else
144     std::ostringstream message;                << 139     { 
145     message << "Invalid run condition." << G4e << 140       std::ostringstream message;
146             << "Proposed step is negative; hst << 141       message << "Invalid run condition." << G4endl
147             << "Requested step cannot be negat << 142               << "Proposed step is negative; hstep = " << hstep << "." << G4endl
148     G4Exception("G4MagInt_Driver::AccurateAdva << 143               << "Requested step cannot be negative! Aborting event.";
149                 "GeomField0003", EventMustBeAb << 144       G4Exception("G4MagInt_Driver::AccurateAdvance()", 
150     return false;                              << 145                   "GeomField0003", EventMustBeAborted, message);
                                                   >> 146       return false;
                                                   >> 147     }
151   }                                               148   }
152                                                   149 
153   y_current.DumpToArray( ystart );                150   y_current.DumpToArray( ystart );
154                                                   151 
155   startCurveLength= y_current.GetCurveLength()    152   startCurveLength= y_current.GetCurveLength();
156   x1= startCurveLength;                           153   x1= startCurveLength; 
157   x2= x1 + hstep;                                 154   x2= x1 + hstep;
158                                                   155 
159   if ( (hinitial > 0.0) && (hinitial < hstep)     156   if ( (hinitial > 0.0) && (hinitial < hstep)
160     && (hinitial > perMillion * hstep) )          157     && (hinitial > perMillion * hstep) )
161   {                                               158   {
162      h = hinitial;                                159      h = hinitial;
163   }                                               160   }
164   else  //  Initial Step size "h" defaults to     161   else  //  Initial Step size "h" defaults to the full interval
165   {                                               162   {
166      h = hstep;                                   163      h = hstep;
167   }                                               164   }
168                                                   165 
169   x = x1;                                         166   x = x1;
170                                                   167 
171   for ( i=0; i<nvar; ++i)  { y[i] = ystart[i];    168   for ( i=0; i<nvar; ++i)  { y[i] = ystart[i]; }
172                                                   169 
173   G4bool lastStep= false;                         170   G4bool lastStep= false;
174   nstp = 1;                                       171   nstp = 1;
175                                                   172 
176   do                                              173   do
177   {                                               174   {
178     G4ThreeVector StartPos( y[0], y[1], y[2] )    175     G4ThreeVector StartPos( y[0], y[1], y[2] );
179                                                   176 
180 #ifdef G4DEBUG_FIELD                              177 #ifdef G4DEBUG_FIELD
181     G4double xSubStepStart= x;                    178     G4double xSubStepStart= x; 
182     for (i=0; i<nvar; ++i)  { ySubStepStart[i]    179     for (i=0; i<nvar; ++i)  { ySubStepStart[i] = y[i]; }
183     yFldTrkStart.LoadFromArray(y, fNoIntegrati    180     yFldTrkStart.LoadFromArray(y, fNoIntegrationVariables);
184     yFldTrkStart.SetCurveLength(x);               181     yFldTrkStart.SetCurveLength(x);
185 #endif                                            182 #endif
186                                                   183 
187     pIntStepper->RightHandSide( y, dydx );        184     pIntStepper->RightHandSide( y, dydx );
188     ++fNoTotalSteps;                              185     ++fNoTotalSteps;
189                                                   186 
190     // Perform the Integration                    187     // Perform the Integration
191     //                                            188     //      
192     if( h > fMinimumStep )                        189     if( h > fMinimumStep )
193     {                                             190     { 
194       OneGoodStep(y,dydx,x,h,eps,hdid,hnext) ;    191       OneGoodStep(y,dydx,x,h,eps,hdid,hnext) ;
195       //--------------------------------------    192       //--------------------------------------
                                                   >> 193       lastStepSucceeded = (hdid == h);   
196 #ifdef G4DEBUG_FIELD                              194 #ifdef G4DEBUG_FIELD
197       if (dbg>2)                                  195       if (dbg>2)
198       {                                           196       {
199          // PrintStatus( ySubStepStart, xSubSt << 197         PrintStatus( ySubStepStart, xSubStepStart, y, x, h,  nstp); // Only
200         G4DriverReporter::PrintStatus( ySubSte << 
201       }                                           198       }
202 #endif                                            199 #endif
203     }                                             200     }
204     else                                          201     else
205     {                                             202     {
206       G4FieldTrack yFldTrk( G4ThreeVector(0,0,    203       G4FieldTrack yFldTrk( G4ThreeVector(0,0,0), 
207                             G4ThreeVector(0,0,    204                             G4ThreeVector(0,0,0), 0., 0., 0., 0. );
208       G4double dchord_step, dyerr, dyerr_len;     205       G4double dchord_step, dyerr, dyerr_len;   // What to do with these ?
209       yFldTrk.LoadFromArray(y, fNoIntegrationV    206       yFldTrk.LoadFromArray(y, fNoIntegrationVariables); 
210       yFldTrk.SetCurveLength( x );                207       yFldTrk.SetCurveLength( x );
211                                                   208 
212       QuickAdvance( yFldTrk, dydx, h, dchord_s    209       QuickAdvance( yFldTrk, dydx, h, dchord_step, dyerr_len );
213       //--------------------------------------    210       //-----------------------------------------------------
214                                                   211 
215       yFldTrk.DumpToArray(y);                     212       yFldTrk.DumpToArray(y);    
216                                                   213 
217 #ifdef G4FLD_STATS                                214 #ifdef G4FLD_STATS
218       ++fNoSmallSteps;                            215       ++fNoSmallSteps; 
219       if ( dyerr_len > fDyerr_max )  { fDyerr_    216       if ( dyerr_len > fDyerr_max )  { fDyerr_max = dyerr_len; }
220       fDyerrPos_smTot += dyerr_len;               217       fDyerrPos_smTot += dyerr_len;
221       fSumH_sm += h;  // Length total for 'sma    218       fSumH_sm += h;  // Length total for 'small' steps
222       if (nstp==1)  { ++fNoInitialSmallSteps;  << 219       if (nstp<=1)  { ++fNoInitialSmallSteps; }
223 #endif                                            220 #endif
224 #ifdef G4DEBUG_FIELD                              221 #ifdef G4DEBUG_FIELD
225       if (dbg>1)                                  222       if (dbg>1)
226       {                                           223       {
227         if(fNoSmallSteps<2) { PrintStatus(ySub    224         if(fNoSmallSteps<2) { PrintStatus(ySubStepStart, x1, y, x, h, -nstp); }
228         G4cout << "Another sub-min step, no "     225         G4cout << "Another sub-min step, no " << fNoSmallSteps 
229                << " of " << fNoTotalSteps << "    226                << " of " << fNoTotalSteps << " this time " << nstp << G4endl; 
230         PrintStatus( ySubStepStart, x1, y, x,     227         PrintStatus( ySubStepStart, x1, y, x, h,  nstp);   // Only this
231         G4cout << " dyerr= " << dyerr_len << "    228         G4cout << " dyerr= " << dyerr_len << " relative = " << dyerr_len / h 
232                << " epsilon= " << eps << " hst    229                << " epsilon= " << eps << " hstep= " << hstep 
233                << " h= " << h << " hmin= " <<     230                << " h= " << h << " hmin= " << fMinimumStep << G4endl;
234       }                                           231       }
235 #endif                                            232 #endif        
236       if( h == 0.0 )                              233       if( h == 0.0 )
237       {                                           234       { 
238         G4Exception("G4MagInt_Driver::Accurate    235         G4Exception("G4MagInt_Driver::AccurateAdvance()",
239                     "GeomField0003", FatalExce    236                     "GeomField0003", FatalException,
240                     "Integration Step became Z    237                     "Integration Step became Zero!"); 
241       }                                           238       }
242       dyerr = dyerr_len / h;                      239       dyerr = dyerr_len / h;
243       hdid = h;                                   240       hdid = h;
244       x += hdid;                                  241       x += hdid;
245                                                   242 
246       // Compute suggested new step               243       // Compute suggested new step
247       hnext = ComputeNewStepSize( dyerr/eps, h    244       hnext = ComputeNewStepSize( dyerr/eps, h);
                                                   >> 245 
                                                   >> 246       // .. hnext= ComputeNewStepSize_WithinLimits( dyerr/eps, h);
                                                   >> 247       lastStepSucceeded = (dyerr<= eps);
248     }                                             248     }
249                                                   249 
                                                   >> 250     if (lastStepSucceeded)  { ++noFullIntegr; }
                                                   >> 251     else                    { ++noSmallIntegr; }
                                                   >> 252 
250     G4ThreeVector EndPos( y[0], y[1], y[2] );     253     G4ThreeVector EndPos( y[0], y[1], y[2] );
251                                                   254 
252 #ifdef  G4DEBUG_FIELD                             255 #ifdef  G4DEBUG_FIELD
253     if( (dbg>0) && (dbg<=2) && (nstp>nStpPr))     256     if( (dbg>0) && (dbg<=2) && (nstp>nStpPr))
254     {                                             257     {
255       if( nstp==nStpPr )  { G4cout << "***** M    258       if( nstp==nStpPr )  { G4cout << "***** Many steps ****" << G4endl; }
256       G4cout << "MagIntDrv: " ;                   259       G4cout << "MagIntDrv: " ; 
257       G4cout << "hdid="  << std::setw(12) << h    260       G4cout << "hdid="  << std::setw(12) << hdid  << " "
258              << "hnext=" << std::setw(12) << h    261              << "hnext=" << std::setw(12) << hnext << " " 
259              << "hstep=" << std::setw(12) << h    262              << "hstep=" << std::setw(12) << hstep << " (requested) " 
260              << G4endl;                           263              << G4endl;
261       PrintStatus( ystart, x1, y, x, h, (nstp=    264       PrintStatus( ystart, x1, y, x, h, (nstp==nStpPr) ? -nstp: nstp); 
262     }                                             265     }
263 #endif                                            266 #endif
264                                                   267 
265     // Check the endpoint                         268     // Check the endpoint
266     G4double endPointDist= (EndPos-StartPos).m    269     G4double endPointDist= (EndPos-StartPos).mag(); 
267     if ( endPointDist >= hdid*(1.+perMillion)     270     if ( endPointDist >= hdid*(1.+perMillion) )
268     {                                             271     {
269       ++fNoBadSteps;                              272       ++fNoBadSteps;
270                                                   273 
271       // Issue a warning only for gross differ    274       // Issue a warning only for gross differences -
272       // we understand how small difference oc    275       // we understand how small difference occur.
273       if ( endPointDist >= hdid*(1.+perThousan    276       if ( endPointDist >= hdid*(1.+perThousand) )
274       {                                           277       { 
275 #ifdef G4DEBUG_FIELD                              278 #ifdef G4DEBUG_FIELD
276         if (dbg)                                  279         if (dbg)
277         {                                         280         {
278           WarnEndPointTooFar ( endPointDist, h    281           WarnEndPointTooFar ( endPointDist, hdid, eps, dbg ); 
279           G4cerr << "  Total steps:  bad " <<     282           G4cerr << "  Total steps:  bad " << fNoBadSteps
280                  << " current h= " << hdid <<  << 283                  << " good " << noGoodSteps << " current h= " << hdid
                                                   >> 284                  << G4endl;
281           PrintStatus( ystart, x1, y, x, hstep    285           PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);  
282         }                                         286         }
283         ++no_warnings;                         << 
284 #endif                                            287 #endif
                                                   >> 288         ++no_warnings;
285       }                                           289       }
286     }                                             290     }
                                                   >> 291     else
                                                   >> 292     {
                                                   >> 293       ++noGoodSteps;
                                                   >> 294     } 
                                                   >> 295 // #endif
287                                                   296 
288     //  Avoid numerous small last steps           297     //  Avoid numerous small last steps
289     if( (h < eps * hstep) || (h < fSmallestFra    298     if( (h < eps * hstep) || (h < fSmallestFraction * startCurveLength) )
290     {                                             299     {
291       // No more integration -- the next step     300       // No more integration -- the next step will not happen
292       lastStep = true;                            301       lastStep = true;  
293     }                                             302     }
294     else                                          303     else
295     {                                             304     {
296       // Check the proposed next stepsize         305       // Check the proposed next stepsize
297       if(std::fabs(hnext) <= Hmin())              306       if(std::fabs(hnext) <= Hmin())
298       {                                           307       {
299 #ifdef  G4DEBUG_FIELD                             308 #ifdef  G4DEBUG_FIELD
300         // If simply a very small interval is     309         // If simply a very small interval is being integrated, do not warn
301         if( (x < x2 * (1-eps) ) &&        // T    310         if( (x < x2 * (1-eps) ) &&        // The last step can be small: OK
302             (std::fabs(hstep) > Hmin()) ) // a    311             (std::fabs(hstep) > Hmin()) ) // and if we are asked, it's OK
303         {                                         312         {
304           if(dbg>0)                               313           if(dbg>0)
305           {                                       314           {
306             WarnSmallStepSize( hnext, hstep, h    315             WarnSmallStepSize( hnext, hstep, h, x-x1, nstp );  
307             PrintStatus( ystart, x1, y, x, hst    316             PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);
308           }                                       317           }
309           ++no_warnings;                          318           ++no_warnings;
310         }                                         319         }
311 #endif                                            320 #endif
312         // Make sure that the next step is at     321         // Make sure that the next step is at least Hmin.
313         h = Hmin();                               322         h = Hmin();
314       }                                           323       }
315       else                                        324       else
316       {                                           325       {
317         h = hnext;                                326         h = hnext;
318       }                                           327       }
319                                                   328 
320       //  Ensure that the next step does not o    329       //  Ensure that the next step does not overshoot
321       if ( x+h > x2 )                             330       if ( x+h > x2 )
322       {                // When stepsize oversh    331       {                // When stepsize overshoots, decrease it!
323         h = x2 - x ;   // Must cope with diffi    332         h = x2 - x ;   // Must cope with difficult rounding-error
324       }                // issues if hstep << x    333       }                // issues if hstep << x2
325                                                   334 
326       if ( h == 0.0 )                             335       if ( h == 0.0 )
327       {                                           336       {
328         // Cannot progress - accept this as la    337         // Cannot progress - accept this as last step - by default
329         lastStep = true;                          338         lastStep = true;
330 #ifdef G4DEBUG_FIELD                              339 #ifdef G4DEBUG_FIELD
331         if (dbg>2)                                340         if (dbg>2)
332         {                                         341         {
333           int prec= G4cout.precision(12);         342           int prec= G4cout.precision(12); 
334           G4cout << "Warning: G4MagIntegratorD    343           G4cout << "Warning: G4MagIntegratorDriver::AccurateAdvance"
335                  << G4endl                        344                  << G4endl
336                  << "  Integration step 'h' be    345                  << "  Integration step 'h' became "
337                  << h << " due to roundoff. "     346                  << h << " due to roundoff. " << G4endl
338                  << " Calculated as difference    347                  << " Calculated as difference of x2= "<< x2 << " and x=" << x
339                  << "  Forcing termination of     348                  << "  Forcing termination of advance." << G4endl;
340           G4cout.precision(prec);                 349           G4cout.precision(prec);
341         }                                         350         }          
342 #endif                                            351 #endif
343       }                                           352       }
344     }                                             353     }
345   } while ( ((++nstp)<=fMaxNoSteps) && (x < x2 << 354   } while ( ((nstp++)<=fMaxNoSteps) && (x < x2) && (!lastStep) );
346   // Loop checking, 07.10.2016, J. Apostolakis    355   // Loop checking, 07.10.2016, J. Apostolakis
347                                                   356 
348      // Have we reached the end ?                 357      // Have we reached the end ?
349      // --> a better test might be x-x2 > an_e    358      // --> a better test might be x-x2 > an_epsilon
350                                                   359 
351   succeeded = (x>=x2);  // If it was a "forced    360   succeeded = (x>=x2);  // If it was a "forced" last step
352                                                   361 
353   for (i=0; i<nvar; ++i)  { yEnd[i] = y[i]; }     362   for (i=0; i<nvar; ++i)  { yEnd[i] = y[i]; }
354                                                   363 
355   // Put back the values.                         364   // Put back the values.
356   y_current.LoadFromArray( yEnd, fNoIntegratio    365   y_current.LoadFromArray( yEnd, fNoIntegrationVariables );
357   y_current.SetCurveLength( x );                  366   y_current.SetCurveLength( x );
358                                                   367 
359   if(nstp > fMaxNoSteps)                          368   if(nstp > fMaxNoSteps)
360   {                                               369   {
                                                   >> 370     ++no_warnings;
361     succeeded = false;                            371     succeeded = false;
362 #ifdef G4DEBUG_FIELD                              372 #ifdef G4DEBUG_FIELD
363     ++no_warnings;                             << 
364     if (dbg)                                      373     if (dbg)
365     {                                             374     {
366       WarnTooManyStep( x1, x2, x );  //  Issue    375       WarnTooManyStep( x1, x2, x );  //  Issue WARNING
367       PrintStatus( yEnd, x1, y, x, hstep, -nst    376       PrintStatus( yEnd, x1, y, x, hstep, -nstp);
368     }                                             377     }
369 #endif                                            378 #endif
370   }                                               379   }
371                                                   380 
372 #ifdef G4DEBUG_FIELD                              381 #ifdef G4DEBUG_FIELD
373   if( dbg && no_warnings )                        382   if( dbg && no_warnings )
374   {                                               383   {
375     G4cerr << "G4MagIntegratorDriver exit stat    384     G4cerr << "G4MagIntegratorDriver exit status: no-steps " << nstp << G4endl;
376     PrintStatus( yEnd, x1, y, x, hstep, nstp);    385     PrintStatus( yEnd, x1, y, x, hstep, nstp);
377   }                                               386   }
378 #endif                                            387 #endif
379                                                   388 
380   return succeeded;                               389   return succeeded;
381 }  // end of AccurateAdvance .................    390 }  // end of AccurateAdvance ...........................
382                                                   391 
383 // -------------------------------------------    392 // ---------------------------------------------------------
384                                                   393 
385 void                                              394 void
386 G4MagInt_Driver::WarnSmallStepSize( G4double h    395 G4MagInt_Driver::WarnSmallStepSize( G4double hnext, G4double hstep, 
387                                     G4double h    396                                     G4double h, G4double xDone,
388                                     G4int nstp    397                                     G4int nstp)
389 {                                                 398 {
390   static G4ThreadLocal G4int noWarningsIssued     399   static G4ThreadLocal G4int noWarningsIssued = 0;
391   const  G4int maxNoWarnings = 10;   // Number    400   const  G4int maxNoWarnings = 10;   // Number of verbose warnings
392   std::ostringstream message;                     401   std::ostringstream message;
393   if( (noWarningsIssued < maxNoWarnings) || fV    402   if( (noWarningsIssued < maxNoWarnings) || fVerboseLevel > 10 )
394   {                                               403   {
395     message << "The stepsize for the next iter    404     message << "The stepsize for the next iteration, " << hnext
396             << ", is too small - in Step numbe    405             << ", is too small - in Step number " << nstp << "." << G4endl
397             << "The minimum for the driver is     406             << "The minimum for the driver is " << Hmin()  << G4endl
398             << "Requested integr. length was "    407             << "Requested integr. length was " << hstep << " ." << G4endl
399             << "The size of this sub-step was     408             << "The size of this sub-step was " << h     << " ." << G4endl
400             << "The integrations has already g    409             << "The integrations has already gone " << xDone;
401   }                                               410   }
402   else                                            411   else
403   {                                               412   {
404     message << "Too small 'next' step " << hne    413     message << "Too small 'next' step " << hnext
405             << ", step-no: " << nstp << G4endl    414             << ", step-no: " << nstp << G4endl
406             << ", this sub-step: " << h           415             << ", this sub-step: " << h     
407             << ",  req_tot_len: " << hstep        416             << ",  req_tot_len: " << hstep 
408             << ", done: " << xDone << ", min:     417             << ", done: " << xDone << ", min: " << Hmin();
409   }                                               418   }
410   G4Exception("G4MagInt_Driver::WarnSmallStepS    419   G4Exception("G4MagInt_Driver::WarnSmallStepSize()", "GeomField1001",
411               JustWarning, message);              420               JustWarning, message);
412   ++noWarningsIssued;                             421   ++noWarningsIssued;
413 }                                                 422 }
414                                                   423 
415 // -------------------------------------------    424 // ---------------------------------------------------------
416                                                   425 
417 void                                              426 void
418 G4MagInt_Driver::WarnTooManyStep( G4double x1s    427 G4MagInt_Driver::WarnTooManyStep( G4double x1start, 
419                                   G4double x2e    428                                   G4double x2end, 
420                                   G4double xCu    429                                   G4double xCurrent )
421 {                                                 430 {
422    std::ostringstream message;                    431    std::ostringstream message;
423    message << "The number of steps used in the    432    message << "The number of steps used in the Integration driver"
424            << " (Runge-Kutta) is too many." <<    433            << " (Runge-Kutta) is too many." << G4endl
425            << "Integration of the interval was    434            << "Integration of the interval was not completed !" << G4endl
426            << "Only a " << (xCurrent-x1start)*    435            << "Only a " << (xCurrent-x1start)*100/(x2end-x1start)
427            << " % fraction of it was done.";      436            << " % fraction of it was done.";
428    G4Exception("G4MagInt_Driver::WarnTooManySt    437    G4Exception("G4MagInt_Driver::WarnTooManyStep()", "GeomField1001",
429                JustWarning, message);             438                JustWarning, message);
430 }                                                 439 }
431                                                   440 
432 // -------------------------------------------    441 // ---------------------------------------------------------
433                                                   442 
434 void                                              443 void
435 G4MagInt_Driver::WarnEndPointTooFar (G4double     444 G4MagInt_Driver::WarnEndPointTooFar (G4double endPointDist, 
436                                      G4double     445                                      G4double h , 
437                                      G4double     446                                      G4double eps,
438                                      G4int        447                                      G4int    dbg)
439 {                                                 448 {
440   static G4ThreadLocal G4double maxRelError =     449   static G4ThreadLocal G4double maxRelError = 0.0;
441   G4bool isNewMax, prNewMax;                      450   G4bool isNewMax, prNewMax;
442                                                   451 
443   isNewMax = endPointDist > (1.0 + maxRelError    452   isNewMax = endPointDist > (1.0 + maxRelError) * h;
444   prNewMax = endPointDist > (1.0 + 1.05 * maxR    453   prNewMax = endPointDist > (1.0 + 1.05 * maxRelError) * h;
445   if( isNewMax ) { maxRelError= endPointDist /    454   if( isNewMax ) { maxRelError= endPointDist / h - 1.0; }
446                                                   455 
447   if( (dbg != 0) && (h > G4GeometryTolerance:: << 456   if( dbg && (h > G4GeometryTolerance::GetInstance()->GetSurfaceTolerance()) 
448           && ( (dbg>1) || prNewMax || (endPoin    457           && ( (dbg>1) || prNewMax || (endPointDist >= h*(1.+eps) ) ) )
449   {                                               458   { 
450     static G4ThreadLocal G4int noWarnings = 0;    459     static G4ThreadLocal G4int noWarnings = 0;
451     std::ostringstream message;                   460     std::ostringstream message;
452     if( (noWarnings++ < 10) || (dbg>2) )          461     if( (noWarnings++ < 10) || (dbg>2) )
453     {                                             462     {
454       message << "The integration produced an     463       message << "The integration produced an end-point which " << G4endl
455               << "is further from the start-po    464               << "is further from the start-point than the curve length."
456               << G4endl;                          465               << G4endl;
457     }                                             466     }
458     message << "  Distance of endpoints = " <<    467     message << "  Distance of endpoints = " << endPointDist
459             << ", curve length = " << h << G4e    468             << ", curve length = " << h << G4endl
460             << "  Difference (curveLen-endpDis    469             << "  Difference (curveLen-endpDist)= " << (h - endPointDist)
461             << ", relative = " << (h-endPointD    470             << ", relative = " << (h-endPointDist) / h 
462             << ", epsilon =  " << eps;            471             << ", epsilon =  " << eps;
463     G4Exception("G4MagInt_Driver::WarnEndPoint    472     G4Exception("G4MagInt_Driver::WarnEndPointTooFar()", "GeomField1001",
464                 JustWarning, message);            473                 JustWarning, message);
465   }                                               474   }
466 }                                                 475 }
467                                                   476 
468 // -------------------------------------------    477 // ---------------------------------------------------------
469                                                   478 
470 void                                              479 void
471 G4MagInt_Driver::OneGoodStep(      G4double y[    480 G4MagInt_Driver::OneGoodStep(      G4double y[],        // InOut
472                              const G4double dy    481                              const G4double dydx[],
473                                    G4double& x    482                                    G4double& x,         // InOut
474                                    G4double ht    483                                    G4double htry,
475                                    G4double ep    484                                    G4double eps_rel_max,
476                                    G4double& h    485                                    G4double& hdid,      // Out
477                                    G4double& h    486                                    G4double& hnext )    // Out
478                                                   487 
479 // Driver for one Runge-Kutta Step with monito    488 // Driver for one Runge-Kutta Step with monitoring of local truncation error
480 // to ensure accuracy and adjust stepsize. Inp    489 // to ensure accuracy and adjust stepsize. Input are dependent variable
481 // array y[0,...,5] and its derivative dydx[0,    490 // array y[0,...,5] and its derivative dydx[0,...,5] at the
482 // starting value of the independent variable     491 // starting value of the independent variable x . Also input are stepsize
483 // to be attempted htry, and the required accu    492 // to be attempted htry, and the required accuracy eps. On output y and x
484 // are replaced by their new values, hdid is t    493 // are replaced by their new values, hdid is the stepsize that was actually
485 // accomplished, and hnext is the estimated ne    494 // accomplished, and hnext is the estimated next stepsize. 
486 // This is similar to the function rkqs from t    495 // This is similar to the function rkqs from the book:
487 // Numerical Recipes in C: The Art of Scientif    496 // Numerical Recipes in C: The Art of Scientific Computing (NRC), Second
488 // Edition, by William H. Press, Saul A. Teuko    497 // Edition, by William H. Press, Saul A. Teukolsky, William T.
489 // Vetterling, and Brian P. Flannery (Cambridg    498 // Vetterling, and Brian P. Flannery (Cambridge University Press 1992),
490 // 16.2 Adaptive StepSize Control for Runge-Ku    499 // 16.2 Adaptive StepSize Control for Runge-Kutta, p. 719
491                                                   500 
492 {                                                 501 {
493   G4double errmax_sq;                             502   G4double errmax_sq;
494   G4double h, htemp, xnew ;                       503   G4double h, htemp, xnew ;
495                                                   504 
496   G4double yerr[G4FieldTrack::ncompSVEC], ytem    505   G4double yerr[G4FieldTrack::ncompSVEC], ytemp[G4FieldTrack::ncompSVEC];
497                                                   506 
498   h = htry ; // Set stepsize to the initial tr    507   h = htry ; // Set stepsize to the initial trial value
499                                                   508 
500   G4double inv_eps_vel_sq = 1.0 / (eps_rel_max    509   G4double inv_eps_vel_sq = 1.0 / (eps_rel_max*eps_rel_max);
501                                                   510 
502   G4double errpos_sq = 0.0;    // square of di    511   G4double errpos_sq = 0.0;    // square of displacement error
503   G4double errvel_sq = 0.0;    // square of mo    512   G4double errvel_sq = 0.0;    // square of momentum vector difference
504   G4double errspin_sq = 0.0;   // square of sp    513   G4double errspin_sq = 0.0;   // square of spin vector difference
505                                                   514 
                                                   >> 515   static G4ThreadLocal G4int tot_no_trials=0; 
506   const G4int max_trials=100;                     516   const G4int max_trials=100; 
507                                                   517 
508   G4ThreeVector Spin(y[9],y[10],y[11]);           518   G4ThreeVector Spin(y[9],y[10],y[11]);
509   G4double spin_mag2 = Spin.mag2();               519   G4double spin_mag2 = Spin.mag2();
510   G4bool hasSpin = (spin_mag2 > 0.0);             520   G4bool hasSpin = (spin_mag2 > 0.0); 
511                                                   521 
512   for (G4int iter=0; iter<max_trials; ++iter)     522   for (G4int iter=0; iter<max_trials; ++iter)
513   {                                               523   {
                                                   >> 524     ++tot_no_trials;
514     pIntStepper-> Stepper(y,dydx,h,ytemp,yerr)    525     pIntStepper-> Stepper(y,dydx,h,ytemp,yerr); 
515     //            *******                         526     //            *******
516     G4double eps_pos = eps_rel_max * std::max(    527     G4double eps_pos = eps_rel_max * std::max(h, fMinimumStep); 
517     G4double inv_eps_pos_sq = 1.0 / (eps_pos*e    528     G4double inv_eps_pos_sq = 1.0 / (eps_pos*eps_pos); 
518                                                   529 
519     // Evaluate accuracy                          530     // Evaluate accuracy
520     //                                            531     //
521     errpos_sq =  sqr(yerr[0]) + sqr(yerr[1]) +    532     errpos_sq =  sqr(yerr[0]) + sqr(yerr[1]) + sqr(yerr[2]) ;
522     errpos_sq *= inv_eps_pos_sq; // Scale rela    533     errpos_sq *= inv_eps_pos_sq; // Scale relative to required tolerance
523                                                   534 
524     // Accuracy for momentum                      535     // Accuracy for momentum
525     G4double magvel_sq=  sqr(y[3]) + sqr(y[4])    536     G4double magvel_sq=  sqr(y[3]) + sqr(y[4]) + sqr(y[5]) ;
526     G4double sumerr_sq =  sqr(yerr[3]) + sqr(y    537     G4double sumerr_sq =  sqr(yerr[3]) + sqr(yerr[4]) + sqr(yerr[5]) ; 
527     if( magvel_sq > 0.0 )                         538     if( magvel_sq > 0.0 )
528     {                                             539     { 
529        errvel_sq = sumerr_sq / magvel_sq;         540        errvel_sq = sumerr_sq / magvel_sq; 
530     }                                             541     }
531     else                                          542     else
532     {                                             543     {
533        std::ostringstream message;                544        std::ostringstream message;
534        message << "Found case of zero momentum    545        message << "Found case of zero momentum." << G4endl
535                << "- iteration= " << iter << "    546                << "- iteration= " << iter << "; h= " << h;
536        G4Exception("G4MagInt_Driver::OneGoodSt    547        G4Exception("G4MagInt_Driver::OneGoodStep()",
537                    "GeomField1001", JustWarnin    548                    "GeomField1001", JustWarning, message);
538        errvel_sq = sumerr_sq;                     549        errvel_sq = sumerr_sq; 
539     }                                             550     }
540     errvel_sq *= inv_eps_vel_sq;                  551     errvel_sq *= inv_eps_vel_sq;
541     errmax_sq = std::max( errpos_sq, errvel_sq    552     errmax_sq = std::max( errpos_sq, errvel_sq ); // Square of maximum error
542                                                   553 
543     if( hasSpin )                                 554     if( hasSpin )
544     {                                             555     { 
545       // Accuracy for spin                        556       // Accuracy for spin
546       errspin_sq =  ( sqr(yerr[9]) + sqr(yerr[    557       errspin_sq =  ( sqr(yerr[9]) + sqr(yerr[10]) + sqr(yerr[11]) )
547                     /  spin_mag2; // ( sqr(y[9    558                     /  spin_mag2; // ( sqr(y[9]) + sqr(y[10]) + sqr(y[11]) );
548       errspin_sq *= inv_eps_vel_sq;               559       errspin_sq *= inv_eps_vel_sq;
549       errmax_sq = std::max( errmax_sq, errspin    560       errmax_sq = std::max( errmax_sq, errspin_sq ); 
550     }                                             561     }
551                                                   562 
552     if ( errmax_sq <= 1.0 )  { break; } // Ste    563     if ( errmax_sq <= 1.0 )  { break; } // Step succeeded. 
553                                                   564 
554     // Step failed; compute the size of retria    565     // Step failed; compute the size of retrial Step.
555     htemp = GetSafety() * h * std::pow( errmax    566     htemp = GetSafety() * h * std::pow( errmax_sq, 0.5*GetPshrnk() );
556                                                   567 
557     if (htemp >= 0.1*h)  { h = htemp; }  // Tr    568     if (htemp >= 0.1*h)  { h = htemp; }  // Truncation error too large,
558     else  { h = 0.1*h; }                 // re    569     else  { h = 0.1*h; }                 // reduce stepsize, but no more
559                                          // th    570                                          // than a factor of 10
560     xnew = x + h;                                 571     xnew = x + h;
561     if(xnew == x)                                 572     if(xnew == x)
562     {                                             573     {
563       std::ostringstream message;                 574       std::ostringstream message;
564       message << "Stepsize underflow in Steppe    575       message << "Stepsize underflow in Stepper !" << G4endl
565               << "- Step's start x=" << x << "    576               << "- Step's start x=" << x << " and end x= " << xnew 
566               << " are equal !! " << G4endl       577               << " are equal !! " << G4endl
567               << "  Due to step-size= " << h      578               << "  Due to step-size= " << h 
568               << ". Note that input step was "    579               << ". Note that input step was " << htry;
569       G4Exception("G4MagInt_Driver::OneGoodSte    580       G4Exception("G4MagInt_Driver::OneGoodStep()",
570                   "GeomField1001", JustWarning    581                   "GeomField1001", JustWarning, message);
571       break;                                      582       break;
572     }                                             583     }
573   }                                               584   }
574                                                   585 
575   // Compute size of next Step                    586   // Compute size of next Step
576   if (errmax_sq > errcon*errcon)                  587   if (errmax_sq > errcon*errcon)
577   {                                               588   { 
578     hnext = GetSafety()*h*std::pow(errmax_sq,     589     hnext = GetSafety()*h*std::pow(errmax_sq, 0.5*GetPgrow());
579   }                                               590   }
580   else                                            591   else
581   {                                               592   {
582     hnext = max_stepping_increase*h ; // No mo    593     hnext = max_stepping_increase*h ; // No more than a factor of 5 increase
583   }                                               594   }
584   x += (hdid = h);                                595   x += (hdid = h);
585                                                   596 
586   for(G4int k=0; k<fNoIntegrationVariables; ++    597   for(G4int k=0; k<fNoIntegrationVariables; ++k) { y[k] = ytemp[k]; }
587                                                   598 
588   return;                                         599   return;
589 }                                                 600 }
590                                                   601 
591 //--------------------------------------------    602 //----------------------------------------------------------------------
592                                                   603 
593 // QuickAdvance just tries one Step - it does     604 // QuickAdvance just tries one Step - it does not ensure accuracy
594 //                                                605 //
595 G4bool G4MagInt_Driver::QuickAdvance(G4FieldTr    606 G4bool G4MagInt_Driver::QuickAdvance(G4FieldTrack& y_posvel,    // INOUT
596                                const G4double     607                                const G4double      dydx[],  
597                                      G4double     608                                      G4double      hstep,       // In
598                                      G4double&    609                                      G4double&     dchord_step,
599                                      G4double&    610                                      G4double&     dyerr_pos_sq,
600                                      G4double&    611                                      G4double&     dyerr_mom_rel_sq )  
601 {                                                 612 {
602   G4Exception("G4MagInt_Driver::QuickAdvance()    613   G4Exception("G4MagInt_Driver::QuickAdvance()", "GeomField0001",
603               FatalException, "Not yet impleme    614               FatalException, "Not yet implemented."); 
604                                                   615 
605   // Use the parameters of this method, to ple    616   // Use the parameters of this method, to please compiler
606   //                                              617   //
607   dchord_step = dyerr_pos_sq = hstep * hstep *    618   dchord_step = dyerr_pos_sq = hstep * hstep * dydx[0]; 
608   dyerr_mom_rel_sq = y_posvel.GetPosition().ma    619   dyerr_mom_rel_sq = y_posvel.GetPosition().mag2();
609   return true;                                    620   return true;
610 }                                                 621 }
611                                                   622 
612 //--------------------------------------------    623 //----------------------------------------------------------------------
613                                                   624 
614 G4bool G4MagInt_Driver::QuickAdvance(G4FieldTr    625 G4bool G4MagInt_Driver::QuickAdvance(G4FieldTrack& y_posvel,    // INOUT
615                                const G4double     626                                const G4double      dydx[],  
616                                      G4double     627                                      G4double      hstep,       // In
617                                      G4double&    628                                      G4double&     dchord_step,
618                                      G4double&    629                                      G4double&     dyerr )
619 {                                                 630 {
620   G4double dyerr_pos_sq, dyerr_mom_rel_sq;        631   G4double dyerr_pos_sq, dyerr_mom_rel_sq;  
621   G4double yerr_vec[G4FieldTrack::ncompSVEC],     632   G4double yerr_vec[G4FieldTrack::ncompSVEC],
622            yarrin[G4FieldTrack::ncompSVEC], ya    633            yarrin[G4FieldTrack::ncompSVEC], yarrout[G4FieldTrack::ncompSVEC]; 
623   G4double s_start;                               634   G4double s_start;
624   G4double dyerr_mom_sq, vel_mag_sq, inv_vel_m    635   G4double dyerr_mom_sq, vel_mag_sq, inv_vel_mag_sq;
625                                                   636 
                                                   >> 637   static G4ThreadLocal G4int no_call = 0; 
                                                   >> 638   ++no_call; 
                                                   >> 639 
626   // Move data into array                         640   // Move data into array
627   y_posvel.DumpToArray( yarrin );      //  yar    641   y_posvel.DumpToArray( yarrin );      //  yarrin  <== y_posvel 
628   s_start = y_posvel.GetCurveLength();            642   s_start = y_posvel.GetCurveLength();
629                                                   643 
630   // Do an Integration Step                       644   // Do an Integration Step
631   pIntStepper-> Stepper(yarrin, dydx, hstep, y    645   pIntStepper-> Stepper(yarrin, dydx, hstep, yarrout, yerr_vec) ; 
632                                                   646 
633   // Estimate curve-chord distance                647   // Estimate curve-chord distance
634   dchord_step= pIntStepper-> DistChord();         648   dchord_step= pIntStepper-> DistChord();
635                                                   649 
636   // Put back the values.  yarrout ==> y_posve    650   // Put back the values.  yarrout ==> y_posvel
637   y_posvel.LoadFromArray( yarrout, fNoIntegrat    651   y_posvel.LoadFromArray( yarrout, fNoIntegrationVariables );
638   y_posvel.SetCurveLength( s_start + hstep );     652   y_posvel.SetCurveLength( s_start + hstep );
639                                                   653 
640 #ifdef  G4DEBUG_FIELD                             654 #ifdef  G4DEBUG_FIELD
641   if(fVerboseLevel>2)                             655   if(fVerboseLevel>2)
642   {                                               656   {
643     G4cout << "G4MagIntDrv: Quick Advance" <<     657     G4cout << "G4MagIntDrv: Quick Advance" << G4endl;
644     PrintStatus( yarrin, s_start, yarrout, s_s    658     PrintStatus( yarrin, s_start, yarrout, s_start+hstep, hstep,  1); 
645   }                                               659   }
646 #endif                                            660 #endif
647                                                   661 
648   // A single measure of the error                662   // A single measure of the error   
649   //      TO-DO :  account for  energy,  spin,    663   //      TO-DO :  account for  energy,  spin, ... ? 
650   vel_mag_sq   = ( sqr(yarrout[3])+sqr(yarrout    664   vel_mag_sq   = ( sqr(yarrout[3])+sqr(yarrout[4])+sqr(yarrout[5]) );
651   inv_vel_mag_sq = 1.0 / vel_mag_sq;              665   inv_vel_mag_sq = 1.0 / vel_mag_sq; 
652   dyerr_pos_sq = ( sqr(yerr_vec[0])+sqr(yerr_v    666   dyerr_pos_sq = ( sqr(yerr_vec[0])+sqr(yerr_vec[1])+sqr(yerr_vec[2]));
653   dyerr_mom_sq = ( sqr(yerr_vec[3])+sqr(yerr_v    667   dyerr_mom_sq = ( sqr(yerr_vec[3])+sqr(yerr_vec[4])+sqr(yerr_vec[5]));
654   dyerr_mom_rel_sq = dyerr_mom_sq * inv_vel_ma    668   dyerr_mom_rel_sq = dyerr_mom_sq * inv_vel_mag_sq;
655                                                   669 
656   // Calculate also the change in the momentum    670   // Calculate also the change in the momentum squared also ???
657   // G4double veloc_square = y_posvel.GetVeloc    671   // G4double veloc_square = y_posvel.GetVelocity().mag2();
658   // ...                                          672   // ...
659                                                   673 
660 #ifdef RETURN_A_NEW_STEP_LENGTH                   674 #ifdef RETURN_A_NEW_STEP_LENGTH
661   // The following step cannot be done here be    675   // The following step cannot be done here because "eps" is not known.
662   dyerr_len = std::sqrt( dyerr_len_sq );          676   dyerr_len = std::sqrt( dyerr_len_sq ); 
663   dyerr_len_sq /= eps ;                           677   dyerr_len_sq /= eps ;
664                                                   678 
665   // Look at the velocity deviation ?             679   // Look at the velocity deviation ?
666   //  sqr(yerr_vec[3])+sqr(yerr_vec[4])+sqr(ye    680   //  sqr(yerr_vec[3])+sqr(yerr_vec[4])+sqr(yerr_vec[5]));
667                                                   681 
668   // Set suggested new step                       682   // Set suggested new step
669   hstep = ComputeNewStepSize( dyerr_len, hstep    683   hstep = ComputeNewStepSize( dyerr_len, hstep);
670 #endif                                            684 #endif
671                                                   685 
672   if( dyerr_pos_sq > ( dyerr_mom_rel_sq * sqr(    686   if( dyerr_pos_sq > ( dyerr_mom_rel_sq * sqr(hstep) ) )
673   {                                               687   {
674     dyerr = std::sqrt(dyerr_pos_sq);              688     dyerr = std::sqrt(dyerr_pos_sq);
675   }                                               689   }
676   else                                            690   else
677   {                                               691   {
678     // Scale it to the current step size - for    692     // Scale it to the current step size - for now
679     dyerr = std::sqrt(dyerr_mom_rel_sq) * hste    693     dyerr = std::sqrt(dyerr_mom_rel_sq) * hstep;
680   }                                               694   }
681                                                   695 
682   return true;                                    696   return true;
683 }                                                 697 }
684                                                   698 
685 // -------------------------------------------    699 // --------------------------------------------------------------------------
686                                                   700 
687 #ifdef QUICK_ADV_ARRAY_IN_AND_OUT                 701 #ifdef QUICK_ADV_ARRAY_IN_AND_OUT
688 G4bool  G4MagInt_Driver::QuickAdvance(G4double    702 G4bool  G4MagInt_Driver::QuickAdvance(G4double  yarrin[],    // In
689                                 const G4double    703                                 const G4double  dydx[],  
690                                       G4double    704                                       G4double  hstep,       // In
691                                       G4double    705                                       G4double  yarrout[],
692                                       G4double    706                                       G4double& dchord_step,
693                                       G4double    707                                       G4double& dyerr )      // In length
694 {                                                 708 {
695   G4Exception("G4MagInt_Driver::QuickAdvance()    709   G4Exception("G4MagInt_Driver::QuickAdvance()", "GeomField0001",
696               FatalException, "Not yet impleme    710               FatalException, "Not yet implemented.");
697   dyerr = dchord_step = hstep * yarrin[0] * dy    711   dyerr = dchord_step = hstep * yarrin[0] * dydx[0];
698   yarrout[0]= yarrin[0];                          712   yarrout[0]= yarrin[0];
699 }                                                 713 }
700 #endif                                            714 #endif 
701                                                   715 
702 // -------------------------------------------    716 // --------------------------------------------------------------------------
703                                                   717 
704 // This method computes new step sizes - but d    718 // This method computes new step sizes - but does not limit changes to
705 // within  certain factors                        719 // within  certain factors
706 //                                                720 // 
707 G4double G4MagInt_Driver::                        721 G4double G4MagInt_Driver::
708 ComputeNewStepSize_WithoutReductionLimit(G4dou << 722 ComputeNewStepSize(G4double  errMaxNorm,    // max error  (normalised)
709                    G4double  hstepCurrent)  //    723                    G4double  hstepCurrent)  // current step size
710 {                                                 724 {
711   G4double hnew;                                  725   G4double hnew;
712                                                   726 
713   // Compute size of next Step for a failed st    727   // Compute size of next Step for a failed step
714   if(errMaxNorm > 1.0 )                           728   if(errMaxNorm > 1.0 )
715   {                                               729   {
716     // Step failed; compute the size of retria    730     // Step failed; compute the size of retrial Step.
717     hnew = GetSafety()*hstepCurrent*std::pow(e    731     hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPshrnk()) ;
718   }                                               732   }
719   else if(errMaxNorm > 0.0 )                      733   else if(errMaxNorm > 0.0 )
720   {                                               734   {
721     // Compute size of next Step for a success    735     // Compute size of next Step for a successful step
722     hnew = GetSafety()*hstepCurrent*std::pow(e    736     hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPgrow()) ;
723   }                                               737   }
724   else                                            738   else
725   {                                               739   {
726     // if error estimate is zero (possible) or    740     // if error estimate is zero (possible) or negative (dubious)
727     hnew = max_stepping_increase * hstepCurren    741     hnew = max_stepping_increase * hstepCurrent; 
728   }                                               742   }
729                                                   743 
730   return hnew;                                    744   return hnew;
731 }                                                 745 }
732                                                   746 
733 // -------------------------------------------    747 // ---------------------------------------------------------------------------
734                                                   748 
735 G4double                                       << 
736 G4MagInt_Driver::ComputeNewStepSize(           << 
737                           G4double  errMaxNorm << 
738                           G4double  hstepCurre << 
739 {                                              << 
740    // Legacy behaviour:                        << 
741    return ComputeNewStepSize_WithoutReductionL << 
742    // 'Improved' behaviour - at least more con << 
743    // return ComputeNewStepSize_WithinLimits(  << 
744 }                                              << 
745                                                << 
746 // This method computes new step sizes limitin    749 // This method computes new step sizes limiting changes within certain factors
747 //                                                750 // 
748 // It shares its logic with AccurateAdvance.      751 // It shares its logic with AccurateAdvance.
749 // They are kept separate currently for optimi    752 // They are kept separate currently for optimisation.
750 //                                                753 //
751 G4double                                          754 G4double 
752 G4MagInt_Driver::ComputeNewStepSize_WithinLimi    755 G4MagInt_Driver::ComputeNewStepSize_WithinLimits( 
753                           G4double  errMaxNorm    756                           G4double  errMaxNorm,    // max error  (normalised)
754                           G4double  hstepCurre    757                           G4double  hstepCurrent)  // current step size
755 {                                                 758 {
756   G4double hnew;                                  759   G4double hnew;
757                                                   760 
758   // Compute size of next Step for a failed st    761   // Compute size of next Step for a failed step
759   if (errMaxNorm > 1.0 )                          762   if (errMaxNorm > 1.0 )
760   {                                               763   {
761     // Step failed; compute the size of retria    764     // Step failed; compute the size of retrial Step.
762     hnew = GetSafety()*hstepCurrent*std::pow(e    765     hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPshrnk()) ;
763                                                   766   
764     if (hnew < max_stepping_decrease*hstepCurr    767     if (hnew < max_stepping_decrease*hstepCurrent)
765     {                                             768     {
766       hnew = max_stepping_decrease*hstepCurren    769       hnew = max_stepping_decrease*hstepCurrent ;
767                          // reduce stepsize, b    770                          // reduce stepsize, but no more
768                          // than this factor (    771                          // than this factor (value= 1/10)
769     }                                             772     }
770   }                                               773   }
771   else                                            774   else
772   {                                               775   {
773     // Compute size of next Step for a success    776     // Compute size of next Step for a successful step
774     if (errMaxNorm > errcon)                      777     if (errMaxNorm > errcon)
775      { hnew = GetSafety()*hstepCurrent*std::po    778      { hnew = GetSafety()*hstepCurrent*std::pow(errMaxNorm,GetPgrow()); }
776     else  // No more than a factor of 5 increa    779     else  // No more than a factor of 5 increase
777      { hnew = max_stepping_increase * hstepCur    780      { hnew = max_stepping_increase * hstepCurrent; }
778   }                                               781   }
779   return hnew;                                    782   return hnew;
780 }                                                 783 }
781                                                   784 
782 // -------------------------------------------    785 // ---------------------------------------------------------------------------
783                                                   786 
784 void G4MagInt_Driver::PrintStatus( const G4dou    787 void G4MagInt_Driver::PrintStatus( const G4double* StartArr,  
785                                          G4dou    788                                          G4double  xstart,
786                                    const G4dou    789                                    const G4double* CurrentArr, 
787                                          G4dou    790                                          G4double  xcurrent,
788                                          G4dou    791                                          G4double  requestStep, 
789                                          G4int    792                                          G4int     subStepNo )
790   // Potentially add as arguments:                793   // Potentially add as arguments:  
791   //                                 <dydx>       794   //                                 <dydx>           - as Initial Force
792   //                                 stepTaken    795   //                                 stepTaken(hdid)  - last step taken
793   //                                 nextStep     796   //                                 nextStep (hnext) - proposal for size
794 {                                                 797 {
795    G4FieldTrack  StartFT(G4ThreeVector(0,0,0),    798    G4FieldTrack  StartFT(G4ThreeVector(0,0,0),
796                  G4ThreeVector(0,0,0), 0., 0.,    799                  G4ThreeVector(0,0,0), 0., 0., 0., 0. );
797    G4FieldTrack  CurrentFT (StartFT);             800    G4FieldTrack  CurrentFT (StartFT);
798                                                   801 
799    StartFT.LoadFromArray( StartArr, fNoIntegra    802    StartFT.LoadFromArray( StartArr, fNoIntegrationVariables); 
800    StartFT.SetCurveLength( xstart);               803    StartFT.SetCurveLength( xstart);
801    CurrentFT.LoadFromArray( CurrentArr, fNoInt    804    CurrentFT.LoadFromArray( CurrentArr, fNoIntegrationVariables); 
802    CurrentFT.SetCurveLength( xcurrent );          805    CurrentFT.SetCurveLength( xcurrent );
803                                                   806 
804    PrintStatus(StartFT, CurrentFT, requestStep    807    PrintStatus(StartFT, CurrentFT, requestStep, subStepNo ); 
805 }                                                 808 }
806                                                   809 
807 // -------------------------------------------    810 // ---------------------------------------------------------------------------
808                                                   811 
809 void G4MagInt_Driver::PrintStatus(const G4Fiel    812 void G4MagInt_Driver::PrintStatus(const G4FieldTrack& StartFT,
810                                   const G4Fiel    813                                   const G4FieldTrack& CurrentFT, 
811                                         G4doub    814                                         G4double      requestStep, 
812                                         G4int     815                                         G4int         subStepNo)
813 {                                                 816 {
814     G4int verboseLevel= fVerboseLevel;            817     G4int verboseLevel= fVerboseLevel;
815     const G4int noPrecision = 5;                  818     const G4int noPrecision = 5;
816     G4long oldPrec= G4cout.precision(noPrecisi << 819     G4int oldPrec= G4cout.precision(noPrecision);
817     // G4cout.setf(ios_base::fixed,ios_base::f    820     // G4cout.setf(ios_base::fixed,ios_base::floatfield);
818                                                   821 
819     const G4ThreeVector StartPosition=       S    822     const G4ThreeVector StartPosition=       StartFT.GetPosition();
820     const G4ThreeVector StartUnitVelocity=   S    823     const G4ThreeVector StartUnitVelocity=   StartFT.GetMomentumDir();
821     const G4ThreeVector CurrentPosition=     C    824     const G4ThreeVector CurrentPosition=     CurrentFT.GetPosition();
822     const G4ThreeVector CurrentUnitVelocity= C    825     const G4ThreeVector CurrentUnitVelocity= CurrentFT.GetMomentumDir();
823                                                   826 
824     G4double  DotStartCurrentVeloc= StartUnitV    827     G4double  DotStartCurrentVeloc= StartUnitVelocity.dot(CurrentUnitVelocity);
825                                                   828 
826     G4double step_len= CurrentFT.GetCurveLengt    829     G4double step_len= CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
827     G4double subStepSize = step_len;              830     G4double subStepSize = step_len;
828                                                   831      
829     if( (subStepNo <= 1) || (verboseLevel > 3)    832     if( (subStepNo <= 1) || (verboseLevel > 3) )
830     {                                             833     {
831        subStepNo = - subStepNo;        // To a    834        subStepNo = - subStepNo;        // To allow printing banner
832                                                   835 
833        G4cout << std::setw( 6)  << " " << std:    836        G4cout << std::setw( 6)  << " " << std::setw( 25)
834               << " G4MagInt_Driver: Current Po    837               << " G4MagInt_Driver: Current Position  and  Direction" << " "
835               << G4endl;                          838               << G4endl; 
836        G4cout << std::setw( 5) << "Step#" << "    839        G4cout << std::setw( 5) << "Step#" << " "
837               << std::setw( 7) << "s-curve" <<    840               << std::setw( 7) << "s-curve" << " "
838               << std::setw( 9) << "X(mm)" << "    841               << std::setw( 9) << "X(mm)" << " "
839               << std::setw( 9) << "Y(mm)" << "    842               << std::setw( 9) << "Y(mm)" << " "  
840               << std::setw( 9) << "Z(mm)" << "    843               << std::setw( 9) << "Z(mm)" << " "
841               << std::setw( 8) << " N_x " << "    844               << std::setw( 8) << " N_x " << " "
842               << std::setw( 8) << " N_y " << "    845               << std::setw( 8) << " N_y " << " "
843               << std::setw( 8) << " N_z " << "    846               << std::setw( 8) << " N_z " << " "
844               << std::setw( 8) << " N^2-1 " <<    847               << std::setw( 8) << " N^2-1 " << " "
845               << std::setw(10) << " N(0).N " <    848               << std::setw(10) << " N(0).N " << " "
846               << std::setw( 7) << "KinEner " <    849               << std::setw( 7) << "KinEner " << " "
847               << std::setw(12) << "Track-l" <<    850               << std::setw(12) << "Track-l" << " "   // Add the Sub-step ??
848               << std::setw(12) << "Step-len" <    851               << std::setw(12) << "Step-len" << " " 
849               << std::setw(12) << "Step-len" <    852               << std::setw(12) << "Step-len" << " " 
850               << std::setw( 9) << "ReqStep" <<    853               << std::setw( 9) << "ReqStep" << " "  
851               << G4endl;                          854               << G4endl;
852     }                                             855     }
853                                                   856 
854     if( (subStepNo <= 0) )                        857     if( (subStepNo <= 0) )
855     {                                             858     {
856       PrintStat_Aux( StartFT,  requestStep, 0.    859       PrintStat_Aux( StartFT,  requestStep, 0., 
857                        0,        0.0,             860                        0,        0.0,         1.0);
858     }                                             861     }
859                                                   862 
860     if( verboseLevel <= 3 )                       863     if( verboseLevel <= 3 )
861     {                                             864     {
862       G4cout.precision(noPrecision);              865       G4cout.precision(noPrecision);
863       PrintStat_Aux( CurrentFT, requestStep, s    866       PrintStat_Aux( CurrentFT, requestStep, step_len, 
864                      subStepNo, subStepSize, D    867                      subStepNo, subStepSize, DotStartCurrentVeloc );
865     }                                             868     }
866                                                   869 
867     G4cout.precision(oldPrec);                    870     G4cout.precision(oldPrec);
868 }                                                 871 }
869                                                   872 
870 // -------------------------------------------    873 // ---------------------------------------------------------------------------
871                                                   874 
872 void G4MagInt_Driver::PrintStat_Aux(const G4Fi    875 void G4MagInt_Driver::PrintStat_Aux(const G4FieldTrack& aFieldTrack,
873                                           G4do    876                                           G4double      requestStep, 
874                                           G4do    877                                           G4double      step_len,
875                                           G4in    878                                           G4int         subStepNo,
876                                           G4do    879                                           G4double      subStepSize,
877                                           G4do    880                                           G4double      dotVeloc_StartCurr)
878 {                                                 881 {
879     const G4ThreeVector Position = aFieldTrack    882     const G4ThreeVector Position = aFieldTrack.GetPosition();
880     const G4ThreeVector UnitVelocity = aFieldT    883     const G4ThreeVector UnitVelocity = aFieldTrack.GetMomentumDir();
881                                                   884  
882     if( subStepNo >= 0)                           885     if( subStepNo >= 0)
883     {                                             886     {
884        G4cout << std::setw( 5) << subStepNo <<    887        G4cout << std::setw( 5) << subStepNo << " ";
885     }                                             888     }
886     else                                          889     else
887     {                                             890     {
888        G4cout << std::setw( 5) << "Start" << "    891        G4cout << std::setw( 5) << "Start" << " ";
889     }                                             892     }
890     G4double curveLen= aFieldTrack.GetCurveLen    893     G4double curveLen= aFieldTrack.GetCurveLength();
891     G4cout << std::setw( 7) << curveLen;          894     G4cout << std::setw( 7) << curveLen;
892     G4cout << std::setw( 9) << Position.x() <<    895     G4cout << std::setw( 9) << Position.x() << " "
893            << std::setw( 9) << Position.y() <<    896            << std::setw( 9) << Position.y() << " "
894            << std::setw( 9) << Position.z() <<    897            << std::setw( 9) << Position.z() << " "
895            << std::setw( 8) << UnitVelocity.x(    898            << std::setw( 8) << UnitVelocity.x() << " "
896            << std::setw( 8) << UnitVelocity.y(    899            << std::setw( 8) << UnitVelocity.y() << " "
897            << std::setw( 8) << UnitVelocity.z(    900            << std::setw( 8) << UnitVelocity.z() << " ";
898     G4long oldprec= G4cout.precision(3);       << 901     G4int oldprec= G4cout.precision(3);
899     G4cout << std::setw( 8) << UnitVelocity.ma    902     G4cout << std::setw( 8) << UnitVelocity.mag2()-1.0 << " ";
900     G4cout.precision(6);                          903     G4cout.precision(6);
901     G4cout << std::setw(10) << dotVeloc_StartC    904     G4cout << std::setw(10) << dotVeloc_StartCurr << " ";
902     G4cout.precision(oldprec);                    905     G4cout.precision(oldprec);
903     G4cout << std::setw( 7) << aFieldTrack.Get    906     G4cout << std::setw( 7) << aFieldTrack.GetKineticEnergy();
904     G4cout << std::setw(12) << step_len << " "    907     G4cout << std::setw(12) << step_len << " ";
905                                                   908 
906     static G4ThreadLocal G4double oldCurveLeng    909     static G4ThreadLocal G4double oldCurveLength = 0.0;
907     static G4ThreadLocal G4double oldSubStepLe    910     static G4ThreadLocal G4double oldSubStepLength = 0.0;
908     static G4ThreadLocal G4int oldSubStepNo =     911     static G4ThreadLocal G4int oldSubStepNo = -1;
909                                                   912 
910     G4double subStep_len = 0.0;                   913     G4double subStep_len = 0.0;
911     if( curveLen > oldCurveLength )               914     if( curveLen > oldCurveLength )
912     {                                             915     {
913       subStep_len= curveLen - oldCurveLength;     916       subStep_len= curveLen - oldCurveLength;
914     }                                             917     }
915     else if (subStepNo == oldSubStepNo)           918     else if (subStepNo == oldSubStepNo)
916     {                                             919     {
917       subStep_len= oldSubStepLength;              920       subStep_len= oldSubStepLength;
918     }                                             921     }
919     oldCurveLength= curveLen;                     922     oldCurveLength= curveLen;
920     oldSubStepLength= subStep_len;                923     oldSubStepLength= subStep_len;
921                                                   924 
922     G4cout << std::setw(12) << subStep_len <<     925     G4cout << std::setw(12) << subStep_len << " "; 
923     G4cout << std::setw(12) << subStepSize <<     926     G4cout << std::setw(12) << subStepSize << " "; 
924     if( requestStep != -1.0 )                     927     if( requestStep != -1.0 )
925     {                                             928     {
926       G4cout << std::setw( 9) << requestStep <    929       G4cout << std::setw( 9) << requestStep << " ";
927     }                                             930     }
928     else                                          931     else
929     {                                             932     {
930        G4cout << std::setw( 9) << " InitialSte    933        G4cout << std::setw( 9) << " InitialStep " << " ";
931     }                                             934     }
932     G4cout << G4endl;                             935     G4cout << G4endl;
933 }                                                 936 }
934                                                   937 
935 // -------------------------------------------    938 // ---------------------------------------------------------------------------
936                                                   939 
937 void G4MagInt_Driver::PrintStatisticsReport()     940 void G4MagInt_Driver::PrintStatisticsReport()
938 {                                                 941 {
939   G4int noPrecBig = 6;                            942   G4int noPrecBig = 6;
940   G4long oldPrec = G4cout.precision(noPrecBig) << 943   G4int oldPrec = G4cout.precision(noPrecBig);
941                                                   944 
942   G4cout << "G4MagInt_Driver Statistics of ste    945   G4cout << "G4MagInt_Driver Statistics of steps undertaken. " << G4endl;
943   G4cout << "G4MagInt_Driver: Number of Steps:    946   G4cout << "G4MagInt_Driver: Number of Steps: "
944          << " Total= " <<  fNoTotalSteps          947          << " Total= " <<  fNoTotalSteps
945          << " Bad= "   <<  fNoBadSteps            948          << " Bad= "   <<  fNoBadSteps 
946          << " Small= " <<  fNoSmallSteps          949          << " Small= " <<  fNoSmallSteps 
947          << " Non-initial small= " << (fNoSmal    950          << " Non-initial small= " << (fNoSmallSteps-fNoInitialSmallSteps)
948          << G4endl;                               951          << G4endl;
949  G4cout.precision(oldPrec);                       952  G4cout.precision(oldPrec);
950 }                                                 953 }
951                                                   954  
952 // -------------------------------------------    955 // ---------------------------------------------------------------------------
953                                                   956 
954 void G4MagInt_Driver::SetSmallestFraction(G4do    957 void G4MagInt_Driver::SetSmallestFraction(G4double newFraction)
955 {                                                 958 {
956   if( (newFraction > 1.e-16) && (newFraction <    959   if( (newFraction > 1.e-16) && (newFraction < 1e-8) )
957   {                                               960   {
958     fSmallestFraction= newFraction;               961     fSmallestFraction= newFraction;
959   }                                               962   }
960   else                                            963   else
961   {                                               964   { 
962     std::ostringstream message;                   965     std::ostringstream message;
963     message << "Smallest Fraction not changed.    966     message << "Smallest Fraction not changed. " << G4endl
964             << "  Proposed value was " << newF    967             << "  Proposed value was " << newFraction << G4endl
965             << "  Value must be between 1.e-8     968             << "  Value must be between 1.e-8 and 1.e-16";
966     G4Exception("G4MagInt_Driver::SetSmallestF    969     G4Exception("G4MagInt_Driver::SetSmallestFraction()",
967                 "GeomField1001", JustWarning,     970                 "GeomField1001", JustWarning, message);
968   }                                               971   }
969 }                                                 972 }
970                                                   973 
971 void G4MagInt_Driver::                            974 void G4MagInt_Driver::
972 GetDerivatives(const G4FieldTrack& y_curr, G4d    975 GetDerivatives(const G4FieldTrack& y_curr, G4double* dydx) const
973 {                                                 976 {
974     G4double ytemp[G4FieldTrack::ncompSVEC];      977     G4double ytemp[G4FieldTrack::ncompSVEC];
975     y_curr.DumpToArray(ytemp);                    978     y_curr.DumpToArray(ytemp);
976     pIntStepper->RightHandSide(ytemp, dydx);      979     pIntStepper->RightHandSide(ytemp, dydx);
977       // Avoid virtual call for GetStepper        980       // Avoid virtual call for GetStepper
978       // Was: GetStepper()->ComputeRightHandSi    981       // Was: GetStepper()->ComputeRightHandSide(ytemp, dydx);
979 }                                                 982 }
980                                                   983 
981 void G4MagInt_Driver::GetDerivatives(const G4F    984 void G4MagInt_Driver::GetDerivatives(const G4FieldTrack& track,
982                                      G4double     985                                      G4double dydx[],
983                                      G4double     986                                      G4double field[]) const
984 {                                                 987 {
985     G4double ytemp[G4FieldTrack::ncompSVEC];      988     G4double ytemp[G4FieldTrack::ncompSVEC];
986     track.DumpToArray(ytemp);                     989     track.DumpToArray(ytemp);
987     pIntStepper->RightHandSide(ytemp, dydx, fi    990     pIntStepper->RightHandSide(ytemp, dydx, field);
988 }                                                 991 }
989                                                   992 
990 G4EquationOfMotion* G4MagInt_Driver::GetEquati    993 G4EquationOfMotion* G4MagInt_Driver::GetEquationOfMotion()
991 {                                                 994 {
992     return pIntStepper->GetEquationOfMotion();    995     return pIntStepper->GetEquationOfMotion();
993 }                                                 996 }
994                                                   997 
995 void G4MagInt_Driver::SetEquationOfMotion(G4Eq    998 void G4MagInt_Driver::SetEquationOfMotion(G4EquationOfMotion *equation)
996 {                                                 999 {
997     pIntStepper->SetEquationOfMotion(equation)    1000     pIntStepper->SetEquationOfMotion(equation);
998 }                                                 1001 }
999                                                   1002 
1000 const G4MagIntegratorStepper* G4MagInt_Driver    1003 const G4MagIntegratorStepper* G4MagInt_Driver::GetStepper() const
1001 {                                                1004 {
1002     return pIntStepper;                          1005     return pIntStepper;
1003 }                                                1006 }
1004                                                  1007 
1005 G4MagIntegratorStepper* G4MagInt_Driver::GetS    1008 G4MagIntegratorStepper* G4MagInt_Driver::GetStepper()
1006 {                                                1009 {
1007     return pIntStepper;                          1010     return pIntStepper;
1008 }                                                1011 }
1009                                                  1012 
1010 void G4MagInt_Driver::                           1013 void G4MagInt_Driver::
1011 RenewStepperAndAdjust(G4MagIntegratorStepper*    1014 RenewStepperAndAdjust(G4MagIntegratorStepper* pItsStepper)
1012 {                                                1015 {  
1013     pIntStepper = pItsStepper;                << 1016     pIntStepper = pItsStepper; 
1014     ReSetParameters();                           1017     ReSetParameters();
1015 }                                             << 
1016                                               << 
1017 void G4MagInt_Driver::StreamInfo( std::ostrea << 
1018 {                                             << 
1019     os << "State of G4MagInt_Driver: " << std << 
1020     os << "  Max number of Steps = " << fMaxN << 
1021        << "    (base # = " << fMaxStepBase << << 
1022     os << "  Safety factor       = " << safet << 
1023     os << "  Power - shrink      = " << pshrn << 
1024     os << "  Power - grow        = " << pgrow << 
1025     os << "  threshold (errcon)  = " << errco << 
1026                                               << 
1027     os << "    fMinimumStep =      " << fMini << 
1028     os << "    Smallest Fraction = " << fSmal << 
1029                                               << 
1030     os << "    No Integrat Vars  = " << fNoIn << 
1031     os << "    Min No Vars       = " << fMinN << 
1032     os << "    Num-Vars          = " << fNoVa << 
1033                                               << 
1034     os << "    verbose level     = " << fVerb << 
1035     os << "    Reintegrates      = " << DoesR << 
1036 }                                             << 
1037                                               << 
1038 void PrintInfo( const G4MagInt_Driver & magDr << 
1039 {                                             << 
1040     os << "State of G4MagInt_Driver: " << std << 
1041     os << "  Max number of Steps = " << magDr << 
1042     //   << "    (base # = " << magDrv.fMaxSt << 
1043     os << "  Safety factor       = " << magDr << 
1044     os << "  Power - shrink      = " << magDr << 
1045     os << "  Power - grow        = " << magDr << 
1046     os << "  threshold (errcon)  = " << magDr << 
1047                                               << 
1048     os << "    fMinimumStep =      " << magDr << 
1049     os << "    Smallest Fraction = " << magDr << 
1050                                               << 
1051     /*****                                    << 
1052     os << "    No Integrat Vars  = " << magDr << 
1053     os << "    Min No Vars       = " << magDr << 
1054     os << "    Num-Vars          = " << magDr << 
1055     *****/                                    << 
1056     os << "    verbose level     = " << magDr << 
1057     os << "    Reintegrates      = " << magDr << 
1058 }                                                1018 }
1059                                                  1019