Go to the documentation of this file.
   18 #include "libmesh/adaptive_time_solver.h" 
   19 #include "libmesh/diff_system.h" 
   20 #include "libmesh/numeric_vector.h" 
   30     target_tolerance(1.e-3),
 
   35     global_tolerance(true)
 
   72     std::unique_ptr<NumericVector<Number>>(
core_time_solver->old_local_nonlinear_solution.get());
 
   94   old_nonlinear_soln = nonlinear_solution;
 
  
Manages consistently variables, degrees of freedom, and coefficient vectors.
 
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.
 
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.
 
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...
 
virtual Real error_order() const override
This method is passed on to the core_time_solver.
 
sys_type & _system
A reference to the system we are solving.
 
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.
 
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
 
bool first_solve
A bool that will be true the first time solve() is called, and false thereafter.
 
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
 
virtual bool element_residual(bool get_jacobian, DiffContext &) override
This method is passed on to the core_time_solver.
 
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
Serial vector of _system.get_vector("_old_nonlinear_solution")
 
virtual void init() override
The initialization function.
 
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
 
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.
 
const NumericVector< Number > & get_vector(const std::string &vec_name) const
 
virtual ~AdaptiveTimeSolver()
Destructor.