|
PolyFEM
|
Backward Differential Formulas. More...
#include <BDF.hpp>
Public Member Functions | |
| BDF (const int order=1) | |
| void | set_parameters (const json ¶ms) override |
| Set the number of steps parameters from a json object. | |
| void | update_quantities (const Eigen::VectorXd &x) override |
| Update the time integration quantities (i.e., \(x\), \(v\), and \(a\)). | |
| Eigen::VectorXd | x_tilde () const override |
| Compute the predicted solution to be used in the inertia term \((x-\tilde{x})^TM(x-\tilde{x})\). | |
| Eigen::VectorXd | compute_velocity (const Eigen::VectorXd &x) const override |
| Compute the current velocity given the current solution and using the stored previous solution(s). | |
| Eigen::VectorXd | compute_acceleration (const Eigen::VectorXd &v) const override |
| Compute the current acceleration given the current velocity and using the stored previous velocity(s). | |
| double | acceleration_scaling () const override |
| Compute the acceleration scaling used to scale forces when integrating a second order ODE. | |
| double | dv_dx (const unsigned prev_ti=0) const override |
| Compute the derivative of the velocity with respect to the solution. | |
| double | beta_dt () const |
| Compute \(\beta\Delta t\). | |
| Eigen::VectorXd | weighted_sum_x_prevs () const |
| Compute the weighted sum of the previous solutions. | |
| Eigen::VectorXd | weighted_sum_v_prevs () const |
| Compute the weighted sum of the previous velocities. | |
Public Member Functions inherited from polyfem::time_integrator::ImplicitTimeIntegrator | |
| ImplicitTimeIntegrator () | |
| virtual | ~ImplicitTimeIntegrator ()=default |
| 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\). | |
| 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 const std::vector< double > & | alphas (const int i) |
Retrieve the alphas used for BDF with i steps. | |
| static double | betas (const int i) |
Retrieve the value of beta used for BDF with i steps. | |
Static Public Member Functions inherited from polyfem::time_integrator::ImplicitTimeIntegrator | |
| 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 | |
| int | max_steps () const override |
| Get the maximum number of steps to use for integration. | |
Protected Member Functions inherited from polyfem::time_integrator::ImplicitTimeIntegrator | |
| 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 | |
| int | max_steps_ = 1 |
| The maximum number of steps to use for integration. | |
Protected Attributes inherited from polyfem::time_integrator::ImplicitTimeIntegrator | |
| 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. | |
Backward Differential Formulas.
\[ x^{t+1} = \left(\sum_{i=0}^{n-1} \alpha_{i} x^{t-i}\right)+ \Delta t \beta v^{t+1}\newline v^{t+1} = \left(\sum_{i=0}^{n-1} \alpha_{i} v^{t-i}\right)+ \Delta t \beta a^{t+1} \]
| polyfem::time_integrator::BDF::BDF | ( | const int | order = 1 | ) |
Definition at line 7 of file BDF.cpp.
References polyfem::log_and_throw_error(), and max_steps_.
|
overridevirtual |
Compute the acceleration scaling used to scale forces when integrating a second order ODE.
\[ \beta^2 \Delta t^2 \]
Implements polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 110 of file BDF.cpp.
References betas(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), and polyfem::time_integrator::ImplicitTimeIntegrator::steps().
Referenced by polyfem::State::solve_transient_navier_stokes().
|
static |
Retrieve the alphas used for BDF with i steps.
| i | number of steps |
Definition at line 21 of file BDF.cpp.
Referenced by dv_dx(), polyfem::State::solve_transient_adjoint(), weighted_sum_v_prevs(), and weighted_sum_x_prevs().
| double polyfem::time_integrator::BDF::beta_dt | ( | ) | const |
Compute \(\beta\Delta t\).
Definition at line 125 of file BDF.cpp.
References betas(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), and polyfem::time_integrator::ImplicitTimeIntegrator::steps().
Referenced by compute_acceleration(), compute_velocity(), and dv_dx().
|
static |
Retrieve the value of beta used for BDF with i steps.
| i | number of steps |
Definition at line 35 of file BDF.cpp.
Referenced by acceleration_scaling(), beta_dt(), polyfem::State::compute_force_jacobian_prev(), polyfem::State::solve_transient_adjoint(), and x_tilde().
|
overridevirtual |
Compute the current acceleration given the current velocity and using the stored previous velocity(s).
\[ a = \frac{v - \sum_{i=0}^{n-1} \alpha_i v^{t-i}}{\beta \Delta t} \]
| v | current velocity |
Implements polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 105 of file BDF.cpp.
References beta_dt(), and weighted_sum_v_prevs().
Referenced by update_quantities().
|
overridevirtual |
Compute the current velocity given the current solution and using the stored previous solution(s).
\[ v = \frac{x - \sum_{i=0}^{n-1} \alpha_i x^{t-i}}{\beta \Delta t} \]
| x | current solution vector |
Implements polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 100 of file BDF.cpp.
References beta_dt(), weighted_sum_x_prevs(), and x.
Referenced by update_quantities().
|
overridevirtual |
Compute the derivative of the velocity with respect to the solution.
\[ \frac{\partial v}{\partial x} = \frac{1}{\beta \Delta t} \]
\[ \frac{\partial v}{\partial x^{t-i}} = \frac{-\alpha_i}{\beta \Delta t} \]
\[ \frac{\partial v}{\partial x^{t-n}} = 0 \]
| prev_ti | index of the previous solution to use (0 -> current; 1 -> previous; 2 -> second previous; etc.) |
Implements polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 116 of file BDF.cpp.
References alphas(), beta_dt(), and polyfem::time_integrator::ImplicitTimeIntegrator::steps().
|
inlineoverrideprotectedvirtual |
Get the maximum number of steps to use for integration.
Reimplemented from polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 95 of file BDF.hpp.
References max_steps_.
Referenced by update_quantities().
|
overridevirtual |
Set the number of steps parameters from a json object.
| params | json containing {"steps": 1} |
Reimplemented from polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 14 of file BDF.cpp.
References polyfem::log_and_throw_error(), and max_steps_.
Referenced by polyfem::State::solve_transient_navier_stokes().
|
overridevirtual |
Update the time integration quantities (i.e., \(x\), \(v\), and \(a\)).
| x | new solution vector |
Implements polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 75 of file BDF.cpp.
References polyfem::time_integrator::ImplicitTimeIntegrator::a_prevs_, compute_acceleration(), compute_velocity(), max_steps(), polyfem::time_integrator::ImplicitTimeIntegrator::steps(), polyfem::time_integrator::ImplicitTimeIntegrator::v_prevs_, x, and polyfem::time_integrator::ImplicitTimeIntegrator::x_prevs_.
Referenced by polyfem::State::solve_transient_navier_stokes().
| Eigen::VectorXd polyfem::time_integrator::BDF::weighted_sum_v_prevs | ( | ) | const |
Compute the weighted sum of the previous velocities.
\[ \sum_{i=0}^{n-1} \alpha_i v^{t-i} \]
Definition at line 62 of file BDF.cpp.
References alphas(), polyfem::time_integrator::ImplicitTimeIntegrator::steps(), polyfem::time_integrator::ImplicitTimeIntegrator::v_prev(), and polyfem::time_integrator::ImplicitTimeIntegrator::v_prevs_.
Referenced by polyfem::State::cache_transient_adjoint_quantities(), compute_acceleration(), and x_tilde().
| Eigen::VectorXd polyfem::time_integrator::BDF::weighted_sum_x_prevs | ( | ) | const |
Compute the weighted sum of the previous solutions.
\[ \sum_{i=0}^{n-1} \alpha_i x^{t-i} \]
Definition at line 49 of file BDF.cpp.
References alphas(), polyfem::time_integrator::ImplicitTimeIntegrator::steps(), polyfem::time_integrator::ImplicitTimeIntegrator::x_prev(), and polyfem::time_integrator::ImplicitTimeIntegrator::x_prevs_.
Referenced by compute_velocity(), polyfem::State::solve_transient_navier_stokes(), and x_tilde().
|
overridevirtual |
Compute the predicted solution to be used in the inertia term \((x-\tilde{x})^TM(x-\tilde{x})\).
\[ \tilde{x} = \left(\sum_{i=0}^{n-1} \alpha_i x^{t-i}\right) + \beta \Delta t \left(\sum_{i=0}^{n-1} \alpha_i v^{t-i}\right) \]
Implements polyfem::time_integrator::ImplicitTimeIntegrator.
Definition at line 95 of file BDF.cpp.
References betas(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), polyfem::time_integrator::ImplicitTimeIntegrator::steps(), weighted_sum_v_prevs(), and weighted_sum_x_prevs().
|
protected |
The maximum number of steps to use for integration.
Definition at line 98 of file BDF.hpp.
Referenced by BDF(), max_steps(), and set_parameters().