Go to the documentation of this file.
20 #ifndef LIBMESH_ADAPTIVE_TIME_SOLVER_H
21 #define LIBMESH_ADAPTIVE_TIME_SOLVER_H
24 #include "libmesh/system_norm.h"
25 #include "libmesh/first_order_unsteady_solver.h"
69 virtual void init()
override;
71 virtual void reinit()
override;
73 virtual void solve()
override = 0;
103 virtual std::unique_ptr<DiffSolver> &
diff_solver()
override;
109 virtual std::unique_ptr<LinearSolver<Number>> &
linear_solver()
override;
219 #endif // LIBMESH_ADAPTIVE_TIME_SOLVER_H
Manages consistently variables, degrees of freedom, and coefficient vectors.
virtual void solve() override=0
This method solves for the solution at the next timestep.
Generic class from which first order UnsteadySolvers should subclass.
virtual Real calculate_norm(System &, NumericVector< Number > &)
A helper function to calculate error norms.
virtual bool nonlocal_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
The libMesh namespace provides an interface to certain functionality in the library.
AdaptiveTimeSolver(sys_type &s)
Constructor.
Real upper_tolerance
This tolerance is the maximum relative error between an exact time integration and a single time step...
This class provides a specific system class.
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver() override
An implicit linear solver to use for adjoint and sensitivity problems.
virtual void advance_timestep() override
This method advances the solution to the next timestep, after a solve() has been performed.
FirstOrderUnsteadySolver Parent
The parent class.
std::unique_ptr< UnsteadySolver > core_time_solver
This object is used to take timesteps.
SystemNorm component_norm
Error calculations are done in this norm, DISCRETE_L2 by default.
Real last_deltat
We need to store the value of the last deltat used, so that advance_timestep() will increment the sys...
bool global_tolerance
This flag, which is true by default, grows (shrinks) the timestep based on the expected global accura...
Real max_growth
Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old de...
virtual Real error_order() const override
This method is passed on to the core_time_solver.
virtual bool side_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
This class provides all data required for a physics package (e.g.
This class wraps another UnsteadySolver derived class, and compares the results of timestepping with ...
Real max_deltat
Do not allow the adaptive time solver to select deltat > max_deltat.
std::vector< float > component_scale
If component_norms is non-empty, each variable's contribution to the error of a system will also be s...
virtual bool element_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
virtual void init() override
The initialization function.
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
Real min_deltat
Do not allow the adaptive time solver to select deltat < min_deltat.
virtual void reinit() override
The reinitialization function.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual std::unique_ptr< DiffSolver > & diff_solver() override
An implicit linear or nonlinear solver to use at each timestep.
Real target_tolerance
This tolerance is the target relative error between an exact time integration and a single time step ...
virtual ~AdaptiveTimeSolver()
Destructor.