libMesh
|
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. More...
#include <twostep_time_solver.h>
Public Types | |
typedef AdaptiveTimeSolver | Parent |
The parent class. More... | |
typedef DifferentiableSystem | sys_type |
The type of system. More... | |
Public Member Functions | |
TwostepTimeSolver (sys_type &s) | |
Constructor. More... | |
~TwostepTimeSolver () | |
Destructor. More... | |
virtual void | solve () override |
This method solves for the solution at the next timestep. More... | |
virtual std::pair< unsigned int, Real > | adjoint_solve (const QoISet &qoi_indices) override |
This method solves for the adjoint solution at the next adjoint timestep (or a steady state adjoint solve) More... | |
virtual void | integrate_qoi_timestep () override |
A method to integrate the system::QoI functionals. More... | |
virtual void | integrate_adjoint_sensitivity (const QoISet &qois, const ParameterVector ¶meter_vector, SensitivityData &sensitivities) override |
A method to integrate the adjoint sensitivity w.r.t a given parameter vector. More... | |
virtual void | integrate_adjoint_refinement_error_estimate (AdjointRefinementEstimator &adjoint_refinement_error_estimator, ErrorVector &QoI_elementwise_error) override |
A method to compute the adjoint refinement error estimate at the current timestep. More... | |
virtual void | init () override |
The initialization function. More... | |
virtual void | reinit () override |
The reinitialization function. More... | |
virtual void | advance_timestep () override |
This method advances the solution to the next timestep, after a solve() has been performed. More... | |
virtual void | adjoint_advance_timestep () override |
This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. More... | |
virtual void | retrieve_timestep () override |
This method retrieves all the stored solutions at the current system.time. More... | |
virtual Real | last_completed_timestep_size () override |
Returns system.deltat if fixed timestep solver is used, the complete timestep size (sum of all substeps) if the adaptive time solver is used. More... | |
virtual Real | error_order () const override |
This method is passed on to the core_time_solver. More... | |
virtual bool | element_residual (bool get_jacobian, DiffContext &) override |
This method is passed on to the core_time_solver. More... | |
virtual bool | side_residual (bool get_jacobian, DiffContext &) override |
This method is passed on to the core_time_solver. More... | |
virtual bool | nonlocal_residual (bool get_jacobian, DiffContext &) override |
This method is passed on to the core_time_solver. More... | |
virtual std::unique_ptr< DiffSolver > & | diff_solver () override |
An implicit linear or nonlinear solver to use at each timestep. More... | |
virtual std::unique_ptr< LinearSolver< Number > > & | linear_solver () override |
An implicit linear solver to use for adjoint and sensitivity problems. More... | |
virtual unsigned int | time_order () const override |
virtual void | init_adjoints () override |
Add adjoint vectors and old_adjoint_vectors as per the indices of QoISet. More... | |
virtual void | init_data () override |
The data initialization function. More... | |
void | update () |
Number | old_nonlinear_solution (const dof_id_type global_dof_number) const |
virtual Real | du (const SystemNorm &norm) const override |
Computes the size of ||u^{n+1} - u^{n}|| in some norm. More... | |
virtual bool | is_steady () const override |
This is not a steady-state solver. More... | |
void | set_first_adjoint_step (bool first_adjoint_step_setting) |
A setter for the first_adjoint_step boolean. More... | |
void | set_first_solve (bool first_solve_setting) |
virtual void | before_timestep () |
This method is for subclasses or users to override to do arbitrary processing between timesteps. More... | |
const sys_type & | system () const |
sys_type & | system () |
void | set_solution_history (const SolutionHistory &_solution_history) |
A setter function users will employ if they need to do something other than save no solution history. More... | |
SolutionHistory & | get_solution_history () |
A getter function that returns a reference to the solution history object owned by TimeSolver. More... | |
bool | is_adjoint () const |
Accessor for querying whether we need to do a primal or adjoint solve. More... | |
void | set_is_adjoint (bool _is_adjoint_value) |
Accessor for setting whether we need to do a primal or adjoint solve. More... | |
Static Public Member Functions | |
static std::string | get_info () |
Gets a string containing the reference information. More... | |
static void | print_info (std::ostream &out_stream=libMesh::out) |
Prints the reference information, by default to libMesh::out . More... | |
static unsigned int | n_objects () |
Prints the number of outstanding (created, but not yet destroyed) objects. More... | |
static void | enable_print_counter_info () |
Methods to enable/disable the reference counter output from print_info() More... | |
static void | disable_print_counter_info () |
Public Attributes | |
std::unique_ptr< UnsteadySolver > | core_time_solver |
This object is used to take timesteps. More... | |
SystemNorm | component_norm |
Error calculations are done in this norm, DISCRETE_L2 by default. More... | |
std::vector< float > | 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. More... | |
Real | target_tolerance |
This tolerance is the target relative error between an exact time integration and a single time step output, scaled by deltat. More... | |
Real | upper_tolerance |
This tolerance is the maximum relative error between an exact time integration and a single time step output, scaled by deltat. More... | |
Real | max_deltat |
Do not allow the adaptive time solver to select deltat > max_deltat. More... | |
Real | min_deltat |
Do not allow the adaptive time solver to select deltat < min_deltat. More... | |
Real | max_growth |
Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old deltat. More... | |
Real | completed_timestep_size |
The adaptive time solver's have two notions of deltat. More... | |
bool | global_tolerance |
This flag, which is true by default, grows (shrinks) the timestep based on the expected global accuracy of the timestepping scheme. More... | |
std::shared_ptr< NumericVector< Number > > | old_local_nonlinear_solution |
Serial vector of _system.get_vector("_old_nonlinear_solution") This is a shared_ptr so that it can be shared between different derived class instances, as in e.g. More... | |
bool | quiet |
Print extra debugging information if quiet == false. More... | |
unsigned int | reduce_deltat_on_diffsolver_failure |
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. More... | |
Protected Types | |
typedef bool(DifferentiablePhysics::* | ResFuncType) (bool, DiffContext &) |
Definitions of argument types for use in refactoring subclasses. More... | |
typedef void(DiffContext::* | ReinitFuncType) (Real) |
typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Data structure to log the information. More... | |
Protected Member Functions | |
virtual Real | calculate_norm (System &, NumericVector< Number > &) |
A helper function to calculate error norms. More... | |
void | prepare_accel (DiffContext &context) |
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. More... | |
bool | compute_second_order_eqns (bool compute_jacobian, DiffContext &c) |
If there are second order variables, then we need to compute their residual equations and corresponding Jacobian. More... | |
void | increment_constructor_count (const std::string &name) noexcept |
Increments the construction counter. More... | |
void | increment_destructor_count (const std::string &name) noexcept |
Increments the destruction counter. More... | |
Protected Attributes | |
bool | first_solve |
A bool that will be true the first time solve() is called, and false thereafter. More... | |
bool | first_adjoint_step |
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. More... | |
std::vector< std::unique_ptr< NumericVector< Number > > > | old_adjoints |
A vector of pointers to vectors holding the adjoint solution at the last time step. More... | |
Real | last_step_deltat |
We will need to move the system.time around to ensure that residuals are built with the right deltat and the right time. More... | |
Real | next_step_deltat |
std::unique_ptr< DiffSolver > | _diff_solver |
An implicit linear or nonlinear solver to use at each timestep. More... | |
std::unique_ptr< LinearSolver< Number > > | _linear_solver |
An implicit linear solver to use for adjoint problems. More... | |
sys_type & | _system |
A reference to the system we are solving. More... | |
std::unique_ptr< SolutionHistory > | solution_history |
A std::unique_ptr to a SolutionHistory object. More... | |
Real | last_deltat |
The deltat for the last completed timestep before the current one. More... | |
Static Protected Attributes | |
static Counts | _counts |
Actually holds the data. More... | |
static Threads::atomic< unsigned int > | _n_objects |
The number of objects. More... | |
static Threads::spin_mutex | _mutex |
Mutual exclusion object to enable thread-safe reference counting. More... | |
static bool | _enable_print_counter = true |
Flag to control whether reference count information is printed when print_info is called. More... | |
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 50 of file twostep_time_solver.h.
|
protectedinherited |
Data structure to log the information.
The log is identified by the class name.
Definition at line 119 of file reference_counter.h.
The parent class.
Definition at line 56 of file twostep_time_solver.h.
|
protectedinherited |
Definition at line 327 of file time_solver.h.
|
protectedinherited |
Definitions of argument types for use in refactoring subclasses.
Definition at line 325 of file time_solver.h.
|
inherited |
The type of system.
Definition at line 69 of file time_solver.h.
|
explicit |
Constructor.
Requires a reference to the system to be solved.
Definition at line 40 of file twostep_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver.
|
default |
Destructor.
|
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::UnsteadySolver.
Definition at line 126 of file adaptive_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::DifferentiableSystem::deltat, libMesh::UnsteadySolver::first_adjoint_step, libMesh::System::get_adjoint_solution(), libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::NumericVector< T >::localize(), libMesh::make_range(), libMesh::System::n_qois(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::System::time.
|
overridevirtual |
This method solves for the adjoint solution at the next adjoint timestep (or a steady state adjoint solve)
Implements libMesh::AdaptiveTimeSolver.
Definition at line 298 of file twostep_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::calculate_norm(), libMesh::AdaptiveTimeSolver::completed_timestep_size, libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::DifferentiableSystem::deltat, libMesh::System::get_adjoint_solution(), libMesh::H1, libMesh::TimeSolver::last_deltat, libMesh::make_range(), libMesh::System::n_qois(), libMesh::System::n_vars(), libMesh::out, libMesh::Real, and libMesh::System::time.
|
overridevirtualinherited |
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 89 of file adaptive_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::AdaptiveTimeSolver::completed_timestep_size, libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::UnsteadySolver::first_solve, 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::System::solution, and libMesh::System::time.
Referenced by solve().
|
inlinevirtualinherited |
This method is for subclasses or users to override to do arbitrary processing between timesteps.
Definition at line 205 of file time_solver.h.
|
protectedvirtualinherited |
A helper function to calculate error norms.
Definition at line 225 of file adaptive_time_solver.C.
References libMesh::System::calculate_norm(), and libMesh::AdaptiveTimeSolver::component_norm.
Referenced by 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 \( \dot{u} - v = 0 \), where \( u \) is the second order variable add by the user and \( v \) 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, compute_jacobian(), 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::libmesh_assert(), libMesh::make_range(), 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().
|
overridevirtualinherited |
An implicit linear or nonlinear solver to use at each timestep.
Reimplemented from libMesh::TimeSolver.
Definition at line 211 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver.
|
staticinherited |
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
|
overridevirtualinherited |
Computes the size of ||u^{n+1} - u^{n}|| in some norm.
Implements libMesh::TimeSolver.
Definition at line 348 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::calculate_norm(), libMesh::System::get_vector(), libMesh::TensorTools::norm(), and libMesh::System::solution.
|
overridevirtualinherited |
This method is passed on to the core_time_solver.
Implements libMesh::TimeSolver.
Definition at line 181 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver, and libMesh::libmesh_assert().
|
staticinherited |
Methods to enable/disable the reference counter output from print_info()
Definition at line 94 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
|
overridevirtualinherited |
This method is passed on to the core_time_solver.
Implements libMesh::UnsteadySolver.
Definition at line 172 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver, and libMesh::libmesh_assert().
|
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().
|
inherited |
A getter function that returns a reference to the solution history object owned by TimeSolver.
Definition at line 124 of file time_solver.C.
References libMesh::TimeSolver::solution_history.
Referenced by libMesh::AdaptiveTimeSolver::init().
|
inlineprotectednoexceptinherited |
Increments the construction counter.
Should be called in the constructor of any derived class that will be reference counted.
Definition at line 183 of file reference_counter.h.
References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
|
inlineprotectednoexceptinherited |
Increments the destruction counter.
Should be called in the destructor of any derived class that will be reference counted.
Definition at line 207 of file reference_counter.h.
References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
|
overridevirtualinherited |
The initialization function.
This method is used to initialize internal data structures before a simulation begins.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 54 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::TimeSolver::get_solution_history(), libMesh::libmesh_assert(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::TimeSolver::set_solution_history().
|
overridevirtualinherited |
Add adjoint vectors and old_adjoint_vectors as per the indices of QoISet.
Reimplemented from libMesh::TimeSolver.
Definition at line 68 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::add_vector(), libMesh::GHOSTED, libMesh::TimeSolver::init_adjoints(), libMesh::make_range(), and libMesh::System::n_qois().
|
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 89 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().
|
overridevirtual |
A method to compute the adjoint refinement error estimate at the current timestep.
int_{tstep_start}^{tstep_end} R(u^h,z) dt The user provides an initialized ARefEE object. Fills in an ErrorVector that contains the weighted sum of errors from all the QoIs and can be used to guide AMR. CURRENTLY ONLY SUPPORTED for Backward Euler.
Implements libMesh::AdaptiveTimeSolver.
Definition at line 400 of file twostep_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::System::get_qoi_error_estimate_value(), libMesh::QoISet::has_index(), libMesh::index_range(), libMesh::make_range(), libMesh::System::n_qois(), libMesh::AdjointRefinementEstimator::qoi_set(), and libMesh::System::set_qoi_error_estimate().
|
overridevirtual |
A method to integrate the adjoint sensitivity w.r.t a given parameter vector.
int_{tstep_start}^{tstep_end} dQ/dp dt = int_{tstep_start}^{tstep_end} ( / p) - ( R (u,z) / p ) dt The midpoint rule is used to integrate each substep
Implements libMesh::AdaptiveTimeSolver.
Definition at line 377 of file twostep_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::make_range(), libMesh::ParameterVector::size(), and libMesh::QoISet::size().
|
overridevirtual |
A method to integrate the system::QoI functionals.
Implements libMesh::AdaptiveTimeSolver.
Definition at line 341 of file twostep_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::System::get_qoi_value(), libMesh::make_range(), libMesh::System::n_qois(), and libMesh::System::set_qoi().
|
inlineinherited |
Accessor for querying whether we need to do a primal or adjoint solve.
Definition at line 277 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 194 of file unsteady_solver.h.
|
inlineoverridevirtualinherited |
Returns system.deltat if fixed timestep solver is used, the complete timestep size (sum of all substeps) if the adaptive time solver is used.
Returns the change in system.time, deltat, for the last timestep which was successfully completed. This only returns the outermost step size in the case of nested time solvers. If no time step has yet been successfully completed, then returns system.deltat.
Reimplemented from libMesh::TimeSolver.
Definition at line 91 of file adaptive_time_solver.h.
References libMesh::AdaptiveTimeSolver::completed_timestep_size.
|
overridevirtualinherited |
An implicit linear solver to use for adjoint and sensitivity problems.
Reimplemented from libMesh::TimeSolver.
Definition at line 218 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver.
|
inlinestaticinherited |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 85 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
Referenced by libMesh::LibMeshInit::~LibMeshInit().
|
overridevirtualinherited |
This method is passed on to the core_time_solver.
Implements libMesh::TimeSolver.
Definition at line 201 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver, and libMesh::libmesh_assert().
|
inherited |
Definition at line 337 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(), libMesh::NewmarkSolver::_general_residual(), and libMesh::FEMPhysics::eulerian_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 81 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
Referenced by libMesh::LibMeshInit::~LibMeshInit().
|
overridevirtualinherited |
The reinitialization function.
This method is used to resize internal data vectors after a mesh change.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 79 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver, and libMesh::libmesh_assert().
|
overridevirtualinherited |
This method retrieves all the stored solutions at the current system.time.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 165 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver.
|
inlineinherited |
A setter for the first_adjoint_step boolean.
Needed for nested time solvers.
Definition at line 199 of file unsteady_solver.h.
References libMesh::UnsteadySolver::first_adjoint_step.
|
inlineinherited |
Definition at line 209 of file unsteady_solver.h.
References libMesh::UnsteadySolver::first_solve.
|
inlineinherited |
Accessor for setting whether we need to do a primal or adjoint solve.
Definition at line 284 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 119 of file time_solver.C.
References libMesh::SolutionHistory::clone(), and libMesh::TimeSolver::solution_history.
Referenced by libMesh::AdaptiveTimeSolver::init().
|
overridevirtualinherited |
This method is passed on to the core_time_solver.
Implements libMesh::TimeSolver.
Definition at line 191 of file adaptive_time_solver.C.
References libMesh::AdaptiveTimeSolver::core_time_solver, and libMesh::libmesh_assert().
|
overridevirtual |
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.
Implements libMesh::AdaptiveTimeSolver.
Definition at line 54 of file twostep_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::calculate_norm(), libMesh::NumericVector< T >::clone(), libMesh::AdaptiveTimeSolver::completed_timestep_size, libMesh::AdaptiveTimeSolver::core_time_solver, libMesh::DifferentiableSystem::deltat, libMesh::UnsteadySolver::first_solve, libMesh::System::get_vector(), libMesh::AdaptiveTimeSolver::global_tolerance, libMesh::TimeSolver::last_deltat, libMesh::libmesh_assert(), libMesh::AdaptiveTimeSolver::max_deltat, libMesh::AdaptiveTimeSolver::max_growth, libMesh::AdaptiveTimeSolver::min_deltat, libMesh::out, libMesh::Utility::pow(), libMesh::TimeSolver::quiet, libMesh::Real, libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure, libMesh::System::solution, libMesh::AdaptiveTimeSolver::target_tolerance, libMesh::System::time, and libMesh::AdaptiveTimeSolver::upper_tolerance.
|
inlineinherited |
Definition at line 210 of file time_solver.h.
References libMesh::TimeSolver::_system.
Referenced by libMesh::TimeSolver::adjoint_solve(), libMesh::TimeSolver::reinit(), and libMesh::TimeSolver::solve().
|
inlineinherited |
Definition at line 215 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.
|
inherited |
Definition at line 361 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(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.
|
staticprotectedinherited |
Actually holds the data.
Definition at line 124 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info().
|
protectedinherited |
An implicit linear or nonlinear solver to use at each timestep.
Definition at line 302 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 143 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 307 of file time_solver.h.
Referenced by libMesh::TimeSolver::linear_solver(), and libMesh::TimeSolver::reinit().
|
staticprotectedinherited |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 137 of file reference_counter.h.
|
staticprotectedinherited |
The number of objects.
Print the reference count information when the number returns to 0.
Definition at line 132 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 312 of file time_solver.h.
Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::AdaptiveTimeSolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::adjoint_advance_timestep(), adjoint_solve(), libMesh::UnsteadySolver::adjoint_solve(), libMesh::TimeSolver::adjoint_solve(), libMesh::NewmarkSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::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::UnsteadySolver::init_adjoints(), libMesh::TimeSolver::init_adjoints(), libMesh::SecondOrderUnsteadySolver::init_data(), libMesh::UnsteadySolver::init_data(), libMesh::TimeSolver::init_data(), libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), integrate_adjoint_refinement_error_estimate(), libMesh::EulerSolver::integrate_adjoint_refinement_error_estimate(), integrate_adjoint_sensitivity(), libMesh::SteadySolver::integrate_adjoint_sensitivity(), libMesh::UnsteadySolver::integrate_adjoint_sensitivity(), libMesh::Euler2Solver::integrate_qoi_timestep(), integrate_qoi_timestep(), libMesh::EulerSolver::integrate_qoi_timestep(), libMesh::SteadySolver::integrate_qoi_timestep(), 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(), solve(), libMesh::UnsteadySolver::solve(), libMesh::EigenTimeSolver::solve(), libMesh::TimeSolver::system(), and libMesh::UnsteadySolver::update().
|
inherited |
The adaptive time solver's have two notions of deltat.
The deltat the solver ended up using for the completed timestep. And the deltat the solver determined would be workable for the coming timestep. The latter gets set as system.deltat. We need a variable to save the deltat used for the completed timestep.
Definition at line 205 of file adaptive_time_solver.h.
Referenced by adjoint_solve(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::last_completed_timestep_size(), and solve().
|
inherited |
Error calculations are done in this norm, DISCRETE_L2 by default.
Definition at line 135 of file adaptive_time_solver.h.
Referenced by libMesh::AdaptiveTimeSolver::calculate_norm().
|
inherited |
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 143 of file adaptive_time_solver.h.
|
inherited |
This object is used to take timesteps.
Definition at line 130 of file adaptive_time_solver.h.
Referenced by libMesh::AdaptiveTimeSolver::adjoint_advance_timestep(), adjoint_solve(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::diff_solver(), libMesh::AdaptiveTimeSolver::element_residual(), libMesh::AdaptiveTimeSolver::error_order(), libMesh::AdaptiveTimeSolver::init(), integrate_adjoint_refinement_error_estimate(), integrate_adjoint_sensitivity(), integrate_qoi_timestep(), libMesh::AdaptiveTimeSolver::linear_solver(), libMesh::AdaptiveTimeSolver::nonlocal_residual(), libMesh::AdaptiveTimeSolver::reinit(), libMesh::AdaptiveTimeSolver::retrieve_timestep(), libMesh::AdaptiveTimeSolver::side_residual(), solve(), and 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 226 of file unsteady_solver.h.
Referenced by libMesh::AdaptiveTimeSolver::adjoint_advance_timestep(), and libMesh::UnsteadySolver::set_first_adjoint_step().
|
protectedinherited |
A bool that will be true the first time solve() is called, and false thereafter.
Definition at line 220 of file unsteady_solver.h.
Referenced by libMesh::NewmarkSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::UnsteadySolver::set_first_solve(), solve(), and libMesh::UnsteadySolver::solve().
|
inherited |
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 222 of file adaptive_time_solver.h.
Referenced by solve().
|
protectedinherited |
The deltat for the last completed timestep before the current one.
Definition at line 332 of file time_solver.h.
Referenced by adjoint_solve(), libMesh::UnsteadySolver::adjoint_solve(), libMesh::TimeSolver::last_completed_timestep_size(), solve(), and libMesh::UnsteadySolver::solve().
|
protectedinherited |
We will need to move the system.time around to ensure that residuals are built with the right deltat and the right time.
Definition at line 237 of file unsteady_solver.h.
Referenced by libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), and libMesh::EulerSolver::integrate_adjoint_refinement_error_estimate().
|
inherited |
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 183 of file adaptive_time_solver.h.
Referenced by solve().
|
inherited |
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 197 of file adaptive_time_solver.h.
Referenced by solve().
|
inherited |
Do not allow the adaptive time solver to select deltat < min_deltat.
The default value is 0.0.
Definition at line 189 of file adaptive_time_solver.h.
Referenced by solve().
|
protectedinherited |
Definition at line 238 of file unsteady_solver.h.
Referenced by libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), and libMesh::EulerSolver::integrate_adjoint_refinement_error_estimate().
|
protectedinherited |
A vector of pointers to vectors holding the adjoint solution at the last time step.
Definition at line 231 of file unsteady_solver.h.
Referenced by libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), libMesh::EulerSolver::integrate_adjoint_refinement_error_estimate(), and libMesh::UnsteadySolver::UnsteadySolver().
|
inherited |
Serial vector of _system.get_vector("_old_nonlinear_solution") This is a shared_ptr so that it can be shared between different derived class instances, as in e.g.
Definition at line 178 of file unsteady_solver.h.
Referenced by libMesh::AdaptiveTimeSolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::AdaptiveTimeSolver::init(), libMesh::UnsteadySolver::init_data(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::UnsteadySolver::reinit(), libMesh::UnsteadySolver::retrieve_timestep(), and libMesh::UnsteadySolver::update().
|
inherited |
Print extra debugging information if quiet == false.
Definition at line 230 of file time_solver.h.
Referenced by 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 259 of file time_solver.h.
Referenced by 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 319 of file time_solver.h.
Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::TimeSolver::get_solution_history(), libMesh::UnsteadySolver::retrieve_timestep(), and libMesh::TimeSolver::set_solution_history().
|
inherited |
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 160 of file adaptive_time_solver.h.
Referenced by solve().
|
inherited |
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 177 of file adaptive_time_solver.h.
Referenced by solve().