PolyFEM
Loading...
Searching...
No Matches
ImplicitTimeIntegrator.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
4
5#include <Eigen/Core>
6
7#include <map>
8#include <vector>
9#include <deque>
10
12{
15 {
16 public:
18 virtual ~ImplicitTimeIntegrator() = default;
19
22 virtual void set_parameters(const json &params) {}
23
30 virtual void init(const Eigen::MatrixXd &x_prevs, const Eigen::MatrixXd &v_prevs, const Eigen::MatrixXd &a_prevs, double dt);
31
34 virtual void update_quantities(const Eigen::VectorXd &x) = 0;
35
38 virtual Eigen::VectorXd x_tilde() const = 0;
39
43 virtual Eigen::VectorXd compute_velocity(const Eigen::VectorXd &x) const = 0;
44
48 virtual Eigen::VectorXd compute_acceleration(const Eigen::VectorXd &v) const = 0;
49
52 virtual double acceleration_scaling() const = 0;
53
56 virtual double dv_dx(const unsigned prev_ti = 0) const = 0;
57
59 const double &dt() const { return dt_; }
60
63 virtual void save_state(const std::string &state_path) const;
64
68 static std::shared_ptr<ImplicitTimeIntegrator> construct_time_integrator(const json &params);
69
72 static const std::vector<std::string> &get_time_integrator_names();
73
75 const Eigen::VectorXd &x_prev() const { return x_prevs_.front(); }
77 const Eigen::VectorXd &v_prev() const { return v_prevs_.front(); }
79 const Eigen::VectorXd &a_prev() const { return a_prevs_.front(); }
80
82 const std::deque<Eigen::VectorXd> &x_prevs() const { return x_prevs_; }
84 const std::deque<Eigen::VectorXd> &v_prevs() const { return v_prevs_; }
86 const std::deque<Eigen::VectorXd> &a_prevs() const { return a_prevs_; }
87
89 int steps() const { return x_prevs_.size(); }
90
91 protected:
93 virtual int max_steps() const { return 1; }
94
97 double dt_ = 1;
98
100 std::deque<Eigen::VectorXd> x_prevs_;
102 std::deque<Eigen::VectorXd> v_prevs_;
104 std::deque<Eigen::VectorXd> a_prevs_;
105
107 void set_x_prev(const Eigen::VectorXd &x_prev) { x_prevs_.front() = x_prev; }
109 void set_v_prev(const Eigen::VectorXd &v_prev) { v_prevs_.front() = v_prev; }
111 void set_a_prev(const Eigen::VectorXd &a_prev) { a_prevs_.front() = a_prev; }
112 };
113} // namespace polyfem::time_integrator
int x
Implicit time integrator of a second order ODE (equivently a system of coupled first order ODEs).
const std::deque< Eigen::VectorXd > & v_prevs() const
Get the (relevant) history of previous velocity value.
const Eigen::VectorXd & v_prev() const
Get the most recent previous velocity value.
std::deque< Eigen::VectorXd > x_prevs_
Store the necessary previous values of the solution for single or multi-step integration.
static std::shared_ptr< ImplicitTimeIntegrator > construct_time_integrator(const json &params)
Factory method for constructing implicit time integrators from the name of the integrator.
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.
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)...
const std::deque< Eigen::VectorXd > & a_prevs() const
Get the (relevant) history of previous acceleration value.
void set_v_prev(const Eigen::VectorXd &v_prev)
Convenience functions for setting the most recent previous velocity.
virtual int max_steps() const
Get the maximum number of steps to use for integration.
const double & dt() const
Access the time step size.
virtual void update_quantities(const Eigen::VectorXd &x)=0
Update the time integration quantities (i.e., , , and ).
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 void set_parameters(const json &params)
Set the time integrator parameters from a json object.
std::deque< Eigen::VectorXd > a_prevs_
Store the necessary previous values of the acceleration for single or multi-step integration.
const Eigen::VectorXd & a_prev() const
Get the most recent previous acceleration value.
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 , , and .
const std::deque< Eigen::VectorXd > & x_prevs() const
Get the (relevant) history of previous solution value.
int steps() const
Get the current 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.
std::deque< Eigen::VectorXd > v_prevs_
Store the necessary previous values of the velocity for single or multi-step integration.
virtual double dv_dx(const unsigned prev_ti=0) const =0
Compute the derivative of the velocity with respect to the solution.
virtual double acceleration_scaling() const =0
Compute the acceleration scaling used to scale forces when integrating a second order ODE.
static const std::vector< std::string > & get_time_integrator_names()
Get a vector of the names of possible ImplicitTimeIntegrators.
virtual void save_state(const std::string &state_path) const
Save the values of , , and .
virtual Eigen::VectorXd x_tilde() const =0
Compute the predicted solution to be used in the inertia term .
nlohmann::json json
Definition Common.hpp:9