libMesh
|
This class defines a theta-method Euler (defaulting to Backward Euler with theta = 1.0) solver to handle time integration of DifferentiableSystems. More...
#include <euler_solver.h>
Public Types | |
typedef FirstOrderUnsteadySolver | Parent |
The parent class. More... | |
typedef DifferentiableSystem | sys_type |
The type of system. More... | |
Public Member Functions | |
EulerSolver (sys_type &s) | |
Constructor. More... | |
virtual | ~EulerSolver () |
Destructor. More... | |
virtual Real | error_order () const override |
Error convergence order: 2 for Crank-Nicolson, 1 otherwise. More... | |
virtual bool | element_residual (bool request_jacobian, DiffContext &) override |
This method uses the DifferentiablePhysics' element_time_derivative() and element_constraint() to build a full residual on an element. More... | |
virtual bool | side_residual (bool request_jacobian, DiffContext &) override |
This method uses the DifferentiablePhysics' side_time_derivative() and side_constraint() to build a full residual on an element's side. More... | |
virtual bool | nonlocal_residual (bool request_jacobian, DiffContext &) override |
This method uses the DifferentiablePhysics' nonlocal_time_derivative() and nonlocal_constraint() to build a full residual for non-local terms. More... | |
virtual void | integrate_qoi_timestep () override |
A method to integrate the system::QoI functionals. 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 unsigned int | time_order () const override |
virtual void | init () override |
The initialization function. More... | |
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... | |
virtual void | reinit () override |
The reinitialization function. More... | |
virtual void | solve () override |
This method solves for the solution at the next timestep. More... | |
virtual void | advance_timestep () override |
This method advances the solution to the next timestep, after a solve() has been performed. More... | |
void | update () |
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 | 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 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... | |
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 () |
virtual std::unique_ptr< DiffSolver > & | diff_solver () |
An implicit linear or nonlinear solver to use at each timestep. More... | |
virtual std::unique_ptr< LinearSolver< Number > > & | linear_solver () |
An implicit linear solver to use for adjoint and sensitivity problems. More... | |
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... | |
virtual Real | last_completed_timestep_size () |
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... | |
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 | |
Real | theta |
The value for the theta method to employ: 1.0 corresponds to backwards Euler, 0.0 corresponds to forwards Euler, 0.5 corresponds to Crank-Nicolson. 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 bool | _general_residual (bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit, bool compute_second_order_eqns) |
This method is the underlying implementation of the public residual methods. 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 defines a theta-method Euler (defaulting to Backward Euler with theta = 1.0) solver to handle time integration of DifferentiableSystems.
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 43 of file euler_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 49 of file euler_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 32 of file euler_solver.C.
|
virtualdefault |
Destructor.
|
protectedvirtual |
This method is the underlying implementation of the public residual methods.
Definition at line 102 of file euler_solver.C.
References libMesh::TimeSolver::_system, libMesh::DenseVector< T >::add(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::DifferentiableSystem::deltat, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DiffContext::fixed_solution_derivative, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_solution(), libMesh::DiffContext::get_elem_solution_rate(), libMesh::DifferentiableSystem::get_physics(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::FirstOrderUnsteadySolver::prepare_accel(), libMesh::DenseVector< T >::size(), libMesh::DenseVector< T >::swap(), libMesh::DenseMatrix< T >::swap(), theta, and libMesh::System::use_fixed_solution.
Referenced by element_residual(), nonlocal_residual(), and side_residual().
|
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::AdaptiveTimeSolver, and libMesh::NewmarkSolver.
Definition at line 239 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, 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, libMesh::TimeSolver::solution_history, and libMesh::System::time.
|
overridevirtualinherited |
This method solves for the adjoint solution at the next adjoint timestep (or a steady state adjoint solve)
Reimplemented from libMesh::TimeSolver.
Reimplemented in libMesh::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.
Definition at line 228 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, and libMesh::TimeSolver::last_deltat.
|
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::TimeSolver.
Reimplemented in libMesh::AdaptiveTimeSolver, and libMesh::NewmarkSolver.
Definition at line 192 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, 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, libMesh::TimeSolver::solution_history, and libMesh::System::time.
Referenced by libMesh::NewmarkSolver::advance_timestep(), and libMesh::UnsteadySolver::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.
|
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 _general_residual(), libMesh::Euler2Solver::_general_residual(), element_residual(), libMesh::Euler2Solver::element_residual(), nonlocal_residual(), and libMesh::Euler2Solver::nonlocal_residual().
|
inlinevirtualinherited |
An implicit linear or nonlinear solver to use at each timestep.
Reimplemented in libMesh::AdaptiveTimeSolver.
Definition at line 220 of file time_solver.h.
References libMesh::TimeSolver::_diff_solver.
Referenced by libMesh::TimeSolver::adjoint_solve(), adjust_linear_solvers(), libMesh::TimeSolver::init(), libMesh::TimeSolver::init_data(), libMesh::TimeSolver::reinit(), and libMesh::TimeSolver::solve().
|
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.
|
overridevirtual |
This method uses the DifferentiablePhysics' element_time_derivative() and element_constraint() to build a full residual on an element.
What combination it uses will depend on theta.
Implements libMesh::TimeSolver.
Definition at line 53 of file euler_solver.C.
References libMesh::DifferentiablePhysics::_eulerian_time_deriv(), _general_residual(), libMesh::TimeSolver::_system, libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::DifferentiablePhysics::damping_residual(), libMesh::DiffContext::elem_reinit(), libMesh::DifferentiablePhysics::element_constraint(), libMesh::DifferentiableSystem::get_physics(), libMesh::DifferentiablePhysics::have_second_order_vars(), and libMesh::DifferentiablePhysics::mass_residual().
|
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.
|
overridevirtual |
Error convergence order: 2 for Crank-Nicolson, 1 otherwise.
Implements libMesh::UnsteadySolver.
Definition at line 43 of file euler_solver.C.
References theta.
|
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::TimeSolver.
Reimplemented in libMesh::AdaptiveTimeSolver, and libMesh::SecondOrderUnsteadySolver.
Definition at line 60 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::add_vector(), and libMesh::TimeSolver::init().
Referenced by libMesh::SecondOrderUnsteadySolver::init().
|
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. The integration scheme used is consistent with the theta used for the unsteady solution. CURRENTLY ONLY SUPPORTED for Backward Euler.
Implements libMesh::FirstOrderUnsteadySolver.
Definition at line 254 of file euler_solver.C.
References libMesh::TimeSolver::_system, std::abs(), libMesh::NumericVector< T >::clone(), libMesh::DifferentiableSystem::deltat, libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::System::get_adjoint_solution(), libMesh::AdjointRefinementEstimator::get_global_QoI_error_estimate(), libMesh::QoISet::has_index(), libMesh::index_range(), libMesh::UnsteadySolver::last_step_deltat, libMesh::make_range(), libMesh::System::n_qois(), libMesh::UnsteadySolver::next_step_deltat, libMesh::UnsteadySolver::old_adjoints, libMesh::AdjointRefinementEstimator::qoi_set(), libMesh::Real, libMesh::UnsteadySolver::retrieve_timestep(), libMesh::System::set_qoi_error_estimate(), std::sqrt(), libMesh::NumericVector< T >::swap(), theta, libMesh::System::time, libMesh::TOLERANCE, and libMesh::System::update().
|
overridevirtualinherited |
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 trapezoidal rule is used to numerically integrate the timestep.
Reimplemented from libMesh::TimeSolver.
Reimplemented in libMesh::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.
Definition at line 280 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::ImplicitSystem::adjoint_qoi_parameter_sensitivity(), libMesh::DifferentiableSystem::deltat, libMesh::make_range(), libMesh::Real, libMesh::System::remove_vector(), libMesh::UnsteadySolver::retrieve_timestep(), libMesh::ParameterVector::size(), libMesh::QoISet::size(), and libMesh::System::time.
|
overridevirtual |
A method to integrate the system::QoI functionals.
Implements libMesh::FirstOrderUnsteadySolver.
Definition at line 196 of file euler_solver.C.
References libMesh::TimeSolver::_system, libMesh::ExplicitSystem::assemble_qoi(), libMesh::DifferentiableSystem::deltat, libMesh::System::get_qoi_value(), libMesh::make_range(), libMesh::System::n_qois(), libMesh::UnsteadySolver::retrieve_timestep(), libMesh::System::set_qoi(), theta, and libMesh::System::time.
|
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.
|
virtualinherited |
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 in libMesh::AdaptiveTimeSolver.
Definition at line 160 of file time_solver.C.
References libMesh::TimeSolver::last_deltat.
|
inlinevirtualinherited |
An implicit linear solver to use for adjoint and sensitivity problems.
Reimplemented in libMesh::AdaptiveTimeSolver.
Definition at line 225 of file time_solver.h.
References libMesh::TimeSolver::_linear_solver.
Referenced by libMesh::TimeSolver::init(), libMesh::TimeSolver::init_data(), and libMesh::TimeSolver::reinit().
|
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().
|
overridevirtual |
This method uses the DifferentiablePhysics' nonlocal_time_derivative() and nonlocal_constraint() to build a full residual for non-local terms.
What combination it uses will depend on theta.
Implements libMesh::TimeSolver.
Definition at line 85 of file euler_solver.C.
References _general_residual(), libMesh::TimeSolver::_system, libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::DifferentiableSystem::have_second_order_scalar_vars(), libMesh::DifferentiablePhysics::nonlocal_constraint(), libMesh::DifferentiablePhysics::nonlocal_damping_residual(), libMesh::DifferentiablePhysics::nonlocal_mass_residual(), libMesh::DiffContext::nonlocal_reinit(), and libMesh::DifferentiablePhysics::nonlocal_time_derivative().
|
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 _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 _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::TimeSolver.
Reimplemented in libMesh::SecondOrderUnsteadySolver, and libMesh::AdaptiveTimeSolver.
Definition at line 104 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::GHOSTED, libMesh::NumericVector< T >::localize(), libMesh::System::n_dofs(), libMesh::System::n_local_dofs(), libMesh::UnsteadySolver::old_local_nonlinear_solution, libMesh::TimeSolver::reinit(), and libMesh::SERIAL.
Referenced by libMesh::SecondOrderUnsteadySolver::reinit().
|
overridevirtualinherited |
This method retrieves all the stored solutions at the current system.time.
Reimplemented from libMesh::TimeSolver.
Reimplemented in libMesh::SecondOrderUnsteadySolver, and libMesh::AdaptiveTimeSolver.
Definition at line 269 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, libMesh::TimeSolver::solution_history, and libMesh::System::time.
Referenced by libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), integrate_adjoint_refinement_error_estimate(), libMesh::UnsteadySolver::integrate_adjoint_sensitivity(), libMesh::Euler2Solver::integrate_qoi_timestep(), and integrate_qoi_timestep().
|
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().
|
overridevirtual |
This method uses the DifferentiablePhysics' side_time_derivative() and side_constraint() to build a full residual on an element's side.
What combination it uses will depend on theta.
Implements libMesh::TimeSolver.
Definition at line 70 of file euler_solver.C.
References _general_residual(), libMesh::DiffContext::elem_side_reinit(), libMesh::DifferentiablePhysics::side_constraint(), libMesh::DifferentiablePhysics::side_damping_residual(), libMesh::DifferentiablePhysics::side_mass_residual(), and libMesh::DifferentiablePhysics::side_time_derivative().
|
overridevirtualinherited |
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::TimeSolver.
Reimplemented in libMesh::NewmarkSolver, libMesh::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.
Definition at line 127 of file unsteady_solver.C.
References libMesh::TimeSolver::_diff_solver, libMesh::TimeSolver::_system, libMesh::UnsteadySolver::advance_timestep(), libMesh::DifferentiableSystem::deltat, libMesh::DiffSolver::DIVERGED_BACKTRACKING_FAILURE, libMesh::DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS, libMesh::UnsteadySolver::first_solve, libMesh::TimeSolver::last_deltat, libMesh::out, libMesh::TimeSolver::quiet, and libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure.
Referenced by libMesh::NewmarkSolver::solve().
|
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 _general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::AdaptiveTimeSolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::TwostepTimeSolver::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(), 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(), libMesh::TwostepTimeSolver::integrate_adjoint_refinement_error_estimate(), integrate_adjoint_refinement_error_estimate(), libMesh::TwostepTimeSolver::integrate_adjoint_sensitivity(), libMesh::SteadySolver::integrate_adjoint_sensitivity(), libMesh::UnsteadySolver::integrate_adjoint_sensitivity(), libMesh::Euler2Solver::integrate_qoi_timestep(), libMesh::TwostepTimeSolver::integrate_qoi_timestep(), integrate_qoi_timestep(), libMesh::SteadySolver::integrate_qoi_timestep(), 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(), libMesh::TimeSolver::system(), and libMesh::UnsteadySolver::update().
|
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(), libMesh::TwostepTimeSolver::solve(), and libMesh::UnsteadySolver::solve().
|
protectedinherited |
The deltat for the last completed timestep before the current one.
Definition at line 332 of file time_solver.h.
Referenced by libMesh::TwostepTimeSolver::adjoint_solve(), libMesh::UnsteadySolver::adjoint_solve(), libMesh::TimeSolver::last_completed_timestep_size(), libMesh::TwostepTimeSolver::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 integrate_adjoint_refinement_error_estimate().
|
protectedinherited |
Definition at line 238 of file unsteady_solver.h.
Referenced by libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), and 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(), 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 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 259 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 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().
Real libMesh::EulerSolver::theta |
The value for the theta method to employ: 1.0 corresponds to backwards Euler, 0.0 corresponds to forwards Euler, 0.5 corresponds to Crank-Nicolson.
Definition at line 117 of file euler_solver.h.
Referenced by _general_residual(), error_order(), integrate_adjoint_refinement_error_estimate(), and integrate_qoi_timestep().