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 15045 : IterationAdaptiveDT::validParams()
26 : {
27 15045 : InputParameters params = TimeStepper::validParams();
28 15045 : params.addClassDescription("Adjust the timestep based on the number of iterations");
29 15045 : 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 15045 : 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 15045 : 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 15045 : 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 15045 : 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 15045 : params.addParam<Real>(
52 : "max_function_change",
53 : "The absolute value of the maximum change in timestep_limiting_function over a timestep");
54 45135 : params.addParam<bool>("force_step_every_function_point",
55 30090 : false,
56 : "Forces the timestepper to take "
57 : "a step that is consistent with "
58 : "points defined in the function");
59 15045 : 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 15045 : params.addRequiredParam<Real>("dt", "The default timestep size between solves");
65 15045 : params.addParam<std::vector<Real>>("time_t", {}, "The values of t");
66 15045 : params.addParam<std::vector<Real>>("time_dt", {}, "The values of dt");
67 45135 : params.addParam<Real>("growth_factor",
68 30090 : 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 45135 : params.addParam<Real>("cutback_factor",
73 30090 : 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 45135 : params.addParam<bool>("reject_large_step",
79 30090 : false,
80 : "If 'true', time steps that are too large compared to the "
81 : "ideal time step will be rejected and repeated");
82 45135 : params.addRangeCheckedParam<Real>("reject_large_step_threshold",
83 30090 : 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 15045 : params.declareControllable("growth_factor cutback_factor");
91 :
92 15045 : return params;
93 0 : }
94 :
95 394 : IterationAdaptiveDT::IterationAdaptiveDT(const InputParameters & parameters)
96 : : TimeStepper(parameters),
97 : PostprocessorInterface(this),
98 394 : _dt_old(declareRestartableData<Real>("dt_old", 0.0)),
99 394 : _input_dt(getParam<Real>("dt")),
100 394 : _tfunc_last_step(declareRestartableData<bool>("tfunc_last_step", false)),
101 394 : _sync_last_step(declareRestartableData<bool>("sync_last_step", false)),
102 788 : _linear_iteration_ratio(isParamValid("linear_iteration_ratio")
103 394 : ? getParam<unsigned>("linear_iteration_ratio")
104 : : 25), // Default to 25
105 394 : _adaptive_timestepping(false),
106 788 : _pps_value(
107 788 : parameters.get<std::vector<PostprocessorName>>("timestep_limiting_postprocessor").size()),
108 394 : _timestep_limiting_functions(),
109 394 : _piecewise_timestep_limiting_functions(),
110 394 : _piecewise_linear_timestep_limiting_functions(),
111 394 : _times(0),
112 394 : _max_function_change(-1),
113 394 : _force_step_every_function_point(getParam<bool>("force_step_every_function_point")),
114 788 : _post_function_sync_dt(isParamValid("force_step_every_function_point") &&
115 788 : isParamValid("post_function_sync_dt")
116 1182 : ? getParam<Real>("post_function_sync_dt")
117 : : 0.0),
118 394 : _tfunc_times(getParam<std::vector<Real>>("time_t").begin(),
119 788 : getParam<std::vector<Real>>("time_t").end()),
120 394 : _time_ipol(getParam<std::vector<Real>>("time_t"), getParam<std::vector<Real>>("time_dt")),
121 394 : _use_time_ipol(_time_ipol.getSampleSize() > 0),
122 394 : _growth_factor(getParam<Real>("growth_factor")),
123 394 : _cutback_factor(getParam<Real>("cutback_factor")),
124 394 : _outer_its(declareRestartableData<unsigned int>("outer_its", 0)),
125 394 : _nl_its(_outer_its),
126 394 : _l_its(declareRestartableData<unsigned int>("l_its", 0)),
127 394 : _cutback_occurred(declareRestartableData<bool>("cutback_occurred", false)),
128 394 : _at_function_point(false),
129 394 : _reject_large_step(getParam<bool>("reject_large_step")),
130 788 : _large_step_rejection_threshold(getParam<Real>("reject_large_step_threshold"))
131 : {
132 : auto timestep_limiting_postprocessor_names =
133 394 : parameters.get<std::vector<PostprocessorName>>("timestep_limiting_postprocessor");
134 446 : for (size_t i = 0; i < _pps_value.size(); ++i)
135 52 : _pps_value[i] = &getPostprocessorValueByName(timestep_limiting_postprocessor_names[i]);
136 :
137 394 : if (isParamValid("optimal_iterations"))
138 : {
139 300 : _adaptive_timestepping = true;
140 300 : _optimal_iterations = getParam<int>("optimal_iterations");
141 :
142 300 : if (isParamValid("iteration_window"))
143 0 : _iteration_window = getParam<int>("iteration_window");
144 : else
145 300 : _iteration_window = ceil(_optimal_iterations / 5.0);
146 : }
147 : else
148 : {
149 94 : if (isParamValid("iteration_window"))
150 0 : mooseError("'optimal_iterations' must be used for 'iteration_window' to be used");
151 94 : if (isParamValid("linear_iteration_ratio"))
152 0 : mooseError("'optimal_iterations' must be used for 'linear_iteration_ratio' to be used");
153 : }
154 :
155 394 : if (isParamValid("timestep_limiting_function"))
156 84 : _max_function_change =
157 84 : isParamValid("max_function_change") ? getParam<Real>("max_function_change") : -1;
158 : else
159 : {
160 310 : if (isParamValid("max_function_change"))
161 4 : mooseError("'timestep_limiting_function' must be used for 'max_function_change' to be used");
162 306 : 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 386 : 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 386 : }
171 :
172 : void
173 386 : IterationAdaptiveDT::init()
174 : {
175 386 : if (isParamValid("timestep_limiting_function"))
176 : {
177 84 : std::set<Real> times;
178 :
179 84 : const auto tid = isParamValid("_tid") ? getParam<THREAD_ID>("_tid") : 0;
180 264 : for (const auto & name : getParam<std::vector<FunctionName>>("timestep_limiting_function"))
181 : {
182 180 : const auto * func = &_fe_problem.getFunction(name, tid);
183 180 : _timestep_limiting_functions.push_back(func);
184 :
185 180 : const auto * pfunc = dynamic_cast<const PiecewiseBase *>(func);
186 180 : _piecewise_timestep_limiting_functions.push_back(pfunc);
187 :
188 180 : if (pfunc)
189 : {
190 180 : const auto * plfunc = dynamic_cast<const PiecewiseLinear *>(pfunc);
191 180 : _piecewise_linear_timestep_limiting_functions.push_back(plfunc);
192 :
193 180 : const auto ntimes = pfunc->functionSize();
194 768 : for (unsigned int i = 0; i < ntimes; ++i)
195 588 : times.insert(pfunc->domain(i));
196 : }
197 : else
198 0 : mooseError("timestep_limiting_function must be a PiecewiseBase function");
199 : }
200 84 : _times.resize(times.size());
201 84 : 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 84 : }
210 386 : }
211 :
212 : void
213 386 : IterationAdaptiveDT::preExecute()
214 : {
215 386 : TimeStepper::preExecute();
216 :
217 : // Delete all tfunc times that are at or before the begin time
218 431 : while (!_tfunc_times.empty() && _time + _timestep_tolerance >= *_tfunc_times.begin())
219 45 : _tfunc_times.erase(_tfunc_times.begin());
220 386 : }
221 :
222 : Real
223 773 : IterationAdaptiveDT::computeInitialDT()
224 : {
225 : Real dt;
226 773 : if (_tfunc_last_step)
227 : {
228 20 : dt = _time_ipol.sample(_time_old);
229 20 : if (_verbose)
230 20 : _console << "Setting initial dt to value specified by dt function: " << std::setw(9) << dt
231 20 : << std::endl;
232 : }
233 : else
234 : {
235 753 : dt = _input_dt;
236 753 : if (_verbose)
237 198 : _console << "Setting initial dt to input value: " << std::setw(9) << dt << std::endl;
238 : }
239 773 : return dt;
240 : }
241 :
242 : Real
243 3011 : IterationAdaptiveDT::computeDT()
244 : {
245 3011 : Real dt = _dt_old;
246 :
247 3011 : if (_cutback_occurred)
248 : {
249 33 : _cutback_occurred = false;
250 33 : 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 2978 : else if (_tfunc_last_step)
258 : {
259 11 : _sync_last_step = false;
260 11 : dt = _time_ipol.sample(_time_old);
261 :
262 11 : if (_verbose)
263 11 : _console << "Setting dt to value specified by dt function: " << std::setw(9) << dt
264 11 : << std::endl;
265 : }
266 2967 : else if (_sync_last_step)
267 : {
268 275 : _sync_last_step = false;
269 275 : if (_post_function_sync_dt)
270 : {
271 99 : dt = _post_function_sync_dt;
272 :
273 99 : if (_verbose)
274 99 : _console << "Setting dt to 'post_function_sync_dt': " << std::setw(9) << dt << std::endl;
275 : }
276 : else
277 : {
278 176 : dt = _executioner.unconstrainedDT();
279 :
280 176 : if (_verbose)
281 176 : _console << "Setting dt to unconstrained value used before sync: " << std::setw(9) << dt
282 176 : << std::endl;
283 : }
284 : }
285 2692 : else if (_adaptive_timestepping)
286 2200 : computeAdaptiveDT(dt);
287 492 : else if (_use_time_ipol)
288 0 : dt = computeInterpolationDT();
289 : else
290 : {
291 492 : dt *= _growth_factor;
292 492 : if (dt > _dt_old * _growth_factor)
293 0 : dt = _dt_old * _growth_factor;
294 492 : if (_verbose)
295 308 : _console << "Growing dt based on growth factor (" << _growth_factor
296 308 : << ") and previous dt before sync (" << _dt_old << ") : " << std::setw(9) << dt
297 308 : << std::endl;
298 : }
299 :
300 3011 : return dt;
301 : }
302 :
303 : bool
304 3117 : IterationAdaptiveDT::constrainStep(Real & dt)
305 : {
306 3117 : bool at_sync_point = TimeStepper::constrainStep(dt);
307 :
308 : // Use value from computed dt while rejecting the timestep
309 3117 : if (_dt_from_reject)
310 : {
311 11 : dt = *_dt_from_reject;
312 11 : _dt_from_reject.reset();
313 : }
314 : // Otherwise, limit the timestep to the current postprocessor value
315 : else
316 3106 : limitDTToPostprocessorValue(dt);
317 :
318 : // Limit the timestep to limit change in the function
319 3113 : limitDTByFunction(dt);
320 :
321 : // Adjust to the next tfunc time if needed
322 3113 : if (!_tfunc_times.empty() && _time + dt + _timestep_tolerance >= *_tfunc_times.begin())
323 : {
324 22 : dt = *_tfunc_times.begin() - _time;
325 :
326 22 : if (_verbose)
327 11 : _console << "Limiting dt to sync with dt function time: " << std::setw(9)
328 11 : << *_tfunc_times.begin() << " dt: " << std::setw(9) << dt << std::endl;
329 : }
330 :
331 3113 : return at_sync_point;
332 : }
333 :
334 : Real
335 77 : IterationAdaptiveDT::computeFailedDT()
336 : {
337 77 : _cutback_occurred = true;
338 :
339 : // Can't cut back any more
340 77 : if (_dt <= _dt_min)
341 0 : mooseError("Solve failed and timestep already at dtmin, cannot continue!");
342 :
343 77 : 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 77 : _console << "\nSolve failed, cutting timestep." << std::endl;
351 :
352 77 : return _dt * _cutback_factor_at_failure;
353 : }
354 :
355 : bool
356 8745 : IterationAdaptiveDT::converged() const
357 : {
358 8745 : if (!_reject_large_step)
359 8580 : return TimeStepper::converged();
360 :
361 : // the solver has not converged
362 165 : 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 165 : if (_dt == _dt_min || _t_step < 2)
368 33 : return true;
369 :
370 : // This means we haven't tried constraining the latest step yet
371 132 : if (_dt_from_reject)
372 11 : return false;
373 :
374 : // we get what the next time step should be
375 121 : Real dt_test = _dt;
376 121 : limitDTToPostprocessorValue(dt_test);
377 :
378 : // we cannot constrain the time step any further
379 121 : 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 121 : if (dt_test < _dt * _large_step_rejection_threshold)
385 : {
386 11 : _dt_from_reject = dt_test;
387 11 : return false;
388 : }
389 :
390 : // otherwise we move one
391 110 : return true;
392 : }
393 :
394 : void
395 3227 : IterationAdaptiveDT::limitDTToPostprocessorValue(Real & limitedDT) const
396 : {
397 3227 : if (_pps_value.size() != 0 && _t_step > 1)
398 : {
399 477 : Real limiting_pps_value = *_pps_value[0];
400 477 : unsigned int i_min = 0;
401 642 : for (size_t i = 1; i < _pps_value.size(); ++i)
402 165 : if (*_pps_value[i] < limiting_pps_value)
403 : {
404 110 : limiting_pps_value = *_pps_value[i];
405 110 : i_min = i;
406 : }
407 :
408 477 : if (limitedDT > limiting_pps_value)
409 : {
410 224 : 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 220 : limitedDT = std::max(_dt_min, limiting_pps_value);
416 :
417 220 : if (_verbose)
418 0 : _console << "Limiting dt to postprocessor value. dt = " << limitedDT << std::endl;
419 : }
420 : }
421 3223 : }
422 :
423 : void
424 3113 : IterationAdaptiveDT::limitDTByFunction(Real & limitedDT)
425 : {
426 3113 : Real orig_dt = limitedDT;
427 3113 : const auto nfunc = _timestep_limiting_functions.size();
428 3113 : Real restricted_step = std::numeric_limits<Real>::max();
429 :
430 5687 : for (unsigned int j = 0; j < nfunc; ++j)
431 : {
432 : // Limit by function change for piecewise linear functions.
433 2574 : if (_piecewise_linear_timestep_limiting_functions[j] && _max_function_change > 0)
434 : {
435 : const auto current_function_value =
436 341 : _piecewise_linear_timestep_limiting_functions[j]->value(_time, {});
437 :
438 341 : const auto ntimes = _piecewise_linear_timestep_limiting_functions[j]->functionSize();
439 671 : for (std::size_t next_time_index = 1; next_time_index < ntimes; ++next_time_index)
440 : {
441 : const auto next_time =
442 616 : _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 616 : if (next_time + _timestep_tolerance <= _time)
446 187 : continue;
447 :
448 : // Find out how far we can go without exceeding the max function change.
449 : const auto next_function_value =
450 429 : _piecewise_linear_timestep_limiting_functions[j]->range(next_time_index);
451 429 : const auto change = std::abs(next_function_value - current_function_value);
452 429 : if (change > _max_function_change)
453 : {
454 : // Interpolate to find step.
455 176 : restricted_step =
456 176 : std::min(restricted_step, (_max_function_change / change) * (next_time - _time));
457 176 : break;
458 : }
459 :
460 : // Don't keep going if we've already passed the current limited step.
461 253 : if (next_time > _time + limitedDT)
462 110 : break;
463 : }
464 : }
465 2233 : else if (_timestep_limiting_functions[j])
466 : {
467 2233 : const Real old_value = _timestep_limiting_functions[j]->value(_time_old);
468 2233 : Real new_value = _timestep_limiting_functions[j]->value(_time_old + limitedDT);
469 2233 : Real change = std::abs(new_value - old_value);
470 :
471 2233 : 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 3113 : if (restricted_step < limitedDT)
482 110 : limitedDT = std::max(_dt_min, restricted_step);
483 :
484 3113 : _at_function_point = false;
485 3113 : if (_force_step_every_function_point)
486 3311 : for (unsigned int i = 0; i + 1 < _times.size(); ++i)
487 3256 : if (_time >= _times[i] && _time < _times[i + 1])
488 : {
489 594 : if (limitedDT > _times[i + 1] - _time - _timestep_tolerance)
490 : {
491 319 : limitedDT = _times[i + 1] - _time;
492 319 : _at_function_point = true;
493 : }
494 594 : break;
495 : }
496 :
497 3113 : if (_verbose && limitedDT != orig_dt)
498 : {
499 396 : if (_at_function_point)
500 286 : _console << "Limiting dt to match function point. dt = ";
501 : else
502 110 : _console << "Limiting dt to limit change in function. dt = ";
503 :
504 396 : _console << limitedDT << std::endl;
505 : }
506 3113 : }
507 :
508 : void
509 2200 : IterationAdaptiveDT::computeAdaptiveDT(Real & dt, bool allowToGrow, bool allowToShrink)
510 : {
511 2200 : const unsigned int growth_outer_its(
512 2200 : _optimal_iterations > _iteration_window ? _optimal_iterations - _iteration_window : 0);
513 2200 : const unsigned int shrink_outer_its(_optimal_iterations + _iteration_window);
514 4400 : const unsigned int growth_l_its(_optimal_iterations > _iteration_window
515 2200 : ? _linear_iteration_ratio *
516 2013 : (_optimal_iterations - _iteration_window)
517 : : 0);
518 2200 : const unsigned int shrink_l_its(_linear_iteration_ratio *
519 2200 : (_optimal_iterations + _iteration_window));
520 2200 : const std::string ite_type = (_fe_problem.numLinearSystems() == 0) ? "nl" : "solver";
521 :
522 2200 : if (allowToGrow && (_outer_its < growth_outer_its && _l_its < growth_l_its))
523 : {
524 : // Grow the timestep
525 1980 : dt *= _growth_factor;
526 :
527 1980 : if (_verbose)
528 548 : _console << "Growing dt: " + ite_type + " its = " << _outer_its << " < " << growth_outer_its
529 274 : << " && lin its = " << _l_its << " < " << growth_l_its << " old dt: " << std::setw(9)
530 274 : << _dt_old << " new dt: " << std::setw(9) << dt << '\n';
531 : }
532 220 : else if (allowToShrink && (_outer_its > shrink_outer_its || _l_its > shrink_l_its))
533 : {
534 : // Shrink the timestep
535 220 : dt *= _cutback_factor;
536 :
537 220 : if (_verbose)
538 66 : _console << "Shrinking dt: " + ite_type + " its = " << _outer_its << " > " << shrink_outer_its
539 33 : << " || lin its = " << _l_its << " > " << shrink_l_its << " old dt: " << std::setw(9)
540 33 : << _dt_old << " new dt: " << std::setw(9) << dt << '\n';
541 : }
542 :
543 2200 : _console << std::flush;
544 2200 : }
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 77 : IterationAdaptiveDT::rejectStep()
566 : {
567 77 : TimeStepper::rejectStep();
568 77 : }
569 :
570 : void
571 4906 : IterationAdaptiveDT::acceptStep()
572 : {
573 4906 : TimeStepper::acceptStep();
574 :
575 4906 : _tfunc_last_step = false;
576 4917 : while (!_tfunc_times.empty() && _time + _timestep_tolerance >= *_tfunc_times.begin())
577 : {
578 11 : if (std::abs(_time - *_tfunc_times.begin()) <= _timestep_tolerance)
579 11 : _tfunc_last_step = true;
580 :
581 11 : _tfunc_times.erase(_tfunc_times.begin());
582 : }
583 :
584 : // Reset counts
585 4906 : _outer_its = 0;
586 4906 : _l_its = 0;
587 :
588 : // Use the total number of iterations for multi-system
589 9724 : for (const auto i : make_range(_fe_problem.numNonlinearSystems()))
590 4818 : _outer_its += _fe_problem.getNonlinearSystemBase(/*nl_sys=*/i).nNonlinearIterations();
591 : // Add linear iterations for both nonlinear and linear systems for multi-system
592 9724 : for (const auto i : make_range(_fe_problem.numNonlinearSystems()))
593 4818 : _l_its += _fe_problem.getNonlinearSystemBase(/*nl_sys=*/i).nLinearIterations();
594 4994 : for (const auto i : make_range(_fe_problem.numLinearSystems()))
595 : {
596 88 : _l_its += _fe_problem.getLinearSystem(/*nl_sys=*/i).nLinearIterations();
597 88 : _outer_its += _fe_problem.getLinearSystem(/*nl_sys=*/i).nLinearIterations();
598 : }
599 :
600 5225 : if ((_at_function_point || _executioner.atSyncPoint()) &&
601 319 : _dt + _timestep_tolerance < _executioner.unconstrainedDT())
602 : {
603 308 : _dt_old = _fe_problem.dtOld();
604 308 : _sync_last_step = true;
605 :
606 308 : if (_verbose)
607 308 : _console << "Sync point hit in current step, using previous dt for old dt: " << std::setw(9)
608 308 : << _dt_old << std::endl;
609 : }
610 : else
611 4598 : _dt_old = _dt;
612 4906 : }
|