Loading [MathJax]/extensions/tex2jax.js
OGS
TimeLoop.cpp
Go to the documentation of this file.
1
11#include "TimeLoop.h"
12
13#include <range/v3/algorithm/any_of.hpp>
14#include <range/v3/algorithm/contains.hpp>
15
16#include "BaseLib/Error.h"
17#include "BaseLib/RunTime.h"
24#include "ProcessData.h"
25
26namespace
27{
29 std::vector<std::unique_ptr<ProcessLib::ProcessData>> const&
30 per_process_data,
31 double const t)
32{
33 for (auto& process_data : per_process_data)
34 {
35 process_data->process.updateDeactivatedSubdomains(
36 t, process_data->process_id);
37 }
38}
39
40bool isOutputStep(std::vector<ProcessLib::Output> const& outputs,
41 const int timestep, const NumLib::Time& t,
42 const NumLib::Time& end_time)
43{
44 if (end_time == t)
45 {
46 // the last timestep is an output step
47 return true;
48 }
49
50 return ranges::any_of(outputs, [timestep, t](auto const& output)
51 { return output.isOutputStep(timestep, t); });
52}
53
55 int const timestep, NumLib::Time const& t, double const dt,
56 const NumLib::Time& end_time,
57 std::vector<std::unique_ptr<ProcessLib::ProcessData>> const&
58 per_process_data,
59 std::vector<GlobalVector*> const& process_solutions,
60 std::vector<GlobalVector*> const& process_solutions_prev,
61 std::vector<ProcessLib::Output> const& outputs)
62{
63 if (!isOutputStep(outputs, timestep, t, end_time))
64 {
65 return;
66 }
67
68 for (auto& process_data : per_process_data)
69 {
70 auto const process_id = process_data->process_id;
71 auto& pcs = process_data->process;
72
73 pcs.preOutput(t(), dt, process_solutions, process_solutions_prev,
74 process_id);
75 }
76}
77} // namespace
78
79namespace ProcessLib
80{
82 NumLib::Time const& t, double const dt,
83 std::vector<std::unique_ptr<ProcessData>> const& per_process_data,
84 std::vector<GlobalVector*> const& _process_solutions)
85{
86 for (auto& process_data : per_process_data)
87 {
88 auto const process_id = process_data->process_id;
89 auto& pcs = process_data->process;
90 pcs.preTimestep(_process_solutions, t(), dt, process_id);
91 }
92}
93
95 NumLib::Time const& t, double const dt,
96 std::vector<std::unique_ptr<ProcessData>> const& per_process_data,
97 std::vector<GlobalVector*> const& process_solutions,
98 std::vector<GlobalVector*> const& process_solutions_prev)
99{
100 for (auto& process_data : per_process_data)
101 {
102 auto const process_id = process_data->process_id;
103 auto& pcs = process_data->process;
104
105 pcs.computeSecondaryVariable(t(), dt, process_solutions,
106 *process_solutions_prev[process_id],
107 process_id);
108 pcs.postTimestep(process_solutions, process_solutions_prev, t(), dt,
109 process_id);
110 }
111}
112
113template <NumLib::ODESystemTag ODETag>
115 ProcessData& process_data,
117{
118 using Tag = NumLib::NonlinearSolverTag;
119 // A concrete Picard solver
120 using NonlinearSolverPicard = NumLib::NonlinearSolver<Tag::Picard>;
121 // A concrete Newton solver
122 using NonlinearSolverNewton = NumLib::NonlinearSolver<Tag::Newton>;
123
124 if (dynamic_cast<NonlinearSolverPicard*>(&process_data.nonlinear_solver))
125 {
126 // The Picard solver can also work with a Newton-ready ODE,
127 // because the Newton ODESystem derives from the Picard ODESystem.
128 // So no further checks are needed here.
129
130 process_data.tdisc_ode_sys = std::make_unique<
132 process_data.process_id, ode_sys, *process_data.time_disc);
133 }
134 // TODO (naumov) Provide a function to nonlinear_solver to distinguish the
135 // types. Could be handy, because a nonlinear solver could handle both types
136 // like PETScSNES.
137 else if ((dynamic_cast<NonlinearSolverNewton*>(
138 &process_data.nonlinear_solver) != nullptr)
139#ifdef USE_PETSC
140 || (dynamic_cast<NumLib::PETScNonlinearSolver*>(
141 &process_data.nonlinear_solver) != nullptr)
142#endif // USE_PETSC
143 )
144 {
145 // The Newton-Raphson method needs a Newton-ready ODE.
146
148 if (auto* ode_newton = dynamic_cast<ODENewton*>(&ode_sys))
149 {
150 process_data.tdisc_ode_sys = std::make_unique<
152 process_data.process_id, *ode_newton, *process_data.time_disc);
153 }
154 else
155 {
156 OGS_FATAL(
157 "You are trying to solve a non-Newton-ready ODE with the"
158 " Newton-Raphson method. Aborting");
159 }
160 }
161 else
162 {
163 OGS_FATAL("Encountered unknown nonlinear solver type. Aborting");
164 }
165}
166
168{
169 setTimeDiscretizedODESystem(process_data, process_data.process);
170}
171
172std::pair<std::vector<GlobalVector*>, std::vector<GlobalVector*>>
174 NumLib::Time const& t0,
175 std::vector<std::unique_ptr<ProcessData>> const& per_process_data)
176{
177 std::vector<GlobalVector*> process_solutions;
178 std::vector<GlobalVector*> process_solutions_prev;
179
180 for (auto const& process_data : per_process_data)
181 {
182 auto const process_id = process_data->process_id;
183 auto& ode_sys = *process_data->tdisc_ode_sys;
184
185 // append a solution vector of suitable size
186 process_solutions.emplace_back(
188 ode_sys.getMatrixSpecifications(process_id)));
189 process_solutions_prev.emplace_back(
191 ode_sys.getMatrixSpecifications(process_id)));
192 }
193
194 for (auto const& process_data : per_process_data)
195 {
196 auto& pcs = process_data->process;
197 auto const process_id = process_data->process_id;
198 pcs.setInitialConditions(process_solutions, process_solutions_prev,
199 t0(), process_id);
200
201 auto& time_disc = *process_data->time_disc;
202 time_disc.setInitialState(t0()); // push IC
203 }
204
205 return {process_solutions, process_solutions_prev};
206}
207
209 std::vector<std::unique_ptr<ProcessData>> const& per_process_data,
210 std::vector<GlobalVector*> const& process_solutions,
211 std::vector<GlobalVector*> const& process_solutions_prev)
212{
213 for (auto const& process_data : per_process_data)
214 {
215 auto& nonlinear_solver = process_data->nonlinear_solver;
216
217 setEquationSystem(*process_data);
218 nonlinear_solver.calculateNonEquilibriumInitialResiduum(
219 process_solutions, process_solutions_prev,
220 process_data->process_id);
221 }
222}
223
225 std::vector<GlobalVector*>& x, std::vector<GlobalVector*> const& x_prev,
226 std::size_t const timestep, double const t, double const delta_t,
227 ProcessData const& process_data, std::vector<Output> const& outputs)
228{
229 auto& process = process_data.process;
230 int const process_id = process_data.process_id;
231 auto& time_disc = *process_data.time_disc;
232 auto& nonlinear_solver = process_data.nonlinear_solver;
233
234 setEquationSystem(process_data);
235
236 // Note: Order matters!
237 // First advance to the next timestep, then set known solutions at that
238 // time, afterwards pass the right solution vector and time to the
239 // preTimestep() hook.
240
241 time_disc.nextTimestep(t, delta_t);
242
243 auto const post_iteration_callback =
244 [&](int const iteration, bool const converged,
245 std::vector<GlobalVector*> const& x)
246 {
247 // Note: We don't call the postNonLinearSolver(), preOutput(),
248 // computeSecondaryVariable() and postTimestep() hooks here. This might
249 // lead to some inconsistencies in the data compared to regular output.
250 for (auto const& output : outputs)
251 {
252 output.doOutputNonlinearIteration(process, process_id, timestep,
253 NumLib::Time(t), iteration,
254 converged, x);
255 }
256 };
257
258 auto const nonlinear_solver_status =
259 nonlinear_solver.solve(x, x_prev, post_iteration_callback, process_id);
260
261 if (!nonlinear_solver_status.error_norms_met)
262 {
263 return nonlinear_solver_status;
264 }
265
266 process.postNonLinearSolver(x, x_prev, t, delta_t, process_id);
267
268 return nonlinear_solver_status;
269}
270
272 std::vector<Output>&& outputs,
273 std::vector<std::unique_ptr<ProcessData>>&& per_process_data,
274 std::unique_ptr<NumLib::StaggeredCoupling>&& staggered_coupling,
275 const NumLib::Time& start_time, const NumLib::Time& end_time)
276 : _outputs{std::move(outputs)},
277 _per_process_data(std::move(per_process_data)),
278 _start_time(start_time),
279 _end_time(end_time),
280 _staggered_coupling(std::move(staggered_coupling))
281{
282}
283
285 NumLib::TimeStepAlgorithm const& timestep_algorithm,
286 NumLib::Time const& time)
287{
288 // for the first time step we can't compute the changes to the previous
289 // time step
290 if (time == timestep_algorithm.begin())
291 {
292 return false;
293 }
294 return timestep_algorithm.isSolutionErrorComputationNeeded();
295}
296
297std::pair<NumLib::TimeIncrement, bool> TimeLoop::computeTimeStepping(
298 const double prev_dt, NumLib::Time& t, std::size_t& accepted_steps,
299 std::size_t& rejected_steps,
300 std::vector<TimeStepConstraintCallback> const& time_step_constraints)
301{
302 bool all_process_steps_accepted = true;
303 // Get minimum time step size among step sizes of all processes.
304 NumLib::TimeIncrement dt{std::numeric_limits<double>::max()};
305 constexpr double eps = std::numeric_limits<double>::epsilon();
306
307 bool const is_initial_step =
308 std::any_of(_per_process_data.begin(), _per_process_data.end(),
309 [](auto const& ppd) -> bool
310 { return ppd->timestep_current.timeStepNumber() == 0; });
311
312 for (std::size_t i = 0; i < _per_process_data.size(); i++)
313 {
314 auto& ppd = *_per_process_data[i];
315 auto& timestep_algorithm = *ppd.timestep_algorithm.get();
316
317 auto const& x = *_process_solutions[i];
318 auto const& x_prev = *_process_solutions_prev[i];
319
320 const double solution_error =
321 computationOfChangeNeeded(timestep_algorithm, t)
323 x, x_prev,
324 ppd.conv_crit.get() ? ppd.conv_crit->getVectorNormType()
326 : 0.0;
327
328 ppd.timestep_current.setAccepted(
329 ppd.nonlinear_solver_status.error_norms_met);
330
331 auto [previous_step_accepted, timestepper_dt] = timestep_algorithm.next(
332 solution_error, ppd.nonlinear_solver_status.number_iterations,
333 ppd.timestep_previous, ppd.timestep_current);
334
335 if (!previous_step_accepted)
336 {
337 // Not all processes have accepted steps.
338 all_process_steps_accepted = false;
339 }
340
341 if (!ppd.nonlinear_solver_status.error_norms_met)
342 {
343 WARN(
344 "Time step will be rejected due to nonlinear solver "
345 "divergence.");
346 all_process_steps_accepted = false;
347 }
348
349 if (timestepper_dt > eps || t < timestep_algorithm.end())
350 {
351 dt = NumLib::TimeIncrement{std::min(timestepper_dt, dt())};
352 }
353 }
354
355 if (all_process_steps_accepted)
356 {
358 }
359 else
360 {
362 }
363
364 bool last_step_rejected = false;
365 if (!is_initial_step)
366 {
367 if (all_process_steps_accepted)
368 {
369 accepted_steps++;
370 last_step_rejected = false;
371 }
372 else
373 {
374 if (t <= _end_time)
375 {
376 t -= prev_dt;
377 rejected_steps++;
378 last_step_rejected = true;
379 }
380 }
381 }
382
383 // adjust step size considering external communciation_point_calculators
384 for (auto const& time_step_constraint : time_step_constraints)
385 {
387 std::min(dt(), time_step_constraint(t, dt()))};
388 }
389
390 // Check whether the time stepping is stabilized
391 if (std::abs(dt() - prev_dt) < eps)
392 {
393 if (last_step_rejected)
394 {
395 OGS_FATAL(
396 "The new step size of {} is the same as that of the previous "
397 "rejected time step. \nPlease re-run ogs with a proper "
398 "adjustment in the numerical settings, \ne.g. those for time "
399 "stepper, local or global non-linear solver.",
400 dt);
401 }
402 else
403 {
404 DBUG("The time stepping is stabilized with the step size of {}.",
405 dt);
406 }
407 }
408
409 // Reset the time step with the minimum step size, dt
410 // Update the solution of the previous time step.
411 for (std::size_t i = 0; i < _per_process_data.size(); i++)
412 {
413 if (all_process_steps_accepted)
414 {
415 auto& ppd = *_per_process_data[i];
416 NumLib::updateTimeSteps(dt(), ppd.timestep_previous,
417 ppd.timestep_current);
418 auto& timestep_algorithm = ppd.timestep_algorithm;
419 timestep_algorithm->resetCurrentTimeStep(
420 dt(), ppd.timestep_previous, ppd.timestep_current);
421 }
422
423 auto& x = *_process_solutions[i];
424 auto& x_prev = *_process_solutions_prev[i];
425 if (all_process_steps_accepted)
426 {
427 MathLib::LinAlg::copy(x, x_prev); // pushState
428 }
429 else
430 {
431 if (t <= _end_time)
432 {
433 WARN(
434 "Time step {:d} was rejected {:d} times and it will be "
435 "repeated with a reduced step size.",
436 accepted_steps + 1, _repeating_times_of_rejected_step);
437 MathLib::LinAlg::copy(x_prev, x); // popState
438 }
439 }
440 }
441
442 return {dt, last_step_rejected};
443}
444
445std::vector<TimeLoop::TimeStepConstraintCallback>
447 std::vector<double>&& fixed_times) const
448{
449 std::vector<TimeStepConstraintCallback> const time_step_constraints{
450 [fixed_times = std::move(fixed_times)](NumLib::Time const& t, double dt)
451 { return NumLib::possiblyClampDtToNextFixedTime(t, dt, fixed_times); },
452 [this](NumLib::Time const& t, double dt) -> double
453 {
454 if (t < _end_time && _end_time < t + dt)
455 {
456 return _end_time() - t();
457 }
458 return dt;
459 }};
460 return time_step_constraints;
461}
462
465{
466 for (auto const& process_data : _per_process_data)
467 {
468 auto& pcs = process_data->process;
469 for (auto& output : _outputs)
470 {
471 output.addProcess(pcs);
472 }
473
474 setTimeDiscretizedODESystem(*process_data);
475
476 if (auto* conv_crit =
478 process_data->conv_crit.get()))
479 {
480 int const process_id = process_data->process_id;
481 conv_crit->setDOFTable(pcs.getDOFTable(process_id), pcs.getMesh());
482 }
483 }
484
485 // initial solution storage
488
490 {
491 _staggered_coupling->initializeCoupledSolutions(_process_solutions);
492 }
493
494 updateDeactivatedSubdomains(_per_process_data, _start_time());
495
496 auto const time_step_constraints = generateOutputTimeStepConstraints(
498
499 std::tie(_dt, _last_step_rejected) =
501 _rejected_steps, time_step_constraints);
502
503 // Output initial conditions
504 {
507 }
508
511}
512
514{
515 BaseLib::RunTime time_timestep;
516 time_timestep.start();
517
518 _current_time += _dt();
519
520 const std::size_t timesteps = _accepted_steps + 1;
521 // TODO(wenqing): , input option for time unit.
522 INFO("=== Time stepping at step #{:d} and time {} with step size {}",
523 timesteps, _current_time, _dt);
524
525 updateDeactivatedSubdomains(_per_process_data, _current_time());
526
529 INFO("[time] Time step #{:d} took {:g} s.", timesteps,
530 time_timestep.elapsed());
532}
533
535{
536 const double prev_dt = _dt();
537 auto const current_time = _current_time;
538
539 const std::size_t timesteps = _accepted_steps + 1;
540
541 auto const time_step_constraints = generateOutputTimeStepConstraints(
543
544 // _last_step_rejected is also checked in computeTimeStepping.
545 std::tie(_dt, _last_step_rejected) =
547 _rejected_steps, time_step_constraints);
548
550 {
551 outputSolutions(timesteps, current_time(), &Output::doOutput);
552 }
553
554 if (current_time == (_current_time + _dt()))
555 {
556 ERR("Time step size of {} is too small.\n"
557 "Time stepping stops at step {:d} and at time of {}.",
558 _dt, timesteps, _current_time);
559 return false;
560 }
561
563 {
564 return false;
565 }
566
567 return true;
568}
569
571{
572 INFO(
573 "The whole computation of the time stepping took {:d} steps, in which\n"
574 "\t the accepted steps are {:d}, and the rejected steps are {:d}.\n",
576
577 // output last time step
579 {
582 }
583}
584
586 std::size_t const timesteps)
587{
589
590 NumLib::NonlinearSolverStatus nonlinear_solver_status;
591
593 {
594 nonlinear_solver_status =
596 }
597 else
598 {
599 nonlinear_solver_status =
600 solveUncoupledEquationSystems(t, dt, timesteps);
601 }
602
603 // Run post time step only if the last iteration was successful.
604 // Otherwise it runs the risks to get the same errors as in the last
605 // iteration, an exception thrown in assembly, for example.
606 if (nonlinear_solver_status.error_norms_met)
607 {
608 // Later on, the timestep_algorithm might reject the timestep. We assume
609 // that this is a rare case, so still, we call preOutput() here. We
610 // don't expect a large overhead from it.
611 preOutputForAllProcesses(timesteps, t, dt, _end_time, _per_process_data,
613 _outputs);
614
618 }
619 return nonlinear_solver_status.error_norms_met;
620}
621
623 const NumLib::Time& t, const double dt, const std::size_t timestep_id,
624 ProcessData const& process_data, std::vector<GlobalVector*>& x,
625 std::vector<GlobalVector*> const& x_prev,
626 std::vector<Output> const& outputs)
627{
628 BaseLib::RunTime time_timestep_process;
629 time_timestep_process.start();
630
631 auto const nonlinear_solver_status = solveOneTimeStepOneProcess(
632 x, x_prev, timestep_id, t(), dt, process_data, outputs);
633
634 INFO("[time] Solving process #{:d} took {:g} s in time step #{:d}",
635 process_data.process_id, time_timestep_process.elapsed(), timestep_id);
636
637 return nonlinear_solver_status;
638}
639
640static constexpr std::string_view timestepper_cannot_reduce_dt =
641 "Time stepper cannot reduce the time step size further.";
642
644 const NumLib::Time& t, const double dt, const std::size_t timestep_id)
645{
646 NumLib::NonlinearSolverStatus nonlinear_solver_status;
647
648 for (auto const& process_data : _per_process_data)
649 {
650 auto const process_id = process_data->process_id;
651 nonlinear_solver_status = solveMonolithicProcess(
652 t, dt, timestep_id, *process_data, _process_solutions,
654
655 process_data->nonlinear_solver_status = nonlinear_solver_status;
656 if (!nonlinear_solver_status.error_norms_met)
657 {
658 ERR("The nonlinear solver failed in time step #{:d} at t = {} s "
659 "for process #{:d}.",
660 timestep_id, t, process_id);
661
662 if (!process_data->timestep_algorithm->canReduceTimestepSize(
663 process_data->timestep_current,
664 process_data->timestep_previous))
665 {
666 // save unsuccessful solution
667 for (auto const& output : _outputs)
668 {
669 output.doOutputAlways(
670 process_data->process, process_id, timestep_id, t,
671 process_data->nonlinear_solver_status.number_iterations,
672 process_data->nonlinear_solver_status.error_norms_met,
674 }
676 }
677
678 return nonlinear_solver_status;
679 }
680 }
681
682 return nonlinear_solver_status;
683}
684
687 const NumLib::Time& t, const double dt, const std::size_t timestep_id)
688{
689 auto const nonlinear_solver_status =
691 t(), dt, timestep_id, _process_solutions, _process_solutions_prev,
693
694 _last_step_rejected = nonlinear_solver_status.error_norms_met;
695
696 {
697 for (auto const& process_data : _per_process_data)
698 {
699 auto& pcs = process_data->process;
700 int const process_id = process_data->process_id;
701 auto& ode_sys = *process_data->tdisc_ode_sys;
702 pcs.solveReactionEquation(_process_solutions,
703 _process_solutions_prev, t(), dt, ode_sys,
704 process_id);
705 }
706 }
707
708 return nonlinear_solver_status;
709}
710
711template <typename OutputClassMember>
712void TimeLoop::outputSolutions(unsigned timestep, const double t,
713 OutputClassMember output_class_member) const
714{
715 for (auto const& process_data : _per_process_data)
716 {
717 // If nonlinear solver diverged, the solution has already been
718 // saved.
719 if (!process_data->nonlinear_solver_status.error_norms_met)
720 {
721 continue;
722 }
723
724 auto const process_id = process_data->process_id;
725 auto const& pcs = process_data->process;
726
727 for (auto const& output_object : _outputs)
728 {
729 (output_object.*output_class_member)(
730 pcs, process_id, timestep, NumLib::Time(t),
731 process_data->nonlinear_solver_status.number_iterations,
732 process_data->nonlinear_solver_status.error_norms_met,
734 }
735 }
736}
737
749
751 const double dt) const
752{
753 for (auto const& process_data : _per_process_data)
754 {
755 // If nonlinear solver diverged, the solution has already been
756 // saved.
757 if (!process_data->nonlinear_solver_status.error_norms_met)
758 {
759 continue;
760 }
761
762 auto const process_id = process_data->process_id;
763 auto& pcs = process_data->process;
764
765 process_data->time_disc->nextTimestep(t(), dt);
766
767 pcs.preTimestep(_process_solutions, _start_time(), dt, process_id);
768
769 pcs.preOutput(_start_time(), dt, _process_solutions,
770 _process_solutions_prev, process_id);
771
772 // Update secondary variables, which might be uninitialized, before
773 // output.
774 pcs.computeSecondaryVariable(_start_time(), dt, _process_solutions,
775 *_process_solutions_prev[process_id],
776 process_id);
777 }
778}
779} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
Definition of the RunTime class.
Count the running time.
Definition RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32
Interface of time stepping algorithms.
Time begin() const
return the beginning of time steps
virtual bool isSolutionErrorComputationNeeded() const
virtual void releaseVector(GlobalVector const &x)=0
void doOutput(Process const &process, const int process_id, int const timestep, const NumLib::Time &t, int const iteration, bool const converged, std::vector< GlobalVector * > const &xs) const
Definition Output.cpp:338
void doOutputLastTimestep(Process const &process, const int process_id, int const timestep, const NumLib::Time &t, int const iteration, bool const converged, std::vector< GlobalVector * > const &xs) const
Definition Output.cpp:359
NumLib::NonlinearSolverStatus solveUncoupledEquationSystems(const NumLib::Time &t, const double dt, const std::size_t timestep_id)
Member to solver non coupled systems of equations, which can be a single system of equations,...
Definition TimeLoop.cpp:643
void outputLastTimeStep() const
Definition TimeLoop.cpp:570
void preOutputInitialConditions(NumLib::Time const &t, const double dt) const
Definition TimeLoop.cpp:750
NumLib::TimeIncrement _dt
Definition TimeLoop.h:140
NumLib::Time _current_time
Definition TimeLoop.h:137
std::pair< NumLib::TimeIncrement, bool > computeTimeStepping(const double prev_dt, NumLib::Time &t, std::size_t &accepted_steps, std::size_t &rejected_steps, std::vector< TimeStepConstraintCallback > const &time_step_constraints)
Definition TimeLoop.cpp:297
TimeLoop(std::vector< Output > &&outputs, std::vector< std::unique_ptr< ProcessData > > &&per_process_data, std::unique_ptr< NumLib::StaggeredCoupling > &&staggered_coupling, const NumLib::Time &start_time, const NumLib::Time &end_time)
Definition TimeLoop.cpp:271
std::vector< std::unique_ptr< ProcessData > > _per_process_data
Definition TimeLoop.h:133
void outputSolutions(unsigned timestep, const double t, OutputClassMember output_class_member) const
Definition TimeLoop.cpp:712
std::vector< Output > _outputs
Definition TimeLoop.h:132
std::size_t _accepted_steps
Definition TimeLoop.h:138
std::vector< GlobalVector * > _process_solutions
Definition TimeLoop.h:130
std::unique_ptr< NumLib::StaggeredCoupling > _staggered_coupling
Definition TimeLoop.h:144
void initialize()
initialize output, convergence criterion, etc.
Definition TimeLoop.cpp:464
int _repeating_times_of_rejected_step
Definition TimeLoop.h:141
std::size_t _rejected_steps
Definition TimeLoop.h:139
const NumLib::Time _end_time
Definition TimeLoop.h:136
bool preTsNonlinearSolvePostTs(NumLib::Time const &t, double const dt, std::size_t const timesteps)
Definition TimeLoop.cpp:585
NumLib::NonlinearSolverStatus solveCoupledEquationSystemsByStaggeredScheme(const NumLib::Time &t, const double dt, const std::size_t timestep_id)
Member to solver coupled systems of equations by the staggered scheme.
Definition TimeLoop.cpp:686
std::vector< TimeStepConstraintCallback > generateOutputTimeStepConstraints(std::vector< double > &&fixed_times) const
Definition TimeLoop.cpp:446
const NumLib::Time _start_time
Definition TimeLoop.h:135
std::vector< GlobalVector * > _process_solutions_prev
Definition TimeLoop.h:131
NonlinearSolverTag
Tag used to specify which nonlinear solver will be used.
Definition Types.h:20
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
double computeRelativeNorm(VectorType const &x, VectorType const &y, MathLib::VecNormType norm_type)
Definition LinAlg.h:299
void updateTimeSteps(double const dt, TimeStep &previous_timestep, TimeStep &current_timestep)
Definition TimeStep.h:100
double possiblyClampDtToNextFixedTime(Time const &t, double const dt, std::vector< double > const &fixed_output_times)
static constexpr std::string_view timestepper_cannot_reduce_dt
Definition TimeLoop.cpp:640
void setTimeDiscretizedODESystem(ProcessData &process_data, NumLib::ODESystem< ODETag, NumLib::NonlinearSolverTag::Picard > &ode_sys)
Definition TimeLoop.cpp:114
void calculateNonEquilibriumInitialResiduum(std::vector< std::unique_ptr< ProcessData > > const &per_process_data, std::vector< GlobalVector * > const &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev)
Definition TimeLoop.cpp:208
bool computationOfChangeNeeded(NumLib::TimeStepAlgorithm const &timestep_algorithm, NumLib::Time const &time)
Definition TimeLoop.cpp:284
NumLib::NonlinearSolverStatus solveOneTimeStepOneProcess(std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, std::size_t const timestep, double const t, double const delta_t, ProcessData const &process_data, std::vector< Output > const &outputs)
Definition TimeLoop.cpp:224
void preTimestepForAllProcesses(NumLib::Time const &t, double const dt, std::vector< std::unique_ptr< ProcessData > > const &per_process_data, std::vector< GlobalVector * > const &_process_solutions)
Definition TimeLoop.cpp:81
void postTimestepForAllProcesses(NumLib::Time const &t, double const dt, std::vector< std::unique_ptr< ProcessData > > const &per_process_data, std::vector< GlobalVector * > const &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev)
Definition TimeLoop.cpp:94
std::pair< std::vector< GlobalVector * >, std::vector< GlobalVector * > > setInitialConditions(NumLib::Time const &t0, std::vector< std::unique_ptr< ProcessData > > const &per_process_data)
Definition TimeLoop.cpp:173
void setEquationSystem(ProcessData const &process_data)
std::vector< double > calculateUniqueFixedTimesForAllOutputs(std::vector< Output > const &outputs)
Definition Output.cpp:441
static NumLib::NonlinearSolverStatus solveMonolithicProcess(const NumLib::Time &t, const double dt, const std::size_t timestep_id, ProcessData const &process_data, std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, std::vector< Output > const &outputs)
Definition TimeLoop.cpp:622
void preOutputForAllProcesses(int const timestep, NumLib::Time const &t, double const dt, const NumLib::Time &end_time, std::vector< std::unique_ptr< ProcessLib::ProcessData > > const &per_process_data, std::vector< GlobalVector * > const &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, std::vector< ProcessLib::Output > const &outputs)
Definition TimeLoop.cpp:54
bool isOutputStep(std::vector< ProcessLib::Output > const &outputs, const int timestep, const NumLib::Time &t, const NumLib::Time &end_time)
Definition TimeLoop.cpp:40
void updateDeactivatedSubdomains(std::vector< std::unique_ptr< ProcessLib::ProcessData > > const &per_process_data, double const t)
Definition TimeLoop.cpp:28
static NUMLIB_EXPORT VectorProvider & provider
Status of the non-linear solver.
std::unique_ptr< NumLib::TimeDiscretization > time_disc
Definition ProcessData.h:62
NumLib::NonlinearSolverBase & nonlinear_solver
Definition ProcessData.h:58
std::unique_ptr< NumLib::EquationSystem > tdisc_ode_sys
type-erased time-discretized ODE system
Definition ProcessData.h:64