PolyFEM
Loading...
Searching...
No Matches
ImplicitNewmark.cpp
Go to the documentation of this file.
1#include "ImplicitNewmark.hpp"
2
4{
6 {
7 beta_ = params.at("gamma");
8 gamma_ = params.at("beta");
9 }
10
11 void ImplicitNewmark::update_quantities(const Eigen::VectorXd &x)
12 {
13 const Eigen::VectorXd v = compute_velocity(x);
15 set_v_prev(v);
17 }
18
19 Eigen::VectorXd ImplicitNewmark::x_tilde() const
20 {
21 return x_prev() + dt() * (v_prev() + dt() * (0.5 - beta()) * a_prev());
22 }
23
24 Eigen::VectorXd ImplicitNewmark::compute_velocity(const Eigen::VectorXd &x) const
25 {
26 const double c = gamma() / beta();
27 return c / dt() * (x - x_prev()) + (1 - c) * v_prev() + (1 - c / 2) * dt() * a_prev();
28 }
29
30 Eigen::VectorXd ImplicitNewmark::compute_acceleration(const Eigen::VectorXd &v) const
31 {
32 return (v - v_prev() - (1 - gamma()) * dt() * a_prev()) / (gamma() * dt());
33 }
34
36 {
37 return beta() * dt() * dt();
38 }
39
40 double ImplicitNewmark::dv_dx(const unsigned prev_ti) const
41 {
42 // if (i == n_steps - 1)
43 // throw std::runtime_error("dv_dx is not defined for the last step");
44 const double c = gamma() / beta();
45 if (prev_ti == 0)
46 return c / dt();
47 return ((prev_ti == 1 ? (-c / dt()) : 0)
48 + (1 - c) * dv_dx(prev_ti - 1)
49 + (1 - c / 2) * dt() * da_dx(prev_ti - 1));
50 }
51
52 double ImplicitNewmark::da_dx(const unsigned prev_ti) const
53 {
54 // if (i == n_steps - 1)
55 // throw std::runtime_error("da_dx is not defined for the last step");
56 if (prev_ti == 0)
57 return dv_dx(prev_ti) / (gamma() * dt());
58 return (dv_dx(prev_ti)
59 - dv_dx(prev_ti - 1)
60 - (1 - gamma()) * dt() * da_dx(prev_ti - 1))
61 / (gamma() * dt());
62 }
63} // namespace polyfem::time_integrator
int x
double beta() const
parameter for blending accelerations in the solution update.
double acceleration_scaling() const override
Compute the acceleration scaling used to scale forces when integrating a second order ODE.
double da_dx(const unsigned prev_ti=0) const
double gamma_
parameter for blending accelerations in the velocity update.
double beta_
parameter for blending accelerations in the solution update.
Eigen::VectorXd x_tilde() const override
Compute the predicted solution to be used in the inertia term .
void set_parameters(const json &params) override
Set the gamma and beta parameters from a json object.
double gamma() const
parameter for blending accelerations in the velocity update.
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 dv_dx(const unsigned prev_ti=0) const override
Compute the derivative of the velocity with respect to the solution.
void update_quantities(const Eigen::VectorXd &x) override
Update the time integration quantities (i.e., , , and ).
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.
const Eigen::VectorXd & a_prev() const
Get the most recent previous acceleration value.
void set_x_prev(const Eigen::VectorXd &x_prev)
Convenience functions for setting the most recent previous solution.
nlohmann::json json
Definition Common.hpp:9