LCOV - code coverage report
Current view: top level - src/timesteppers - IterationAdaptiveDT.C (source / functions) Hit Total Coverage
Test: idaholab/moose framework: 419b9d Lines: 285 315 90.5 %
Date: 2025-08-08 20:01:16 Functions: 14 15 93.3 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : //* This file is part of the MOOSE framework
       2             : //* https://mooseframework.inl.gov
       3             : //*
       4             : //* All rights reserved, see COPYRIGHT for full restrictions
       5             : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
       6             : //*
       7             : //* Licensed under LGPL 2.1, please see LICENSE for details
       8             : //* https://www.gnu.org/licenses/lgpl-2.1.html
       9             : 
      10             : // MOOSE includes
      11             : #include "IterationAdaptiveDT.h"
      12             : #include "Function.h"
      13             : #include "PiecewiseLinear.h"
      14             : #include "Transient.h"
      15             : #include "NonlinearSystem.h"
      16             : #include "FEProblemBase.h"
      17             : #include "LinearSystem.h"
      18             : 
      19             : #include <limits>
      20             : #include <set>
      21             : 
      22             : registerMooseObject("MooseApp", IterationAdaptiveDT);
      23             : 
      24             : InputParameters
      25       15107 : IterationAdaptiveDT::validParams()
      26             : {
      27       15107 :   InputParameters params = TimeStepper::validParams();
      28       15107 :   params.addClassDescription("Adjust the timestep based on the number of iterations");
      29       15107 :   params.addParam<int>(
      30             :       "optimal_iterations",
      31             :       "The target number of solver outer iterations for adaptive timestepping. "
      32             :       "For a problem using nonlinear systems, the total number of nonlinear iterations is used. "
      33             :       "For a problem using linear systems, the total number of linear iterations is used.");
      34       15107 :   params.addParam<int>("iteration_window",
      35             :                        "Attempt to grow/shrink timestep if the iteration count "
      36             :                        "is below/above 'optimal_iterations plus/minus "
      37             :                        "iteration_window' (default = optimal_iterations/5).");
      38       15107 :   params.addParam<unsigned>("linear_iteration_ratio",
      39             :                             "The ratio of linear to nonlinear iterations "
      40             :                             "to determine target linear iterations and "
      41             :                             "window for adaptive timestepping (default = "
      42             :                             "25)");
      43       15107 :   params.addParam<std::vector<PostprocessorName>>("timestep_limiting_postprocessor",
      44             :                                                   "If specified, a list of postprocessor values "
      45             :                                                   "used as an upper limit for the "
      46             :                                                   "current time step length");
      47       15107 :   params.addParam<std::vector<FunctionName>>(
      48             :       "timestep_limiting_function",
      49             :       "A list of 'PiecewiseBase' type functions used to control the timestep by "
      50             :       "limiting the change in the function over a timestep");
      51       15107 :   params.addParam<Real>(
      52             :       "max_function_change",
      53             :       "The absolute value of the maximum change in timestep_limiting_function over a timestep");
      54       45321 :   params.addParam<bool>("force_step_every_function_point",
      55       30214 :                         false,
      56             :                         "Forces the timestepper to take "
      57             :                         "a step that is consistent with "
      58             :                         "points defined in the function");
      59       15107 :   params.addRangeCheckedParam<Real>(
      60             :       "post_function_sync_dt",
      61             :       "post_function_sync_dt>0",
      62             :       "Timestep to apply after time sync with function point. To be used in "
      63             :       "conjunction with 'force_step_every_function_point'.");
      64       15107 :   params.addRequiredParam<Real>("dt", "The default timestep size between solves");
      65       15107 :   params.addParam<std::vector<Real>>("time_t", {}, "The values of t");
      66       15107 :   params.addParam<std::vector<Real>>("time_dt", {}, "The values of dt");
      67       45321 :   params.addParam<Real>("growth_factor",
      68       30214 :                         2.0,
      69             :                         "Factor to apply to timestep if easy convergence (if "
      70             :                         "'optimal_iterations' is specified) or if recovering "
      71             :                         "from failed solve");
      72       45321 :   params.addParam<Real>("cutback_factor",
      73       30214 :                         0.5,
      74             :                         "Factor to apply to timestep if difficult convergence "
      75             :                         "occurs (if 'optimal_iterations' is specified). "
      76             :                         "For failed solves, use cutback_factor_at_failure");
      77             : 
      78       45321 :   params.addParam<bool>("reject_large_step",
      79       30214 :                         false,
      80             :                         "If 'true', time steps that are too large compared to the "
      81             :                         "ideal time step will be rejected and repeated");
      82       45321 :   params.addRangeCheckedParam<Real>("reject_large_step_threshold",
      83       30214 :                                     0.1,
      84             :                                     "reject_large_step_threshold > 0 "
      85             :                                     "& reject_large_step_threshold < 1",
      86             :                                     "Ratio between the the ideal time step size and the "
      87             :                                     "current time step size below which a time step will "
      88             :                                     "be rejected if 'reject_large_step' is 'true'");
      89             : 
      90       15107 :   params.declareControllable("growth_factor cutback_factor");
      91             : 
      92       15107 :   return params;
      93           0 : }
      94             : 
      95         425 : IterationAdaptiveDT::IterationAdaptiveDT(const InputParameters & parameters)
      96             :   : TimeStepper(parameters),
      97             :     PostprocessorInterface(this),
      98         425 :     _dt_old(declareRestartableData<Real>("dt_old", 0.0)),
      99         425 :     _input_dt(getParam<Real>("dt")),
     100         425 :     _tfunc_last_step(declareRestartableData<bool>("tfunc_last_step", false)),
     101         425 :     _sync_last_step(declareRestartableData<bool>("sync_last_step", false)),
     102         850 :     _linear_iteration_ratio(isParamValid("linear_iteration_ratio")
     103         425 :                                 ? getParam<unsigned>("linear_iteration_ratio")
     104             :                                 : 25), // Default to 25
     105         425 :     _adaptive_timestepping(false),
     106         850 :     _pps_value(
     107         850 :         parameters.get<std::vector<PostprocessorName>>("timestep_limiting_postprocessor").size()),
     108         425 :     _timestep_limiting_functions(),
     109         425 :     _piecewise_timestep_limiting_functions(),
     110         425 :     _piecewise_linear_timestep_limiting_functions(),
     111         425 :     _times(0),
     112         425 :     _max_function_change(-1),
     113         425 :     _force_step_every_function_point(getParam<bool>("force_step_every_function_point")),
     114         850 :     _post_function_sync_dt(isParamValid("force_step_every_function_point") &&
     115         850 :                                    isParamValid("post_function_sync_dt")
     116        1275 :                                ? getParam<Real>("post_function_sync_dt")
     117             :                                : 0.0),
     118         425 :     _tfunc_times(getParam<std::vector<Real>>("time_t").begin(),
     119         850 :                  getParam<std::vector<Real>>("time_t").end()),
     120         425 :     _time_ipol(getParam<std::vector<Real>>("time_t"), getParam<std::vector<Real>>("time_dt")),
     121         425 :     _use_time_ipol(_time_ipol.getSampleSize() > 0),
     122         425 :     _growth_factor(getParam<Real>("growth_factor")),
     123         425 :     _cutback_factor(getParam<Real>("cutback_factor")),
     124         425 :     _outer_its(declareRestartableData<unsigned int>("outer_its", 0)),
     125         425 :     _nl_its(_outer_its),
     126         425 :     _l_its(declareRestartableData<unsigned int>("l_its", 0)),
     127         425 :     _cutback_occurred(declareRestartableData<bool>("cutback_occurred", false)),
     128         425 :     _at_function_point(false),
     129         425 :     _reject_large_step(getParam<bool>("reject_large_step")),
     130         850 :     _large_step_rejection_threshold(getParam<Real>("reject_large_step_threshold"))
     131             : {
     132             :   auto timestep_limiting_postprocessor_names =
     133         425 :       parameters.get<std::vector<PostprocessorName>>("timestep_limiting_postprocessor");
     134         481 :   for (size_t i = 0; i < _pps_value.size(); ++i)
     135          56 :     _pps_value[i] = &getPostprocessorValueByName(timestep_limiting_postprocessor_names[i]);
     136             : 
     137         425 :   if (isParamValid("optimal_iterations"))
     138             :   {
     139         324 :     _adaptive_timestepping = true;
     140         324 :     _optimal_iterations = getParam<int>("optimal_iterations");
     141             : 
     142         324 :     if (isParamValid("iteration_window"))
     143           0 :       _iteration_window = getParam<int>("iteration_window");
     144             :     else
     145         324 :       _iteration_window = ceil(_optimal_iterations / 5.0);
     146             :   }
     147             :   else
     148             :   {
     149         101 :     if (isParamValid("iteration_window"))
     150           0 :       mooseError("'optimal_iterations' must be used for 'iteration_window' to be used");
     151         101 :     if (isParamValid("linear_iteration_ratio"))
     152           0 :       mooseError("'optimal_iterations' must be used for 'linear_iteration_ratio' to be used");
     153             :   }
     154             : 
     155         425 :   if (isParamValid("timestep_limiting_function"))
     156          91 :     _max_function_change =
     157          91 :         isParamValid("max_function_change") ? getParam<Real>("max_function_change") : -1;
     158             :   else
     159             :   {
     160         334 :     if (isParamValid("max_function_change"))
     161           4 :       mooseError("'timestep_limiting_function' must be used for 'max_function_change' to be used");
     162         330 :     if (_force_step_every_function_point)
     163           4 :       mooseError("'timestep_limiting_function' must be used for 'force_step_every_function_point' "
     164             :                  "to be used");
     165             :   }
     166             : 
     167         417 :   if (!isParamValid("force_step_every_function_point") && isParamValid("post_function_sync_dt"))
     168           0 :     paramError("post_function_sync_dt",
     169             :                "Not applicable if 'force_step_every_function_point = false'");
     170         417 : }
     171             : 
     172             : void
     173         417 : IterationAdaptiveDT::init()
     174             : {
     175         417 :   if (isParamValid("timestep_limiting_function"))
     176             :   {
     177          91 :     std::set<Real> times;
     178             : 
     179          91 :     const auto tid = isParamValid("_tid") ? getParam<THREAD_ID>("_tid") : 0;
     180         286 :     for (const auto & name : getParam<std::vector<FunctionName>>("timestep_limiting_function"))
     181             :     {
     182         195 :       const auto * func = &_fe_problem.getFunction(name, tid);
     183         195 :       _timestep_limiting_functions.push_back(func);
     184             : 
     185         195 :       const auto * pfunc = dynamic_cast<const PiecewiseBase *>(func);
     186         195 :       _piecewise_timestep_limiting_functions.push_back(pfunc);
     187             : 
     188         195 :       if (pfunc)
     189             :       {
     190         195 :         const auto * plfunc = dynamic_cast<const PiecewiseLinear *>(pfunc);
     191         195 :         _piecewise_linear_timestep_limiting_functions.push_back(plfunc);
     192             : 
     193         195 :         const auto ntimes = pfunc->functionSize();
     194         832 :         for (unsigned int i = 0; i < ntimes; ++i)
     195         637 :           times.insert(pfunc->domain(i));
     196             :       }
     197             :       else
     198           0 :         mooseError("timestep_limiting_function must be a PiecewiseBase function");
     199             :     }
     200          91 :     _times.resize(times.size());
     201          91 :     std::copy(times.begin(), times.end(), _times.begin());
     202             : 
     203             :     mooseAssert(_timestep_limiting_functions.size() ==
     204             :                     _piecewise_timestep_limiting_functions.size(),
     205             :                 "Timestep limiting function count inconsistency");
     206             :     mooseAssert(_piecewise_timestep_limiting_functions.size() ==
     207             :                     _piecewise_linear_timestep_limiting_functions.size(),
     208             :                 "Timestep limiting function count inconsistency");
     209          91 :   }
     210         417 : }
     211             : 
     212             : void
     213         417 : IterationAdaptiveDT::preExecute()
     214             : {
     215         417 :   TimeStepper::preExecute();
     216             : 
     217             :   // Delete all tfunc times that are at or before the begin time
     218         466 :   while (!_tfunc_times.empty() && _time + _timestep_tolerance >= *_tfunc_times.begin())
     219          49 :     _tfunc_times.erase(_tfunc_times.begin());
     220         417 : }
     221             : 
     222             : Real
     223         840 : IterationAdaptiveDT::computeInitialDT()
     224             : {
     225             :   Real dt;
     226         840 :   if (_tfunc_last_step)
     227             :   {
     228          22 :     dt = _time_ipol.sample(_time_old);
     229          22 :     if (_verbose)
     230          22 :       _console << "Setting initial dt to value specified by dt function: " << std::setw(9) << dt
     231          22 :                << std::endl;
     232             :   }
     233             :   else
     234             :   {
     235         818 :     dt = _input_dt;
     236         818 :     if (_verbose)
     237         216 :       _console << "Setting initial dt to input value: " << std::setw(9) << dt << std::endl;
     238             :   }
     239         840 :   return dt;
     240             : }
     241             : 
     242             : Real
     243        3285 : IterationAdaptiveDT::computeDT()
     244             : {
     245        3285 :   Real dt = _dt_old;
     246             : 
     247        3285 :   if (_cutback_occurred)
     248             :   {
     249          36 :     _cutback_occurred = false;
     250          36 :     if (_adaptive_timestepping)
     251             :     {
     252             :       // Don't allow it to grow this step, but shrink if needed
     253           0 :       bool allowToGrow = false;
     254           0 :       computeAdaptiveDT(dt, allowToGrow);
     255             :     }
     256             :   }
     257        3249 :   else if (_tfunc_last_step)
     258             :   {
     259          12 :     _sync_last_step = false;
     260          12 :     dt = _time_ipol.sample(_time_old);
     261             : 
     262          12 :     if (_verbose)
     263          12 :       _console << "Setting dt to value specified by dt function: " << std::setw(9) << dt
     264          12 :                << std::endl;
     265             :   }
     266        3237 :   else if (_sync_last_step)
     267             :   {
     268         300 :     _sync_last_step = false;
     269         300 :     if (_post_function_sync_dt)
     270             :     {
     271         108 :       dt = _post_function_sync_dt;
     272             : 
     273         108 :       if (_verbose)
     274         108 :         _console << "Setting dt to 'post_function_sync_dt': " << std::setw(9) << dt << std::endl;
     275             :     }
     276             :     else
     277             :     {
     278         192 :       dt = _executioner.unconstrainedDT();
     279             : 
     280         192 :       if (_verbose)
     281         192 :         _console << "Setting dt to unconstrained value used before sync: " << std::setw(9) << dt
     282         192 :                  << std::endl;
     283             :     }
     284             :   }
     285        2937 :   else if (_adaptive_timestepping)
     286        2403 :     computeAdaptiveDT(dt);
     287         534 :   else if (_use_time_ipol)
     288           0 :     dt = computeInterpolationDT();
     289             :   else
     290             :   {
     291         534 :     dt *= _growth_factor;
     292         534 :     if (dt > _dt_old * _growth_factor)
     293           0 :       dt = _dt_old * _growth_factor;
     294         534 :     if (_verbose)
     295         336 :       _console << "Growing dt based on growth factor (" << _growth_factor
     296         336 :                << ") and previous dt before sync (" << _dt_old << ") : " << std::setw(9) << dt
     297         336 :                << std::endl;
     298             :   }
     299             : 
     300        3285 :   return dt;
     301             : }
     302             : 
     303             : bool
     304        3399 : IterationAdaptiveDT::constrainStep(Real & dt)
     305             : {
     306        3399 :   bool at_sync_point = TimeStepper::constrainStep(dt);
     307             : 
     308             :   // Use value from computed dt while rejecting the timestep
     309        3399 :   if (_dt_from_reject)
     310             :   {
     311          12 :     dt = *_dt_from_reject;
     312          12 :     _dt_from_reject.reset();
     313             :   }
     314             :   // Otherwise, limit the timestep to the current postprocessor value
     315             :   else
     316        3387 :     limitDTToPostprocessorValue(dt);
     317             : 
     318             :   // Limit the timestep to limit change in the function
     319        3395 :   limitDTByFunction(dt);
     320             : 
     321             :   // Adjust to the next tfunc time if needed
     322        3395 :   if (!_tfunc_times.empty() && _time + dt + _timestep_tolerance >= *_tfunc_times.begin())
     323             :   {
     324          24 :     dt = *_tfunc_times.begin() - _time;
     325             : 
     326          24 :     if (_verbose)
     327          12 :       _console << "Limiting dt to sync with dt function time: " << std::setw(9)
     328          12 :                << *_tfunc_times.begin() << " dt: " << std::setw(9) << dt << std::endl;
     329             :   }
     330             : 
     331        3395 :   return at_sync_point;
     332             : }
     333             : 
     334             : Real
     335          84 : IterationAdaptiveDT::computeFailedDT()
     336             : {
     337          84 :   _cutback_occurred = true;
     338             : 
     339             :   // Can't cut back any more
     340          84 :   if (_dt <= _dt_min)
     341           0 :     mooseError("Solve failed and timestep already at dtmin, cannot continue!");
     342             : 
     343          84 :   if (_verbose)
     344             :   {
     345           0 :     _console << "\nSolve failed with dt: " << std::setw(9) << _dt
     346           0 :              << "\nRetrying with reduced dt: " << std::setw(9) << _dt * _cutback_factor_at_failure
     347           0 :              << std::endl;
     348             :   }
     349             :   else
     350          84 :     _console << "\nSolve failed, cutting timestep." << std::endl;
     351             : 
     352          84 :   return _dt * _cutback_factor_at_failure;
     353             : }
     354             : 
     355             : bool
     356        9575 : IterationAdaptiveDT::converged() const
     357             : {
     358        9575 :   if (!_reject_large_step)
     359        9394 :     return TimeStepper::converged();
     360             : 
     361             :   // the solver has not converged
     362         181 :   if (!TimeStepper::converged())
     363           0 :     return false;
     364             : 
     365             :   // we are already at dt_min or at the start of the simulation
     366             :   // in which case we can move on to the next step
     367         181 :   if (_dt == _dt_min || _t_step < 2)
     368          36 :     return true;
     369             : 
     370             :   // This means we haven't tried constraining the latest step yet
     371         145 :   if (_dt_from_reject)
     372          12 :     return false;
     373             : 
     374             :   // we get what the next time step should be
     375         133 :   Real dt_test = _dt;
     376         133 :   limitDTToPostprocessorValue(dt_test);
     377             : 
     378             :   // we cannot constrain the time step any further
     379         133 :   if (dt_test == 0)
     380           0 :     return true;
     381             : 
     382             :   // if the time step is much smaller than the current time step
     383             :   // we need to repeat the current iteration with a smaller time step
     384         133 :   if (dt_test < _dt * _large_step_rejection_threshold)
     385             :   {
     386          12 :     _dt_from_reject = dt_test;
     387          12 :     return false;
     388             :   }
     389             : 
     390             :   // otherwise we move one
     391         121 :   return true;
     392             : }
     393             : 
     394             : void
     395        3520 : IterationAdaptiveDT::limitDTToPostprocessorValue(Real & limitedDT) const
     396             : {
     397        3520 :   if (_pps_value.size() != 0 && _t_step > 1)
     398             :   {
     399         521 :     Real limiting_pps_value = *_pps_value[0];
     400         521 :     unsigned int i_min = 0;
     401         701 :     for (size_t i = 1; i < _pps_value.size(); ++i)
     402         180 :       if (*_pps_value[i] < limiting_pps_value)
     403             :       {
     404         120 :         limiting_pps_value = *_pps_value[i];
     405         120 :         i_min = i;
     406             :       }
     407             : 
     408         521 :     if (limitedDT > limiting_pps_value)
     409             :     {
     410         244 :       if (limiting_pps_value < 0)
     411           4 :         mooseWarning(
     412           4 :             "Negative timestep limiting postprocessor '" +
     413           8 :             getParam<std::vector<PostprocessorName>>("timestep_limiting_postprocessor")[i_min] +
     414           8 :             "': " + std::to_string(limiting_pps_value));
     415         240 :       limitedDT = std::max(_dt_min, limiting_pps_value);
     416             : 
     417         240 :       if (_verbose)
     418           0 :         _console << "Limiting dt to postprocessor value. dt = " << limitedDT << std::endl;
     419             :     }
     420             :   }
     421        3516 : }
     422             : 
     423             : void
     424        3395 : IterationAdaptiveDT::limitDTByFunction(Real & limitedDT)
     425             : {
     426        3395 :   Real orig_dt = limitedDT;
     427        3395 :   const auto nfunc = _timestep_limiting_functions.size();
     428        3395 :   Real restricted_step = std::numeric_limits<Real>::max();
     429             : 
     430        6203 :   for (unsigned int j = 0; j < nfunc; ++j)
     431             :   {
     432             :     // Limit by function change for piecewise linear functions.
     433        2808 :     if (_piecewise_linear_timestep_limiting_functions[j] && _max_function_change > 0)
     434             :     {
     435             :       const auto current_function_value =
     436         372 :           _piecewise_linear_timestep_limiting_functions[j]->value(_time, {});
     437             : 
     438         372 :       const auto ntimes = _piecewise_linear_timestep_limiting_functions[j]->functionSize();
     439         732 :       for (std::size_t next_time_index = 1; next_time_index < ntimes; ++next_time_index)
     440             :       {
     441             :         const auto next_time =
     442         672 :             _piecewise_linear_timestep_limiting_functions[j]->domain(next_time_index);
     443             : 
     444             :         // Skip ahead to find time point that is just past the current time.
     445         672 :         if (next_time + _timestep_tolerance <= _time)
     446         204 :           continue;
     447             : 
     448             :         // Find out how far we can go without exceeding the max function change.
     449             :         const auto next_function_value =
     450         468 :             _piecewise_linear_timestep_limiting_functions[j]->range(next_time_index);
     451         468 :         const auto change = std::abs(next_function_value - current_function_value);
     452         468 :         if (change > _max_function_change)
     453             :         {
     454             :           // Interpolate to find step.
     455         192 :           restricted_step =
     456         192 :               std::min(restricted_step, (_max_function_change / change) * (next_time - _time));
     457         192 :           break;
     458             :         }
     459             : 
     460             :         // Don't keep going if we've already passed the current limited step.
     461         276 :         if (next_time > _time + limitedDT)
     462         120 :           break;
     463             :       }
     464             :     }
     465        2436 :     else if (_timestep_limiting_functions[j])
     466             :     {
     467        2436 :       const Real old_value = _timestep_limiting_functions[j]->value(_time_old);
     468        2436 :       Real new_value = _timestep_limiting_functions[j]->value(_time_old + limitedDT);
     469        2436 :       Real change = std::abs(new_value - old_value);
     470             : 
     471        2436 :       if (_max_function_change > 0.0 && change > _max_function_change)
     472             :         do
     473             :         {
     474           0 :           limitedDT /= 2.0;
     475           0 :           new_value = _timestep_limiting_functions[j]->value(_time_old + limitedDT);
     476           0 :           change = std::abs(new_value - old_value);
     477           0 :         } while (change > _max_function_change);
     478             :     }
     479             :   }
     480             : 
     481        3395 :   if (restricted_step < limitedDT)
     482         120 :     limitedDT = std::max(_dt_min, restricted_step);
     483             : 
     484        3395 :   _at_function_point = false;
     485        3395 :   if (_force_step_every_function_point)
     486        3612 :     for (unsigned int i = 0; i + 1 < _times.size(); ++i)
     487        3552 :       if (_time >= _times[i] && _time < _times[i + 1])
     488             :       {
     489         648 :         if (limitedDT > _times[i + 1] - _time - _timestep_tolerance)
     490             :         {
     491         348 :           limitedDT = _times[i + 1] - _time;
     492         348 :           _at_function_point = true;
     493             :         }
     494         648 :         break;
     495             :       }
     496             : 
     497        3395 :   if (_verbose && limitedDT != orig_dt)
     498             :   {
     499         432 :     if (_at_function_point)
     500         312 :       _console << "Limiting dt to match function point. dt = ";
     501             :     else
     502         120 :       _console << "Limiting dt to limit change in function. dt = ";
     503             : 
     504         432 :     _console << limitedDT << std::endl;
     505             :   }
     506        3395 : }
     507             : 
     508             : void
     509        2403 : IterationAdaptiveDT::computeAdaptiveDT(Real & dt, bool allowToGrow, bool allowToShrink)
     510             : {
     511        2403 :   const unsigned int growth_outer_its(
     512        2403 :       _optimal_iterations > _iteration_window ? _optimal_iterations - _iteration_window : 0);
     513        2403 :   const unsigned int shrink_outer_its(_optimal_iterations + _iteration_window);
     514        4806 :   const unsigned int growth_l_its(_optimal_iterations > _iteration_window
     515        2403 :                                       ? _linear_iteration_ratio *
     516        2199 :                                             (_optimal_iterations - _iteration_window)
     517             :                                       : 0);
     518        2403 :   const unsigned int shrink_l_its(_linear_iteration_ratio *
     519        2403 :                                   (_optimal_iterations + _iteration_window));
     520        2403 :   const std::string ite_type = (_fe_problem.numLinearSystems() == 0) ? "nl" : "solver";
     521             : 
     522        2403 :   if (allowToGrow && (_outer_its < growth_outer_its && _l_its < growth_l_its))
     523             :   {
     524             :     // Grow the timestep
     525        2163 :     dt *= _growth_factor;
     526             : 
     527        2163 :     if (_verbose)
     528         598 :       _console << "Growing dt: " + ite_type + " its = " << _outer_its << " < " << growth_outer_its
     529         299 :                << " && lin its = " << _l_its << " < " << growth_l_its << " old dt: " << std::setw(9)
     530         299 :                << _dt_old << " new dt: " << std::setw(9) << dt << '\n';
     531             :   }
     532         240 :   else if (allowToShrink && (_outer_its > shrink_outer_its || _l_its > shrink_l_its))
     533             :   {
     534             :     // Shrink the timestep
     535         240 :     dt *= _cutback_factor;
     536             : 
     537         240 :     if (_verbose)
     538          72 :       _console << "Shrinking dt: " + ite_type + " its = " << _outer_its << " > " << shrink_outer_its
     539          36 :                << " || lin its = " << _l_its << " > " << shrink_l_its << " old dt: " << std::setw(9)
     540          36 :                << _dt_old << " new dt: " << std::setw(9) << dt << '\n';
     541             :   }
     542             : 
     543        2403 :   _console << std::flush;
     544        2403 : }
     545             : 
     546             : Real
     547           0 : IterationAdaptiveDT::computeInterpolationDT()
     548             : {
     549           0 :   Real dt = _time_ipol.sample(_time_old);
     550             : 
     551           0 :   if (dt > _dt_old * _growth_factor)
     552             :   {
     553           0 :     dt = _dt_old * _growth_factor;
     554             : 
     555           0 :     if (_verbose)
     556           0 :       _console << "Growing dt to recover from cutback. "
     557           0 :                << " old dt: " << std::setw(9) << _dt_old << " new dt: " << std::setw(9) << dt
     558           0 :                << std::endl;
     559             :   }
     560             : 
     561           0 :   return dt;
     562             : }
     563             : 
     564             : void
     565         108 : IterationAdaptiveDT::rejectStep()
     566             : {
     567         108 :   TimeStepper::rejectStep();
     568         108 : }
     569             : 
     570             : void
     571        5363 : IterationAdaptiveDT::acceptStep()
     572             : {
     573        5363 :   TimeStepper::acceptStep();
     574             : 
     575        5363 :   _tfunc_last_step = false;
     576        5375 :   while (!_tfunc_times.empty() && _time + _timestep_tolerance >= *_tfunc_times.begin())
     577             :   {
     578          12 :     if (std::abs(_time - *_tfunc_times.begin()) <= _timestep_tolerance)
     579          12 :       _tfunc_last_step = true;
     580             : 
     581          12 :     _tfunc_times.erase(_tfunc_times.begin());
     582             :   }
     583             : 
     584             :   // Reset counts
     585        5363 :   _outer_its = 0;
     586        5363 :   _l_its = 0;
     587             : 
     588             :   // Use the total number of iterations for multi-system
     589       10630 :   for (const auto i : make_range(_fe_problem.numNonlinearSystems()))
     590        5267 :     _outer_its += _fe_problem.getNonlinearSystemBase(/*nl_sys=*/i).nNonlinearIterations();
     591             :   // Add linear iterations for both nonlinear and linear systems for multi-system
     592       10630 :   for (const auto i : make_range(_fe_problem.numNonlinearSystems()))
     593        5267 :     _l_its += _fe_problem.getNonlinearSystemBase(/*nl_sys=*/i).nLinearIterations();
     594        5459 :   for (const auto i : make_range(_fe_problem.numLinearSystems()))
     595             :   {
     596          96 :     _l_its += _fe_problem.getLinearSystem(/*nl_sys=*/i).nLinearIterations();
     597          96 :     _outer_its += _fe_problem.getLinearSystem(/*nl_sys=*/i).nLinearIterations();
     598             :   }
     599             : 
     600        5711 :   if ((_at_function_point || _executioner.atSyncPoint()) &&
     601         348 :       _dt + _timestep_tolerance < _executioner.unconstrainedDT())
     602             :   {
     603         336 :     _dt_old = _fe_problem.dtOld();
     604         336 :     _sync_last_step = true;
     605             : 
     606         336 :     if (_verbose)
     607         336 :       _console << "Sync point hit in current step, using previous dt for old dt: " << std::setw(9)
     608         336 :                << _dt_old << std::endl;
     609             :   }
     610             :   else
     611        5027 :     _dt_old = _dt;
     612        5363 : }

Generated by: LCOV version 1.14