#include <adaptive_time_solver.h>
Public Types | |
typedef FirstOrderUnsteadySolver | Parent |
typedef DifferentiableSystem | sys_type |
Public Member Functions | |
AdaptiveTimeSolver (sys_type &s) | |
virtual | ~AdaptiveTimeSolver () |
virtual void | init () override |
virtual void | reinit () override |
virtual void | solve () override=0 |
virtual void | advance_timestep () override |
virtual Real | error_order () const override |
virtual bool | element_residual (bool get_jacobian, DiffContext &) override |
virtual bool | side_residual (bool get_jacobian, DiffContext &) override |
virtual bool | nonlocal_residual (bool get_jacobian, DiffContext &) override |
virtual std::unique_ptr< DiffSolver > & | diff_solver () override |
virtual std::unique_ptr< LinearSolver< Number > > & | linear_solver () override |
virtual unsigned int | time_order () const override |
virtual void | init_data () override |
virtual void | adjoint_advance_timestep () override |
virtual void | retrieve_timestep () override |
Number | old_nonlinear_solution (const dof_id_type global_dof_number) const |
virtual Real | du (const SystemNorm &norm) const override |
virtual bool | is_steady () const override |
virtual void | before_timestep () |
const sys_type & | system () const |
sys_type & | system () |
void | set_solution_history (const SolutionHistory &_solution_history) |
bool | is_adjoint () const |
void | set_is_adjoint (bool _is_adjoint_value) |
Static Public Member Functions | |
static std::string | get_info () |
static void | print_info (std::ostream &out=libMesh::out) |
static unsigned int | n_objects () |
static void | enable_print_counter_info () |
static void | disable_print_counter_info () |
Public Attributes | |
std::unique_ptr< UnsteadySolver > | core_time_solver |
SystemNorm | component_norm |
std::vector< float > | component_scale |
Real | target_tolerance |
Real | upper_tolerance |
Real | max_deltat |
Real | min_deltat |
Real | max_growth |
bool | global_tolerance |
std::unique_ptr< NumericVector< Number > > | old_local_nonlinear_solution |
bool | quiet |
unsigned int | reduce_deltat_on_diffsolver_failure |
Protected Types | |
typedef bool(DifferentiablePhysics::* | ResFuncType) (bool, DiffContext &) |
typedef void(DiffContext::* | ReinitFuncType) (Real) |
typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
virtual Real | calculate_norm (System &, NumericVector< Number > &) |
void | prepare_accel (DiffContext &context) |
bool | compute_second_order_eqns (bool compute_jacobian, DiffContext &c) |
void | increment_constructor_count (const std::string &name) |
void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
Real | last_deltat |
bool | first_solve |
bool | first_adjoint_step |
std::unique_ptr< DiffSolver > | _diff_solver |
std::unique_ptr< LinearSolver< Number > > | _linear_solver |
sys_type & | _system |
std::unique_ptr< SolutionHistory > | solution_history |
Static Protected Attributes | |
static Counts | _counts |
static Threads::atomic< unsigned int > | _n_objects |
static Threads::spin_mutex | _mutex |
static bool | _enable_print_counter = true |
This class wraps another UnsteadySolver derived class, and compares the results of timestepping with deltat and timestepping with 2*deltat to adjust future timestep lengths.
Currently this class only works on fully coupled Systems
This class is part of the new DifferentiableSystem framework, which is still experimental. Users of this framework should beware of bugs and future API changes.
Definition at line 49 of file adaptive_time_solver.h.
|
protectedinherited |
Data structure to log the information. The log is identified by the class name.
Definition at line 117 of file reference_counter.h.
The parent class
Definition at line 55 of file adaptive_time_solver.h.
|
protectedinherited |
Definition at line 273 of file time_solver.h.
|
protectedinherited |
Definitions of argument types for use in refactoring subclasses.
Definition at line 271 of file time_solver.h.
|
inherited |
The type of system
Definition at line 65 of file time_solver.h.
|
explicit |
Constructor. Requires a reference to the system to be solved.
Definition at line 27 of file adaptive_time_solver.C.
References libMesh::UnsteadySolver::old_local_nonlinear_solution.
|
virtual |
Destructor.
Definition at line 47 of file adaptive_time_solver.C.
References libMesh::UnsteadySolver::old_local_nonlinear_solution.
|
overridevirtualinherited |
This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will be done before every UnsteadySolver::adjoint_solve().
Reimplemented from libMesh::TimeSolver.
Reimplemented in libMesh::NewmarkSolver.
Definition at line 178 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, libMesh::UnsteadySolver::first_adjoint_step, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::NumericVector< T >::localize(), libMesh::UnsteadySolver::old_local_nonlinear_solution, libMesh::TimeSolver::solution_history, and libMesh::System::time.
|
overridevirtual |
This method advances the solution to the next timestep, after a solve() has been performed. Often this will be done after every UnsteadySolver::solve(), but adaptive mesh refinement and/or adaptive time step selection may require some solve() steps to be repeated.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 87 of file adaptive_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::UnsteadySolver::first_solve, libMesh::System::get_vector(), last_deltat, libMesh::System::solution, and libMesh::System::time.
|
inlinevirtualinherited |
This method is for subclasses or users to override to do arbitrary processing between timesteps
Definition at line 167 of file time_solver.h.
|
protectedvirtual |
A helper function to calculate error norms
Definition at line 155 of file adaptive_time_solver.C.
References libMesh::System::calculate_norm(), and component_norm.
Referenced by libMesh::TwostepTimeSolver::solve().
|
protectedinherited |
If there are second order variables, then we need to compute their residual equations and corresponding Jacobian. The residual equation will simply be , where is the second order variable add by the user and is the variable added by the time-solver as the "velocity" variable.
Definition at line 32 of file first_order_unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_residual(), libMesh::DiffContext::get_elem_solution_derivative(), libMesh::DiffContext::get_elem_solution_rate_derivative(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), libMesh::DifferentiableSystem::get_second_order_dot_var(), libMesh::DiffContext::get_system(), libMesh::FEMContext::interior_rate(), libMesh::FEMContext::interior_value(), libMesh::DifferentiablePhysics::is_second_order_var(), libMesh::QBase::n_points(), libMesh::DiffContext::n_vars(), libMesh::Variable::type(), and libMesh::System::variable().
Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::EulerSolver::nonlocal_residual(), and libMesh::Euler2Solver::nonlocal_residual().
|
overridevirtual |
An implicit linear or nonlinear solver to use at each timestep.
Reimplemented from libMesh::TimeSolver.
Definition at line 141 of file adaptive_time_solver.C.
References core_time_solver.
|
staticinherited |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
Referenced by libMesh::LibMeshInit::LibMeshInit().
|
overridevirtualinherited |
Computes the size of ||u^{n+1} - u^{n}|| in some norm.
Implements libMesh::TimeSolver.
Definition at line 227 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::calculate_norm(), libMesh::System::get_vector(), and libMesh::System::solution.
|
overridevirtual |
This method is passed on to the core_time_solver
Implements libMesh::TimeSolver.
Definition at line 111 of file adaptive_time_solver.C.
References core_time_solver.
|
staticinherited |
Methods to enable/disable the reference counter output from print_info()
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
|
overridevirtual |
This method is passed on to the core_time_solver
Implements libMesh::UnsteadySolver.
Definition at line 102 of file adaptive_time_solver.C.
References core_time_solver.
|
staticinherited |
Gets a string containing the reference information.
Definition at line 47 of file reference_counter.C.
References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().
Referenced by libMesh::ReferenceCounter::print_info().
|
inlineprotectedinherited |
Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 181 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
|
inlineprotectedinherited |
Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 194 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
|
overridevirtual |
The initialization function. This method is used to initialize internal data structures before a simulation begins.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 58 of file adaptive_time_solver.C.
References core_time_solver, and libMesh::UnsteadySolver::old_local_nonlinear_solution.
|
overridevirtualinherited |
The data initialization function. This method is used to initialize internal data structures after the underlying System has been initialized
Reimplemented from libMesh::TimeSolver.
Reimplemented in libMesh::SecondOrderUnsteadySolver.
Definition at line 55 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::GHOSTED, libMesh::TimeSolver::init_data(), libMesh::System::n_dofs(), libMesh::System::n_local_dofs(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::SERIAL.
Referenced by libMesh::SecondOrderUnsteadySolver::init_data().
|
inlineinherited |
Accessor for querying whether we need to do a primal or adjoint solve
Definition at line 233 of file time_solver.h.
References libMesh::TimeSolver::_is_adjoint.
Referenced by libMesh::FEMSystem::build_context().
|
inlineoverridevirtualinherited |
This is not a steady-state solver.
Implements libMesh::TimeSolver.
Definition at line 154 of file unsteady_solver.h.
|
overridevirtual |
An implicit linear solver to use for adjoint and sensitivity problems.
Reimplemented from libMesh::TimeSolver.
Definition at line 148 of file adaptive_time_solver.C.
References core_time_solver.
|
inlinestaticinherited |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 83 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
|
overridevirtual |
This method is passed on to the core_time_solver
Implements libMesh::TimeSolver.
Definition at line 131 of file adaptive_time_solver.C.
References core_time_solver.
|
inherited |
Definition at line 216 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::n_dofs(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.
Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), and libMesh::NewmarkSolver::_general_residual().
|
protectedinherited |
If there are second order variables in the system, then we also prepare the accel for those variables so the user can treat them as such.
Definition at line 25 of file first_order_unsteady_solver.C.
References libMesh::DiffContext::elem_solution_accel_derivative, libMesh::DiffContext::get_elem_solution_accel(), libMesh::DiffContext::get_elem_solution_rate(), and libMesh::DiffContext::get_elem_solution_rate_derivative().
Referenced by libMesh::EulerSolver::_general_residual(), and libMesh::Euler2Solver::_general_residual().
|
staticinherited |
Prints the reference information, by default to libMesh::out
.
Definition at line 87 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
|
overridevirtual |
The reinitialization function. This method is used to resize internal data vectors after a mesh change.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 77 of file adaptive_time_solver.C.
References core_time_solver.
|
overridevirtualinherited |
This method retrieves all the stored solutions at the current system.time
Reimplemented from libMesh::TimeSolver.
Reimplemented in libMesh::SecondOrderUnsteadySolver.
Definition at line 204 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::NumericVector< T >::localize(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::TimeSolver::solution_history.
|
inlineinherited |
Accessor for setting whether we need to do a primal or adjoint solve
Definition at line 240 of file time_solver.h.
References libMesh::TimeSolver::_is_adjoint.
Referenced by libMesh::DifferentiableSystem::adjoint_solve(), libMesh::FEMSystem::postprocess(), and libMesh::DifferentiableSystem::solve().
|
inherited |
A setter function users will employ if they need to do something other than save no solution history
Definition at line 97 of file time_solver.C.
References libMesh::SolutionHistory::clone(), and libMesh::TimeSolver::solution_history.
|
overridevirtual |
This method is passed on to the core_time_solver
Implements libMesh::TimeSolver.
Definition at line 121 of file adaptive_time_solver.C.
References core_time_solver.
|
overridepure virtual |
This method solves for the solution at the next timestep. Usually we will only need to solve one (non)linear system per timestep, but more complex subclasses may override this.
Reimplemented from libMesh::UnsteadySolver.
Implemented in libMesh::TwostepTimeSolver.
|
inlineinherited |
Definition at line 172 of file time_solver.h.
References libMesh::TimeSolver::_system.
Referenced by libMesh::TimeSolver::reinit(), and libMesh::TimeSolver::solve().
|
inlineinherited |
Definition at line 177 of file time_solver.h.
References libMesh::TimeSolver::_system.
|
inlineoverridevirtualinherited |
For example, EulerSolver will have time_order()
= 1 and NewmarkSolver will have time_order()
= 2.
Implements libMesh::UnsteadySolver.
Definition at line 90 of file first_order_unsteady_solver.h.
|
staticprotectedinherited |
Actually holds the data.
Definition at line 122 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
|
protectedinherited |
An implicit linear or nonlinear solver to use at each timestep.
Definition at line 248 of file time_solver.h.
Referenced by libMesh::NewmarkSolver::compute_initial_accel(), libMesh::TimeSolver::diff_solver(), and libMesh::UnsteadySolver::solve().
|
staticprotectedinherited |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 141 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
|
protectedinherited |
An implicit linear solver to use for adjoint problems.
Definition at line 253 of file time_solver.h.
Referenced by libMesh::TimeSolver::linear_solver().
|
staticprotectedinherited |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 135 of file reference_counter.h.
|
staticprotectedinherited |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 130 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
|
protectedinherited |
A reference to the system we are solving.
Definition at line 258 of file time_solver.h.
Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::NewmarkSolver::advance_timestep(), advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::NewmarkSolver::compute_initial_accel(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::UnsteadySolver::du(), libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::EigenTimeSolver::element_residual(), libMesh::SecondOrderUnsteadySolver::init(), libMesh::UnsteadySolver::init(), libMesh::TimeSolver::init(), libMesh::EigenTimeSolver::init(), libMesh::SecondOrderUnsteadySolver::init_data(), libMesh::UnsteadySolver::init_data(), libMesh::TimeSolver::init_data(), libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::SecondOrderUnsteadySolver::old_solution_accel(), libMesh::SecondOrderUnsteadySolver::old_solution_rate(), libMesh::NewmarkSolver::project_initial_accel(), libMesh::SecondOrderUnsteadySolver::project_initial_rate(), libMesh::SecondOrderUnsteadySolver::reinit(), libMesh::UnsteadySolver::reinit(), libMesh::TimeSolver::reinit(), libMesh::UnsteadySolver::retrieve_timestep(), libMesh::EigenTimeSolver::side_residual(), libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), libMesh::EigenTimeSolver::solve(), and libMesh::TimeSolver::system().
SystemNorm libMesh::AdaptiveTimeSolver::component_norm |
Error calculations are done in this norm, DISCRETE_L2 by default.
Definition at line 119 of file adaptive_time_solver.h.
Referenced by calculate_norm().
std::vector<float> libMesh::AdaptiveTimeSolver::component_scale |
If component_norms is non-empty, each variable's contribution to the error of a system will also be scaled by component_scale[var], unless component_scale is empty in which case all variables will be weighted equally.
Definition at line 127 of file adaptive_time_solver.h.
std::unique_ptr<UnsteadySolver> libMesh::AdaptiveTimeSolver::core_time_solver |
This object is used to take timesteps
Definition at line 114 of file adaptive_time_solver.h.
Referenced by diff_solver(), element_residual(), error_order(), init(), linear_solver(), nonlocal_residual(), reinit(), side_residual(), libMesh::TwostepTimeSolver::solve(), and libMesh::TwostepTimeSolver::TwostepTimeSolver().
|
protectedinherited |
A bool that will be true the first time adjoint_advance_timestep() is called, (when the primal solution is to be used to set adjoint boundary conditions) and false thereafter
Definition at line 168 of file unsteady_solver.h.
Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep().
|
protectedinherited |
A bool that will be true the first time solve() is called, and false thereafter
Definition at line 162 of file unsteady_solver.h.
Referenced by libMesh::NewmarkSolver::advance_timestep(), advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::TwostepTimeSolver::solve(), and libMesh::UnsteadySolver::solve().
bool libMesh::AdaptiveTimeSolver::global_tolerance |
This flag, which is true by default, grows (shrinks) the timestep based on the expected global accuracy of the timestepping scheme. Global in this sense means the cumulative final-time accuracy of the scheme. For example, the backward Euler scheme's truncation error is locally of order 2, so that after N timesteps of size deltat, the result is first-order accurate. If you set this to false, you can grow (shrink) your timestep based on the local accuracy rather than the global accuracy of the core TimeSolver.
Definition at line 198 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
|
protected |
We need to store the value of the last deltat used, so that advance_timestep() will increment the system time correctly.
Definition at line 207 of file adaptive_time_solver.h.
Referenced by advance_timestep(), and libMesh::TwostepTimeSolver::solve().
Real libMesh::AdaptiveTimeSolver::max_deltat |
Do not allow the adaptive time solver to select deltat > max_deltat. If you use the default max_deltat=0.0, then deltat is unlimited.
Definition at line 167 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
Real libMesh::AdaptiveTimeSolver::max_growth |
Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old deltat. If you use the default max_growth=0.0, then the deltat growth is unlimited.
Definition at line 181 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
Real libMesh::AdaptiveTimeSolver::min_deltat |
Do not allow the adaptive time solver to select deltat < min_deltat. The default value is 0.0.
Definition at line 173 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
|
inherited |
Serial vector of _system.get_vector("_old_nonlinear_solution")
Definition at line 138 of file unsteady_solver.h.
Referenced by AdaptiveTimeSolver(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), init(), libMesh::UnsteadySolver::init_data(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::UnsteadySolver::reinit(), libMesh::UnsteadySolver::retrieve_timestep(), and ~AdaptiveTimeSolver().
|
inherited |
Print extra debugging information if quiet == false.
Definition at line 192 of file time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), and libMesh::EigenTimeSolver::solve().
|
inherited |
This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat and let the DiffSolver repeat the latest failed solve with a reduced timestep.
Definition at line 221 of file time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve(), and libMesh::UnsteadySolver::solve().
|
protectedinherited |
A std::unique_ptr to a SolutionHistory object. Default is NoSolutionHistory, which the user can override by declaring a different kind of SolutionHistory in the application
Definition at line 265 of file time_solver.h.
Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::UnsteadySolver::retrieve_timestep(), and libMesh::TimeSolver::set_solution_history().
Real libMesh::AdaptiveTimeSolver::target_tolerance |
This tolerance is the target relative error between an exact time integration and a single time step output, scaled by deltat. integrator, scaled by deltat. If the estimated error exceeds or undershoots the target error tolerance, future timesteps will be run with deltat shrunk or grown to compensate.
The default value is 1.0e-2; obviously users should select their own tolerance.
If a negative target_tolerance is specified, then its absolute value is used to scale the estimated error from the first simulation time step, and this becomes the target tolerance of all future time steps.
Definition at line 144 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
Real libMesh::AdaptiveTimeSolver::upper_tolerance |
This tolerance is the maximum relative error between an exact time integration and a single time step output, scaled by deltat. If this error tolerance is exceeded by the estimated error of the current time step, that time step will be repeated with a smaller deltat.
If you use the default upper_tolerance=0.0, then the current time step will not be repeated regardless of estimated error.
If a negative upper_tolerance is specified, then its absolute value is used to scale the estimated error from the first simulation time step, and this becomes the upper tolerance of all future time steps.
Definition at line 161 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().