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 iteration, std::vector<GlobalVector*> const& x)
245 {
246 // Note: We don't call the postNonLinearSolver(), preOutput(),
247 // computeSecondaryVariable() and postTimestep() hooks here. This might
248 // lead to some inconsistencies in the data compared to regular output.
249 for (auto const& output : outputs)
250 {
251 output.doOutputNonlinearIteration(process, process_id, timestep,
252 NumLib::Time(t), iteration, x);
253 }
254 };
255
256 auto const nonlinear_solver_status =
257 nonlinear_solver.solve(x, x_prev, post_iteration_callback, process_id);
258
259 if (!nonlinear_solver_status.error_norms_met)
260 {
261 return nonlinear_solver_status;
262 }
263
264 process.postNonLinearSolver(x, x_prev, t, delta_t, process_id);
265
266 return nonlinear_solver_status;
267}
268
270 std::vector<Output>&& outputs,
271 std::vector<std::unique_ptr<ProcessData>>&& per_process_data,
272 std::unique_ptr<NumLib::StaggeredCoupling>&& staggered_coupling,
273 const NumLib::Time& start_time, const NumLib::Time& end_time)
274 : _outputs{std::move(outputs)},
275 _per_process_data(std::move(per_process_data)),
276 _start_time(start_time),
277 _end_time(end_time),
278 _staggered_coupling(std::move(staggered_coupling))
279{
280}
281
283 NumLib::TimeStepAlgorithm const& timestep_algorithm,
284 NumLib::Time const& time)
285{
286 // for the first time step we can't compute the changes to the previous
287 // time step
288 if (time == timestep_algorithm.begin())
289 {
290 return false;
291 }
292 return timestep_algorithm.isSolutionErrorComputationNeeded();
293}
294
295std::pair<NumLib::TimeIncrement, bool> TimeLoop::computeTimeStepping(
296 const double prev_dt, NumLib::Time& t, std::size_t& accepted_steps,
297 std::size_t& rejected_steps,
298 std::vector<TimeStepConstraintCallback> const& time_step_constraints)
299{
300 bool all_process_steps_accepted = true;
301 // Get minimum time step size among step sizes of all processes.
302 NumLib::TimeIncrement dt{std::numeric_limits<double>::max()};
303 constexpr double eps = std::numeric_limits<double>::epsilon();
304
305 bool const is_initial_step =
306 std::any_of(_per_process_data.begin(), _per_process_data.end(),
307 [](auto const& ppd) -> bool
308 { return ppd->timestep_current.timeStepNumber() == 0; });
309
310 for (std::size_t i = 0; i < _per_process_data.size(); i++)
311 {
312 auto& ppd = *_per_process_data[i];
313 auto& timestep_algorithm = *ppd.timestep_algorithm.get();
314
315 auto const& x = *_process_solutions[i];
316 auto const& x_prev = *_process_solutions_prev[i];
317
318 const double solution_error =
319 computationOfChangeNeeded(timestep_algorithm, t)
321 x, x_prev,
322 ppd.conv_crit.get() ? ppd.conv_crit->getVectorNormType()
324 : 0.0;
325
326 ppd.timestep_current.setAccepted(
327 ppd.nonlinear_solver_status.error_norms_met);
328
329 auto [previous_step_accepted, timestepper_dt] = timestep_algorithm.next(
330 solution_error, ppd.nonlinear_solver_status.number_iterations,
331 ppd.timestep_previous, ppd.timestep_current);
332
333 if (!previous_step_accepted)
334 {
335 // Not all processes have accepted steps.
336 all_process_steps_accepted = false;
337 }
338
339 if (!ppd.nonlinear_solver_status.error_norms_met)
340 {
341 WARN(
342 "Time step will be rejected due to nonlinear solver "
343 "divergence.");
344 all_process_steps_accepted = false;
345 }
346
347 if (timestepper_dt > eps || t < timestep_algorithm.end())
348 {
349 dt = NumLib::TimeIncrement{std::min(timestepper_dt, dt())};
350 }
351 }
352
353 if (all_process_steps_accepted)
354 {
356 }
357 else
358 {
360 }
361
362 bool last_step_rejected = false;
363 if (!is_initial_step)
364 {
365 if (all_process_steps_accepted)
366 {
367 accepted_steps++;
368 last_step_rejected = false;
369 }
370 else
371 {
372 if (t <= _end_time)
373 {
374 t -= prev_dt;
375 rejected_steps++;
376 last_step_rejected = true;
377 }
378 }
379 }
380
381 // adjust step size considering external communciation_point_calculators
382 for (auto const& time_step_constraint : time_step_constraints)
383 {
385 std::min(dt(), time_step_constraint(t, dt()))};
386 }
387
388 // Check whether the time stepping is stabilized
389 if (std::abs(dt() - prev_dt) < eps)
390 {
391 if (last_step_rejected)
392 {
393 OGS_FATAL(
394 "The new step size of {} is the same as that of the previous "
395 "rejected time step. \nPlease re-run ogs with a proper "
396 "adjustment in the numerical settings, \ne.g. those for time "
397 "stepper, local or global non-linear solver.",
398 dt);
399 }
400 else
401 {
402 DBUG("The time stepping is stabilized with the step size of {}.",
403 dt);
404 }
405 }
406
407 // Reset the time step with the minimum step size, dt
408 // Update the solution of the previous time step.
409 for (std::size_t i = 0; i < _per_process_data.size(); i++)
410 {
411 if (all_process_steps_accepted)
412 {
413 auto& ppd = *_per_process_data[i];
414 NumLib::updateTimeSteps(dt(), ppd.timestep_previous,
415 ppd.timestep_current);
416 auto& timestep_algorithm = ppd.timestep_algorithm;
417 timestep_algorithm->resetCurrentTimeStep(
418 dt(), ppd.timestep_previous, ppd.timestep_current);
419 }
420
421 auto& x = *_process_solutions[i];
422 auto& x_prev = *_process_solutions_prev[i];
423 if (all_process_steps_accepted)
424 {
425 MathLib::LinAlg::copy(x, x_prev); // pushState
426 }
427 else
428 {
429 if (t <= _end_time)
430 {
431 WARN(
432 "Time step {:d} was rejected {:d} times and it will be "
433 "repeated with a reduced step size.",
434 accepted_steps + 1, _repeating_times_of_rejected_step);
435 MathLib::LinAlg::copy(x_prev, x); // popState
436 }
437 }
438 }
439
440 return {dt, last_step_rejected};
441}
442
443std::vector<TimeLoop::TimeStepConstraintCallback>
445 std::vector<double>&& fixed_times) const
446{
447 std::vector<TimeStepConstraintCallback> const time_step_constraints{
448 [fixed_times = std::move(fixed_times)](NumLib::Time const& t, double dt)
449 { return NumLib::possiblyClampDtToNextFixedTime(t, dt, fixed_times); },
450 [this](NumLib::Time const& t, double dt) -> double
451 {
452 if (t < _end_time && _end_time < t + dt)
453 {
454 return _end_time() - t();
455 }
456 return dt;
457 }};
458 return time_step_constraints;
459}
460
463{
464 for (auto const& process_data : _per_process_data)
465 {
466 auto& pcs = process_data->process;
467 for (auto& output : _outputs)
468 {
469 output.addProcess(pcs);
470 }
471
472 setTimeDiscretizedODESystem(*process_data);
473
474 if (auto* conv_crit =
476 process_data->conv_crit.get()))
477 {
478 int const process_id = process_data->process_id;
479 conv_crit->setDOFTable(pcs.getDOFTable(process_id), pcs.getMesh());
480 }
481 }
482
483 // initial solution storage
486
488 {
489 _staggered_coupling->initializeCoupledSolutions(_process_solutions);
490 }
491
492 updateDeactivatedSubdomains(_per_process_data, _start_time());
493
494 auto const time_step_constraints = generateOutputTimeStepConstraints(
496
497 std::tie(_dt, _last_step_rejected) =
499 _rejected_steps, time_step_constraints);
500
501 // Output initial conditions
502 {
505 }
506
509}
510
512{
513 BaseLib::RunTime time_timestep;
514 time_timestep.start();
515
516 _current_time += _dt();
517
518 const std::size_t timesteps = _accepted_steps + 1;
519 // TODO(wenqing): , input option for time unit.
520 INFO("=== Time stepping at step #{:d} and time {} with step size {}",
521 timesteps, _current_time, _dt);
522
523 updateDeactivatedSubdomains(_per_process_data, _current_time());
524
527 INFO("[time] Time step #{:d} took {:g} s.", timesteps,
528 time_timestep.elapsed());
530}
531
533{
534 const double prev_dt = _dt();
535 auto const current_time = _current_time;
536
537 const std::size_t timesteps = _accepted_steps + 1;
538
539 auto const time_step_constraints = generateOutputTimeStepConstraints(
541
542 // _last_step_rejected is also checked in computeTimeStepping.
543 std::tie(_dt, _last_step_rejected) =
545 _rejected_steps, time_step_constraints);
546
548 {
549 outputSolutions(timesteps, current_time(), &Output::doOutput);
550 }
551
552 if (current_time == (_current_time + _dt()))
553 {
554 ERR("Time step size of {} is too small.\n"
555 "Time stepping stops at step {:d} and at time of {}.",
556 _dt, timesteps, _current_time);
557 return false;
558 }
559
561 {
562 return false;
563 }
564
565 return true;
566}
567
569{
570 INFO(
571 "The whole computation of the time stepping took {:d} steps, in which\n"
572 "\t the accepted steps are {:d}, and the rejected steps are {:d}.\n",
574
575 // output last time step
577 {
580 }
581}
582
584 std::size_t const timesteps)
585{
587
588 NumLib::NonlinearSolverStatus nonlinear_solver_status;
589
591 {
592 nonlinear_solver_status =
594 }
595 else
596 {
597 nonlinear_solver_status =
598 solveUncoupledEquationSystems(t, dt, timesteps);
599 }
600
601 // Run post time step only if the last iteration was successful.
602 // Otherwise it runs the risks to get the same errors as in the last
603 // iteration, an exception thrown in assembly, for example.
604 if (nonlinear_solver_status.error_norms_met)
605 {
606 // Later on, the timestep_algorithm might reject the timestep. We assume
607 // that this is a rare case, so still, we call preOutput() here. We
608 // don't expect a large overhead from it.
609 preOutputForAllProcesses(timesteps, t, dt, _end_time, _per_process_data,
611 _outputs);
612
616 }
617 return nonlinear_solver_status.error_norms_met;
618}
619
621 const NumLib::Time& t, const double dt, const std::size_t timestep_id,
622 ProcessData const& process_data, std::vector<GlobalVector*>& x,
623 std::vector<GlobalVector*> const& x_prev,
624 std::vector<Output> const& outputs)
625{
626 BaseLib::RunTime time_timestep_process;
627 time_timestep_process.start();
628
629 auto const nonlinear_solver_status = solveOneTimeStepOneProcess(
630 x, x_prev, timestep_id, t(), dt, process_data, outputs);
631
632 INFO("[time] Solving process #{:d} took {:g} s in time step #{:d}",
633 process_data.process_id, time_timestep_process.elapsed(), timestep_id);
634
635 return nonlinear_solver_status;
636}
637
638static constexpr std::string_view timestepper_cannot_reduce_dt =
639 "Time stepper cannot reduce the time step size further.";
640
642 const NumLib::Time& t, const double dt, const std::size_t timestep_id)
643{
644 NumLib::NonlinearSolverStatus nonlinear_solver_status;
645
646 for (auto const& process_data : _per_process_data)
647 {
648 auto const process_id = process_data->process_id;
649 nonlinear_solver_status = solveMonolithicProcess(
650 t, dt, timestep_id, *process_data, _process_solutions,
652
653 process_data->nonlinear_solver_status = nonlinear_solver_status;
654 if (!nonlinear_solver_status.error_norms_met)
655 {
656 ERR("The nonlinear solver failed in time step #{:d} at t = {} s "
657 "for process #{:d}.",
658 timestep_id, t, process_id);
659
660 if (!process_data->timestep_algorithm->canReduceTimestepSize(
661 process_data->timestep_current,
662 process_data->timestep_previous))
663 {
664 // save unsuccessful solution
665 for (auto const& output : _outputs)
666 {
667 output.doOutputAlways(
668 process_data->process, process_id, timestep_id, t,
669 process_data->nonlinear_solver_status.number_iterations,
671 }
673 }
674
675 return nonlinear_solver_status;
676 }
677 }
678
679 return nonlinear_solver_status;
680}
681
684 const NumLib::Time& t, const double dt, const std::size_t timestep_id)
685{
686 auto const nonlinear_solver_status =
688 t(), dt, timestep_id, _process_solutions, _process_solutions_prev,
690
691 _last_step_rejected = nonlinear_solver_status.error_norms_met;
692
693 {
694 for (auto const& process_data : _per_process_data)
695 {
696 auto& pcs = process_data->process;
697 int const process_id = process_data->process_id;
698 auto& ode_sys = *process_data->tdisc_ode_sys;
699 pcs.solveReactionEquation(_process_solutions,
700 _process_solutions_prev, t(), dt, ode_sys,
701 process_id);
702 }
703 }
704
705 return nonlinear_solver_status;
706}
707
708template <typename OutputClassMember>
709void TimeLoop::outputSolutions(unsigned timestep, const double t,
710 OutputClassMember output_class_member) const
711{
712 for (auto const& process_data : _per_process_data)
713 {
714 // If nonlinear solver diverged, the solution has already been
715 // saved.
716 if (!process_data->nonlinear_solver_status.error_norms_met)
717 {
718 continue;
719 }
720
721 auto const process_id = process_data->process_id;
722 auto const& pcs = process_data->process;
723
724 for (auto const& output_object : _outputs)
725 {
726 (output_object.*output_class_member)(
727 pcs, process_id, timestep, NumLib::Time(t),
728 process_data->nonlinear_solver_status.number_iterations,
730 }
731 }
732}
733
745
747 const double dt) const
748{
749 for (auto const& process_data : _per_process_data)
750 {
751 // If nonlinear solver diverged, the solution has already been
752 // saved.
753 if (!process_data->nonlinear_solver_status.error_norms_met)
754 {
755 continue;
756 }
757
758 auto const process_id = process_data->process_id;
759 auto& pcs = process_data->process;
760
761 process_data->time_disc->nextTimestep(t(), dt);
762
763 pcs.preTimestep(_process_solutions, _start_time(), dt, process_id);
764
765 pcs.preOutput(_start_time(), dt, _process_solutions,
766 _process_solutions_prev, process_id);
767
768 // Update secondary variables, which might be uninitialized, before
769 // output.
770 pcs.computeSecondaryVariable(_start_time(), dt, _process_solutions,
771 *_process_solutions_prev[process_id],
772 process_id);
773 }
774}
775} // 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 doOutputLastTimestep(Process const &process, const int process_id, int const timestep, const NumLib::Time &t, int const iteration, std::vector< GlobalVector * > const &xs) const
Definition Output.cpp:355
void doOutput(Process const &process, const int process_id, int const timestep, const NumLib::Time &t, int const iteration, std::vector< GlobalVector * > const &xs) const
Definition Output.cpp:336
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:641
void outputLastTimeStep() const
Definition TimeLoop.cpp:568
void preOutputInitialConditions(NumLib::Time const &t, const double dt) const
Definition TimeLoop.cpp:746
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:295
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:269
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:709
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:462
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:583
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:683
std::vector< TimeStepConstraintCallback > generateOutputTimeStepConstraints(std::vector< double > &&fixed_times) const
Definition TimeLoop.cpp:444
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:638
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:282
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:433
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:620
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