PolyFEM
Loading...
Searching...
No Matches
ImplicitEuler.cpp
Go to the documentation of this file.
1#include "ImplicitEuler.hpp"
2
4{
5 void ImplicitEuler::update_quantities(const Eigen::VectorXd &x)
6 {
7 const Eigen::VectorXd v = compute_velocity(x);
9 set_v_prev(v);
11 }
12
13 Eigen::VectorXd ImplicitEuler::x_tilde() const
14 {
15 return x_prev() + dt() * v_prev();
16 }
17
18 Eigen::VectorXd ImplicitEuler::compute_velocity(const Eigen::VectorXd &x) const
19 {
20 return (x - x_prev()) / dt();
21 }
22
23 Eigen::VectorXd ImplicitEuler::compute_acceleration(const Eigen::VectorXd &v) const
24 {
25 return (v - v_prev()) / dt();
26 }
27
29 {
30 return dt() * dt();
31 }
32
33 double ImplicitEuler::dv_dx(const unsigned prev_ti) const
34 {
35 if (prev_ti > 1)
36 return 0;
37 return (prev_ti == 0 ? 1 : -1) / dt();
38 }
39} // namespace polyfem::time_integrator
int 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).
void update_quantities(const Eigen::VectorXd &x) override
Update the time integration quantities (i.e., , , and ).
double dv_dx(const unsigned prev_ti=0) const override
Compute the derivative of the velocity with respect to the solution.
double acceleration_scaling() const override
Compute the acceleration scaling used to scale forces when integrating a second order ODE.
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)...
Eigen::VectorXd x_tilde() const override
Compute the predicted solution to be used in the inertia term .
const Eigen::VectorXd & v_prev() const
Get the most recent previous velocity value.
const Eigen::VectorXd & x_prev() const
Get the most recent previous solution value.
void set_a_prev(const Eigen::VectorXd &a_prev)
Convenience functions for setting the most recent previous acceleration.
void set_v_prev(const Eigen::VectorXd &v_prev)
Convenience functions for setting the most recent previous velocity.
const double & dt() const
Access the time step size.
void set_x_prev(const Eigen::VectorXd &x_prev)
Convenience functions for setting the most recent previous solution.