53 const G4double G4MagInt_Driver::max_stepping_increase = 5.0;
54 const G4double G4MagInt_Driver::max_stepping_decrease = 0.1;
59 const G4int G4MagInt_Driver::fMaxStepBase = 250;
61 #ifndef G4NO_FIELD_STATISTICS
72 G4int statisticsVerbose)
73 : fSmallestFraction( 1.0e-12 ),
74 fNoIntegrationVariables(numComponents),
76 fNoVars( std::
max( fNoIntegrationVariables, fMinNoVars )),
77 fStatisticsVerboseLevel(statisticsVerbose),
78 fNoTotalSteps(0), fNoBadSteps(0), fNoSmallSteps(0),
79 fNoInitialSmallSteps(0),
80 fDyerr_max(0.0), fDyerr_mx2(0.0),
81 fDyerrPos_smTot(0.0), fDyerrPos_lgTot(0.0), fDyerrVel_lgTot(0.0),
82 fSumH_sm(0.0), fSumH_lg(0.0),
89 fMinimumStep= hminimum;
95 if( (fVerboseLevel > 0) || (fStatisticsVerboseLevel > 1) )
97 G4cout <<
"MagIntDriver version: Accur-Adv: "
98 <<
"invE_nS, QuickAdv-2sqrt with Statistics "
114 if( fStatisticsVerboseLevel > 1 )
138 G4int nstp, i, no_warnings=0;
143 static G4int nStpPr=50;
151 G4bool succeeded =
true, lastStepSucceeded;
155 G4int noFullIntegr=0, noSmallIntegr = 0 ;
157 const G4int nvar= fNoVars;
167 std::ostringstream message;
168 message <<
"Proposed step is zero; hstep = " << hstep <<
" !";
175 std::ostringstream message;
176 message <<
"Invalid run condition." <<
G4endl
177 <<
"Proposed step is negative; hstep = " << hstep <<
"." <<
G4endl
178 <<
"Requested step cannot be negative! Aborting event.";
188 x1= startCurveLength;
191 if ( (hinitial > 0.0) && (hinitial < hstep)
203 for (i=0;i<nvar;i++) { y[i] = ystart[i]; }
214 for (i=0;i<nvar;i++) { ySubStepStart[i] = y[i]; }
227 if( h > fMinimumStep )
231 lastStepSucceeded= (hdid == h);
234 PrintStatus( ySubStepStart, xSubStepStart, y, x, h, nstp);
242 G4double dchord_step, dyerr, dyerr_len;
246 QuickAdvance( yFldTrk, dydx, h, dchord_step, dyerr_len );
253 if ( dyerr_len > fDyerr_max) { fDyerr_max= dyerr_len; }
254 fDyerrPos_smTot += dyerr_len;
256 if (nstp<=1) { fNoInitialSmallSteps++; }
261 if(fNoSmallSteps<2) {
PrintStatus(ySubStepStart, x1, y, x, h, -nstp); }
262 G4cout <<
"Another sub-min step, no " << fNoSmallSteps
263 <<
" of " << fNoTotalSteps <<
" this time " << nstp <<
G4endl;
265 G4cout <<
" dyerr= " << dyerr_len <<
" relative = " << dyerr_len / h
266 <<
" epsilon= " << eps <<
" hstep= " << hstep
267 <<
" h= " << h <<
" hmin= " << fMinimumStep <<
G4endl;
274 "Integration Step became Zero!");
276 dyerr = dyerr_len / h;
284 lastStepSucceeded= (dyerr<=
eps);
287 if (lastStepSucceeded) { noFullIntegr++; }
288 else { noSmallIntegr++; }
293 if( (dbg>0) && (dbg<=2) && (nstp>nStpPr))
295 if( nstp==nStpPr ) {
G4cout <<
"***** Many steps ****" <<
G4endl; }
297 G4cout <<
"hdid=" << std::setw(12) << hdid <<
" "
298 <<
"hnext=" << std::setw(12) << hnext <<
" "
299 <<
"hstep=" << std::setw(12) << hstep <<
" (requested) "
301 PrintStatus( ystart, x1, y, x, h, (nstp==nStpPr) ? -nstp: nstp);
306 G4double endPointDist= (EndPos-StartPos).mag();
319 G4cerr <<
" Total steps: bad " << fNoBadSteps
320 <<
" good " << noGoodSteps <<
" current h= " << hdid
322 PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);
335 if( (h < eps * hstep) || (h < fSmallestFraction * startCurveLength) )
343 if(std::fabs(hnext) <=
Hmin())
347 if( (x < x2 * (1-eps) ) &&
348 (std::fabs(hstep) >
Hmin()) )
353 PrintStatus( ystart, x1, y, x, hstep, no_warnings?nstp:-nstp);
380 G4cout <<
"Warning: G4MagIntegratorDriver::AccurateAdvance"
382 <<
" Integration step 'h' became "
383 << h <<
" due to roundoff. " <<
G4endl
384 <<
" Calculated as difference of x2= "<< x2 <<
" and x=" << x
385 <<
" Forcing termination of advance." <<
G4endl;
391 }
while ( ((nstp++)<=fMaxNoSteps) && (x < x2) && (!lastStep) );
399 for (i=0;i<nvar;i++) { yEnd[i] = y[i]; }
405 if(nstp > fMaxNoSteps)
419 if( dbg && no_warnings )
421 G4cerr <<
"G4MagIntegratorDriver exit status: no-steps " << nstp <<
G4endl;
437 const G4int maxNoWarnings = 10;
438 std::ostringstream message;
439 if( (noWarningsIssued < maxNoWarnings) || fVerboseLevel > 10 )
441 message <<
"The stepsize for the next iteration, " << hnext
442 <<
", is too small - in Step number " << nstp <<
"." <<
G4endl
443 <<
"The minimum for the driver is " <<
Hmin() <<
G4endl
444 <<
"Requested integr. length was " << hstep <<
" ." <<
G4endl
445 <<
"The size of this sub-step was " << h <<
" ." <<
G4endl
446 <<
"The integrations has already gone " << xDone;
450 message <<
"Too small 'next' step " << hnext
451 <<
", step-no: " << nstp <<
G4endl
452 <<
", this sub-step: " << h
453 <<
", req_tot_len: " << hstep
454 <<
", done: " << xDone <<
", min: " <<
Hmin();
456 G4Exception(
"G4MagInt_Driver::WarnSmallStepSize()",
"GeomField1001",
468 std::ostringstream message;
469 message <<
"The number of steps used in the Integration driver"
470 <<
" (Runge-Kutta) is too many." <<
G4endl
471 <<
"Integration of the interval was not completed !" <<
G4endl
472 <<
"Only a " << (xCurrent-x1start)*100/(x2end-x1start)
473 <<
" % fraction of it was done.";
474 G4Exception(
"G4MagInt_Driver::WarnTooManyStep()",
"GeomField1001",
487 G4bool isNewMax, prNewMax;
489 isNewMax = endPointDist > (1.0 + maxRelError) * h;
490 prNewMax = endPointDist > (1.0 + 1.05 * maxRelError) * h;
491 if( isNewMax ) { maxRelError= endPointDist / h - 1.0; }
494 && ( (dbg>1) || prNewMax || (endPointDist >= h*(1.+
eps) ) ) )
497 std::ostringstream message;
498 if( (noWarnings ++ < 10) || (dbg>2) )
500 message <<
"The integration produced an end-point which " <<
G4endl
501 <<
"is further from the start-point than the curve length."
504 message <<
" Distance of endpoints = " << endPointDist
505 <<
", curve length = " << h <<
G4endl
506 <<
" Difference (curveLen-endpDist)= " << (h - endPointDist)
507 <<
", relative = " << (h-endPointDist) / h
508 <<
", epsilon = " <<
eps;
509 G4Exception(
"G4MagInt_Driver::WarnEndPointTooFar()",
"GeomField1001",
546 G4double inv_eps_vel_sq = 1.0 / (eps_rel_max*eps_rel_max);
555 const G4int max_trials=100;
559 G4bool hasSpin= (spin_mag2 > 0.0);
561 for (iter=0; iter<max_trials ;iter++)
564 pIntStepper-> Stepper(y,dydx,h,ytemp,yerr);
567 G4double inv_eps_pos_sq = 1.0 / (eps_pos*eps_pos);
571 errpos_sq =
sqr(yerr[0]) +
sqr(yerr[1]) +
sqr(yerr[2]) ;
572 errpos_sq *= inv_eps_pos_sq;
577 if( magvel_sq > 0.0 ) {
578 errvel_sq = sumerr_sq / magvel_sq;
580 G4cerr <<
"** G4MagIntegrationDriver: found case of zero momentum."
581 <<
" iteration= " << iter <<
" h= " << h <<
G4endl;
582 errvel_sq = sumerr_sq;
584 errvel_sq *= inv_eps_vel_sq;
585 errmax_sq =
std::max( errpos_sq, errvel_sq );
590 errspin_sq = (
sqr(yerr[9]) +
sqr(yerr[10]) +
sqr(yerr[11]) )
592 errspin_sq *= inv_eps_vel_sq;
593 errmax_sq =
std::max( errmax_sq, errspin_sq );
596 if ( errmax_sq <= 1.0 ) {
break; }
601 if (htemp >= 0.1*h) { h = htemp; }
607 G4cerr <<
"G4MagIntegratorDriver::OneGoodStep:" <<
G4endl
608 <<
" Stepsize underflow in Stepper " <<
G4endl ;
609 G4cerr <<
" Step's start x=" << x <<
" and end x= " << xnew
610 <<
" are equal !! " << G4endl
611 <<
" Due to step-size= " << h
612 <<
" . Note that input step was " << htry <<
G4endl;
620 fDyerrPos_lgTot += errpos_sq;
621 fDyerrVel_lgTot += errvel_sq * h * h;
625 if (errmax_sq > errcon*errcon)
631 hnext = max_stepping_increase*h ;
635 for(
G4int k=0;k<fNoIntegrationVariables;k++) { y[k] = ytemp[k]; }
652 G4Exception(
"G4MagInt_Driver::QuickAdvance()",
"GeomField0001",
656 dchord_step = dyerr_pos_sq = hstep * hstep * dydx[0];
670 G4double dyerr_pos_sq, dyerr_mom_rel_sq;
674 G4double dyerr_mom_sq, vel_mag_sq, inv_vel_mag_sq;
684 pIntStepper-> Stepper(yarrin, dydx, hstep, yarrout, yerr_vec) ;
688 dchord_step= pIntStepper-> DistChord();
699 PrintStatus( yarrin, s_start, yarrout, s_start+hstep, hstep, 1);
705 vel_mag_sq = (
sqr(yarrout[3])+
sqr(yarrout[4])+
sqr(yarrout[5]) );
706 inv_vel_mag_sq = 1.0 / vel_mag_sq;
707 dyerr_pos_sq = (
sqr(yerr_vec[0])+
sqr(yerr_vec[1])+
sqr(yerr_vec[2]));
708 dyerr_mom_sq = (
sqr(yerr_vec[3])+
sqr(yerr_vec[4])+
sqr(yerr_vec[5]));
709 dyerr_mom_rel_sq = dyerr_mom_sq * inv_vel_mag_sq;
715 #ifdef RETURN_A_NEW_STEP_LENGTH
717 dyerr_len = std::sqrt( dyerr_len_sq );
718 dyerr_len_sq /=
eps ;
727 if( dyerr_pos_sq > ( dyerr_mom_rel_sq *
sqr(hstep) ) )
729 dyerr = std::sqrt(dyerr_pos_sq);
734 dyerr = std::sqrt(dyerr_mom_rel_sq) * hstep;
742 #ifdef QUICK_ADV_ARRAY_IN_AND_OUT
751 G4Exception(
"G4MagInt_Driver::QuickAdvance()",
"GeomField0001",
753 dyerr = dchord_step = hstep * yarrin[0] * dydx[0];
754 yarrout[0]= yarrin[0];
771 if(errMaxNorm > 1.0 )
775 }
else if(errMaxNorm > 0.0 ) {
780 hnew = max_stepping_increase * hstepCurrent;
801 if (errMaxNorm > 1.0 )
806 if (hnew < max_stepping_decrease*hstepCurrent)
808 hnew = max_stepping_decrease*hstepCurrent ;
816 if (errMaxNorm > errcon)
819 { hnew = max_stepping_increase * hstepCurrent; }
843 CurrentFT.
LoadFromArray( CurrentArr, fNoIntegrationVariables);
846 PrintStatus(StartFT, CurrentFT, requestStep, subStepNo );
857 G4int verboseLevel= fVerboseLevel;
867 G4double DotStartCurrentVeloc= StartUnitVelocity.
dot(CurrentUnitVelocity);
872 if( (subStepNo <= 1) || (verboseLevel > 3) )
874 subStepNo = - subStepNo;
876 G4cout << std::setw( 6) <<
" " << std::setw( 25)
877 <<
" G4MagInt_Driver: Current Position and Direction" <<
" "
879 G4cout << std::setw( 5) <<
"Step#" <<
" "
880 << std::setw( 7) <<
"s-curve" <<
" "
881 << std::setw( 9) <<
"X(mm)" <<
" "
882 << std::setw( 9) <<
"Y(mm)" <<
" "
883 << std::setw( 9) <<
"Z(mm)" <<
" "
884 << std::setw( 8) <<
" N_x " <<
" "
885 << std::setw( 8) <<
" N_y " <<
" "
886 << std::setw( 8) <<
" N_z " <<
" "
887 << std::setw( 8) <<
" N^2-1 " <<
" "
888 << std::setw(10) <<
" N(0).N " <<
" "
889 << std::setw( 7) <<
"KinEner " <<
" "
890 << std::setw(12) <<
"Track-l" <<
" "
891 << std::setw(12) <<
"Step-len" <<
" "
892 << std::setw(12) <<
"Step-len" <<
" "
893 << std::setw( 9) <<
"ReqStep" <<
" "
897 if( (subStepNo <= 0) )
904 if( verboseLevel <= 3 )
906 G4cout.precision(noPrecision);
908 subStepNo, subStepSize, DotStartCurrentVeloc );
924 G4cout.precision(oldPrec);
942 G4cout << std::setw( 5) << subStepNo <<
" ";
946 G4cout << std::setw( 5) <<
"Start" <<
" ";
949 G4cout << std::setw( 7) << curveLen;
950 G4cout << std::setw( 9) << Position.
x() <<
" "
951 << std::setw( 9) << Position.
y() <<
" "
952 << std::setw( 9) << Position.
z() <<
" "
953 << std::setw( 8) << UnitVelocity.
x() <<
" "
954 << std::setw( 8) << UnitVelocity.
y() <<
" "
955 << std::setw( 8) << UnitVelocity.
z() <<
" ";
957 G4cout << std::setw( 8) << UnitVelocity.
mag2()-1.0 <<
" ";
959 G4cout << std::setw(10) << dotVeloc_StartCurr <<
" ";
960 G4cout.precision(oldprec);
962 G4cout << std::setw(12) << step_len <<
" ";
969 if( curveLen > oldCurveLength )
971 subStep_len= curveLen - oldCurveLength;
973 else if (subStepNo == oldSubStepNo)
975 subStep_len= oldSubStepLength;
977 oldCurveLength= curveLen;
978 oldSubStepLength= subStep_len;
980 G4cout << std::setw(12) << subStep_len <<
" ";
981 G4cout << std::setw(12) << subStepSize <<
" ";
982 if( requestStep != -1.0 )
984 G4cout << std::setw( 9) << requestStep <<
" ";
988 G4cout << std::setw( 9) <<
" InitialStep " <<
" ";
1000 G4cout <<
"G4MagInt_Driver Statistics of steps undertaken. " <<
G4endl;
1001 G4cout <<
"G4MagInt_Driver: Number of Steps: "
1002 <<
" Total= " << fNoTotalSteps
1003 <<
" Bad= " << fNoBadSteps
1004 <<
" Small= " << fNoSmallSteps
1005 <<
" Non-initial small= " << (fNoSmallSteps-fNoInitialSmallSteps)
1010 <<
" maximum= " << fDyerr_max
1011 <<
" Sum small= " << fDyerrPos_smTot
1012 <<
" std::sqrt(Sum large^2): pos= " << std::sqrt(fDyerrPos_lgTot)
1013 <<
" vel= " << std::sqrt( fDyerrVel_lgTot )
1014 <<
" Total h-distance: small= " << fSumH_sm
1015 <<
" large= " << fSumH_lg
1019 G4int noPrecSmall=4;
1021 G4cout.precision(noPrecSmall);
1022 G4cout <<
"MIDnums: " << fMinimumStep
1023 <<
" " << fNoTotalSteps
1024 <<
" " << fNoSmallSteps
1025 <<
" " << fNoSmallSteps-fNoInitialSmallSteps
1026 <<
" " << fNoBadSteps
1027 <<
" " << fDyerr_max
1028 <<
" " << fDyerr_mx2
1029 <<
" " << fDyerrPos_smTot
1031 <<
" " << fDyerrPos_lgTot
1032 <<
" " << fDyerrVel_lgTot
1038 G4cout.precision(oldPrec);
1045 if( (newFraction > 1.e-16) && (newFraction < 1e-8) )
1047 fSmallestFraction= newFraction;
1051 G4cerr <<
"Warning: SmallestFraction not changed. " <<
G4endl
1052 <<
" Proposed value was " << newFraction <<
G4endl
1053 <<
" Value must be between 1.e-8 and 1.e-16" <<
G4endl;
static constexpr double perMillion
virtual void ComputeRightHandSide(const G4double y[], G4double dydx[])
G4double GetCurveLength() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetKineticEnergy() const
double dot(const Hep3Vector &) const
void SetCurveLength(G4double nCurve_s)
void RenewStepperAndAdjust(G4MagIntegratorStepper *pItsStepper)
void WarnTooManyStep(G4double x1start, G4double x2end, G4double xCurrent)
void WarnEndPointTooFar(G4double endPointDist, G4double hStepSize, G4double epsilonRelative, G4int debugFlag)
G4double GetSurfaceTolerance() const
const G4ThreeVector & GetMomentumDir() const
static const G4double eps
G4MagInt_Driver(G4double hminimum, G4MagIntegratorStepper *pItsStepper, G4int numberOfComponents=6, G4int statisticsVerbosity=1)
void OneGoodStep(G4double ystart[], const G4double dydx[], G4double &x, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
G4double GetSafety() const
void PrintStatisticsReport()
virtual G4int IntegratorOrder() const =0
G4ThreeVector GetPosition() const
G4GLOB_DLL std::ostream G4cout
void PrintStatus(const G4double *StartArr, G4double xstart, const G4double *CurrentArr, G4double xcurrent, G4double requestStep, G4int subStepNo)
G4double ComputeNewStepSize_WithinLimits(G4double errMaxNorm, G4double hstepCurrent)
G4double GetPgrow() const
void PrintStat_Aux(const G4FieldTrack &aFieldTrack, G4double requestStep, G4double actualStep, G4int subStepNo, G4double subStepSize, G4double dotVelocities)
G4double GetPshrnk() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
G4bool QuickAdvance(G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr)
void DumpToArray(G4double valArr[ncompSVEC]) const
void WarnSmallStepSize(G4double hnext, G4double hstep, G4double h, G4double xDone, G4int noSteps)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void LoadFromArray(const G4double valArr[ncompSVEC], G4int noVarsIntegrated)
void SetSmallestFraction(G4double val)
G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent)
static constexpr double perThousand
G4bool AccurateAdvance(G4FieldTrack &y_current, G4double hstep, G4double eps, G4double hinitial=0.0)
static G4GeometryTolerance * GetInstance()
G4GLOB_DLL std::ostream G4cerr