PolyFEM
|
Implicit time integrator of a second order ODE (equivently a system of coupled first order ODEs). More...
#include <ImplicitTimeIntegrator.hpp>
Public Member Functions | |
ImplicitTimeIntegrator () | |
virtual | ~ImplicitTimeIntegrator ()=default |
virtual void | set_parameters (const json ¶ms) |
Set the time integrator parameters from a json object. | |
virtual void | init (const Eigen::MatrixXd &x_prevs, const Eigen::MatrixXd &v_prevs, const Eigen::MatrixXd &a_prevs, double dt) |
Initialize the time integrator with the previous values for \(x\), \(v\), and \(a\). | |
virtual void | update_quantities (const Eigen::VectorXd &x)=0 |
Update the time integration quantities (i.e., \(x\), \(v\), and \(a\)). | |
virtual Eigen::VectorXd | x_tilde () const =0 |
Compute the predicted solution to be used in the inertia term \((x-\tilde{x})^TM(x-\tilde{x})\). | |
virtual Eigen::VectorXd | compute_velocity (const Eigen::VectorXd &x) const =0 |
Compute the current velocity given the current solution and using the stored previous solution(s). | |
virtual Eigen::VectorXd | compute_acceleration (const Eigen::VectorXd &v) const =0 |
Compute the current acceleration given the current velocity and using the stored previous velocity(s). | |
virtual double | acceleration_scaling () const =0 |
Compute the acceleration scaling used to scale forces when integrating a second order ODE. | |
virtual double | dv_dx (const unsigned prev_ti=0) const =0 |
Compute the derivative of the velocity with respect to the solution. | |
const double & | dt () const |
Access the time step size. | |
virtual void | save_state (const std::string &state_path) const |
Save the values of \(x\), \(v\), and \(a\). | |
const Eigen::VectorXd & | x_prev () const |
Get the most recent previous solution value. | |
const Eigen::VectorXd & | v_prev () const |
Get the most recent previous velocity value. | |
const Eigen::VectorXd & | a_prev () const |
Get the most recent previous acceleration value. | |
const std::deque< Eigen::VectorXd > & | x_prevs () const |
Get the (relevant) history of previous solution value. | |
const std::deque< Eigen::VectorXd > & | v_prevs () const |
Get the (relevant) history of previous velocity value. | |
const std::deque< Eigen::VectorXd > & | a_prevs () const |
Get the (relevant) history of previous acceleration value. | |
int | steps () const |
Get the current number of steps to use for integration. | |
Static Public Member Functions | |
static std::shared_ptr< ImplicitTimeIntegrator > | construct_time_integrator (const json ¶ms) |
Factory method for constructing implicit time integrators from the name of the integrator. | |
static const std::vector< std::string > & | get_time_integrator_names () |
Get a vector of the names of possible ImplicitTimeIntegrators. | |
Protected Member Functions | |
virtual int | max_steps () const |
Get the maximum number of steps to use for integration. | |
void | set_x_prev (const Eigen::VectorXd &x_prev) |
Convenience functions for setting the most recent previous solution. | |
void | set_v_prev (const Eigen::VectorXd &v_prev) |
Convenience functions for setting the most recent previous velocity. | |
void | set_a_prev (const Eigen::VectorXd &a_prev) |
Convenience functions for setting the most recent previous acceleration. | |
Protected Attributes | |
double | dt_ = 1 |
Time step size. | |
std::deque< Eigen::VectorXd > | x_prevs_ |
Store the necessary previous values of the solution for single or multi-step integration. | |
std::deque< Eigen::VectorXd > | v_prevs_ |
Store the necessary previous values of the velocity for single or multi-step integration. | |
std::deque< Eigen::VectorXd > | a_prevs_ |
Store the necessary previous values of the acceleration for single or multi-step integration. | |
Implicit time integrator of a second order ODE (equivently a system of coupled first order ODEs).
Definition at line 14 of file ImplicitTimeIntegrator.hpp.
|
inline |
Definition at line 17 of file ImplicitTimeIntegrator.hpp.
|
virtualdefault |
|
inline |
Get the most recent previous acceleration value.
Definition at line 79 of file ImplicitTimeIntegrator.hpp.
References a_prevs_.
Referenced by polyfem::time_integrator::ImplicitNewmark::compute_acceleration(), polyfem::time_integrator::ImplicitNewmark::compute_velocity(), set_a_prev(), and polyfem::time_integrator::ImplicitNewmark::x_tilde().
|
inline |
Get the (relevant) history of previous acceleration value.
Definition at line 86 of file ImplicitTimeIntegrator.hpp.
References a_prevs_.
Referenced by init(), and save_state().
|
pure virtual |
Compute the acceleration scaling used to scale forces when integrating a second order ODE.
Implemented in polyfem::time_integrator::BDF, polyfem::time_integrator::ImplicitEuler, and polyfem::time_integrator::ImplicitNewmark.
|
pure virtual |
Compute the current acceleration given the current velocity and using the stored previous velocity(s).
v | current velocity |
Implemented in polyfem::time_integrator::BDF, polyfem::time_integrator::ImplicitEuler, and polyfem::time_integrator::ImplicitNewmark.
|
pure virtual |
Compute the current velocity given the current solution and using the stored previous solution(s).
x | current solution |
Implemented in polyfem::time_integrator::BDF, polyfem::time_integrator::ImplicitEuler, and polyfem::time_integrator::ImplicitNewmark.
Referenced by polyfem::solver::RayleighDampingForm::first_derivative_unweighted(), and polyfem::solver::RayleighDampingForm::value_unweighted().
|
static |
Factory method for constructing implicit time integrators from the name of the integrator.
name | name of the type of ImplicitTimeIntegrator to construct |
Definition at line 66 of file ImplicitTimeIntegrator.cpp.
References polyfem::logger(), and polyfem::utils::StringUtils::startswith().
Referenced by polyfem::State::init_linear_solve(), polyfem::State::init_nonlinear_tensor_solve(), polyfem::mesh::LocalRelaxationData< M >::init_solve_data(), and polyfem::State::solve_transient_linear().
|
inline |
Access the time step size.
Definition at line 59 of file ImplicitTimeIntegrator.hpp.
References dt_.
Referenced by polyfem::time_integrator::BDF::acceleration_scaling(), polyfem::time_integrator::ImplicitEuler::acceleration_scaling(), polyfem::time_integrator::ImplicitNewmark::acceleration_scaling(), polyfem::time_integrator::BDF::beta_dt(), polyfem::time_integrator::ImplicitEuler::compute_acceleration(), polyfem::time_integrator::ImplicitNewmark::compute_acceleration(), polyfem::time_integrator::ImplicitEuler::compute_velocity(), polyfem::time_integrator::ImplicitNewmark::compute_velocity(), polyfem::time_integrator::ImplicitNewmark::da_dx(), polyfem::time_integrator::ImplicitEuler::dv_dx(), polyfem::time_integrator::ImplicitNewmark::dv_dx(), init(), polyfem::solver::RayleighDampingForm::stiffness(), polyfem::time_integrator::BDF::x_tilde(), polyfem::time_integrator::ImplicitEuler::x_tilde(), and polyfem::time_integrator::ImplicitNewmark::x_tilde().
|
pure virtual |
Compute the derivative of the velocity with respect to the solution.
prev_ti | index of the previous solution to use (0 -> current; 1 -> previous; 2 -> second previous; etc.) |
Implemented in polyfem::time_integrator::BDF, polyfem::time_integrator::ImplicitEuler, and polyfem::time_integrator::ImplicitNewmark.
Referenced by polyfem::solver::RayleighDampingForm::second_derivative_unweighted(), and polyfem::solver::RayleighDampingForm::value_unweighted().
|
static |
Get a vector of the names of possible ImplicitTimeIntegrators.
Definition at line 95 of file ImplicitTimeIntegrator.cpp.
|
virtual |
Initialize the time integrator with the previous values for \(x\), \(v\), and \(a\).
x_prev | previous value(s) for the solution |
v_prev | previous value(s) for the velocity |
a_prev | previous value(s) for the acceleration |
dt | time step size |
Definition at line 18 of file ImplicitTimeIntegrator.cpp.
References a_prevs(), a_prevs_, dt(), dt_, max_steps(), v_prevs(), v_prevs_, x_prevs(), and x_prevs_.
Referenced by polyfem::State::solve_transient_navier_stokes().
|
inlineprotectedvirtual |
Get the maximum number of steps to use for integration.
Reimplemented in polyfem::time_integrator::BDF.
Definition at line 93 of file ImplicitTimeIntegrator.hpp.
Referenced by init().
|
virtual |
Save the values of \(x\), \(v\), and \(a\).
state_path | path for the output file containing \(x, v, a\) as hdf5 |
Definition at line 44 of file ImplicitTimeIntegrator.cpp.
References a_prevs(), v_prevs(), polyfem::io::write_matrix(), x_prev(), and x_prevs().
|
inlineprotected |
Convenience functions for setting the most recent previous acceleration.
Definition at line 111 of file ImplicitTimeIntegrator.hpp.
References a_prev(), and a_prevs_.
Referenced by polyfem::time_integrator::ImplicitEuler::update_quantities(), and polyfem::time_integrator::ImplicitNewmark::update_quantities().
|
inlinevirtual |
Set the time integrator parameters from a json object.
params | json containing parameters specific to each time integrator |
Reimplemented in polyfem::time_integrator::BDF, and polyfem::time_integrator::ImplicitNewmark.
Definition at line 22 of file ImplicitTimeIntegrator.hpp.
|
inlineprotected |
Convenience functions for setting the most recent previous velocity.
Definition at line 109 of file ImplicitTimeIntegrator.hpp.
References v_prev(), and v_prevs_.
Referenced by polyfem::time_integrator::ImplicitEuler::update_quantities(), and polyfem::time_integrator::ImplicitNewmark::update_quantities().
|
inlineprotected |
Convenience functions for setting the most recent previous solution.
Definition at line 107 of file ImplicitTimeIntegrator.hpp.
References x_prev(), and x_prevs_.
Referenced by polyfem::time_integrator::ImplicitEuler::update_quantities(), and polyfem::time_integrator::ImplicitNewmark::update_quantities().
|
inline |
Get the current number of steps to use for integration.
Definition at line 89 of file ImplicitTimeIntegrator.hpp.
References x_prevs_.
Referenced by polyfem::time_integrator::BDF::acceleration_scaling(), polyfem::time_integrator::BDF::beta_dt(), polyfem::time_integrator::BDF::dv_dx(), polyfem::time_integrator::BDF::update_quantities(), polyfem::time_integrator::BDF::weighted_sum_v_prevs(), polyfem::time_integrator::BDF::weighted_sum_x_prevs(), and polyfem::time_integrator::BDF::x_tilde().
|
pure virtual |
Update the time integration quantities (i.e., \(x\), \(v\), and \(a\)).
x | new solution vector |
Implemented in polyfem::time_integrator::BDF, polyfem::time_integrator::ImplicitEuler, and polyfem::time_integrator::ImplicitNewmark.
|
inline |
Get the most recent previous velocity value.
Definition at line 77 of file ImplicitTimeIntegrator.hpp.
References v_prevs_.
Referenced by polyfem::State::cache_transient_adjoint_quantities(), polyfem::time_integrator::ImplicitEuler::compute_acceleration(), polyfem::time_integrator::ImplicitNewmark::compute_acceleration(), polyfem::time_integrator::ImplicitNewmark::compute_velocity(), set_v_prev(), polyfem::time_integrator::BDF::weighted_sum_v_prevs(), polyfem::time_integrator::ImplicitEuler::x_tilde(), and polyfem::time_integrator::ImplicitNewmark::x_tilde().
|
inline |
Get the (relevant) history of previous velocity value.
Definition at line 84 of file ImplicitTimeIntegrator.hpp.
References v_prevs_.
Referenced by init(), and save_state().
|
inline |
Get the most recent previous solution value.
Definition at line 75 of file ImplicitTimeIntegrator.hpp.
References x_prevs_.
Referenced by polyfem::time_integrator::ImplicitEuler::compute_velocity(), polyfem::time_integrator::ImplicitNewmark::compute_velocity(), save_state(), set_x_prev(), polyfem::time_integrator::BDF::weighted_sum_x_prevs(), polyfem::time_integrator::ImplicitEuler::x_tilde(), and polyfem::time_integrator::ImplicitNewmark::x_tilde().
|
inline |
Get the (relevant) history of previous solution value.
Definition at line 82 of file ImplicitTimeIntegrator.hpp.
References x_prevs_.
Referenced by init(), and save_state().
|
pure virtual |
Compute the predicted solution to be used in the inertia term \((x-\tilde{x})^TM(x-\tilde{x})\).
Implemented in polyfem::time_integrator::BDF, polyfem::time_integrator::ImplicitEuler, and polyfem::time_integrator::ImplicitNewmark.
Referenced by polyfem::solver::InertiaForm::first_derivative_unweighted(), and polyfem::solver::InertiaForm::value_unweighted().
|
protected |
Store the necessary previous values of the acceleration for single or multi-step integration.
Definition at line 104 of file ImplicitTimeIntegrator.hpp.
Referenced by a_prev(), a_prevs(), init(), set_a_prev(), and polyfem::time_integrator::BDF::update_quantities().
|
protected |
Time step size.
Default of one for static sims, this should be set using init().
Definition at line 97 of file ImplicitTimeIntegrator.hpp.
|
protected |
Store the necessary previous values of the velocity for single or multi-step integration.
Definition at line 102 of file ImplicitTimeIntegrator.hpp.
Referenced by init(), set_v_prev(), polyfem::time_integrator::BDF::update_quantities(), v_prev(), v_prevs(), and polyfem::time_integrator::BDF::weighted_sum_v_prevs().
|
protected |
Store the necessary previous values of the solution for single or multi-step integration.
Definition at line 100 of file ImplicitTimeIntegrator.hpp.
Referenced by init(), set_x_prev(), steps(), polyfem::time_integrator::BDF::update_quantities(), polyfem::time_integrator::BDF::weighted_sum_x_prevs(), x_prev(), and x_prevs().