PolyFEM
Loading...
Searching...
No Matches
polyfem::time_integrator::ImplicitNewmark Class Reference

Implicit Newmark-beta method. More...

#include <ImplicitNewmark.hpp>

Inheritance diagram for polyfem::time_integrator::ImplicitNewmark:
[legend]
Collaboration diagram for polyfem::time_integrator::ImplicitNewmark:
[legend]

Public Member Functions

 ImplicitNewmark ()
 
void set_parameters (const json &params) override
 Set the gamma and beta 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 da_dx (const unsigned prev_ti=0) const
 
double beta () const
 \(\beta\) parameter for blending accelerations in the solution update.
 
double gamma () const
 \(\gamma\) parameter for blending accelerations in the velocity update.
 
- 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.
 

Protected Attributes

double beta_ = 0.25
 \(\beta\) parameter for blending accelerations in the solution update.
 
double gamma_ = 0.5
 \(\gamma\) parameter for blending accelerations in the velocity update.
 
- 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.
 

Additional Inherited Members

- Static Public Member Functions inherited from polyfem::time_integrator::ImplicitTimeIntegrator
static std::shared_ptr< ImplicitTimeIntegratorconstruct_time_integrator (const json &params)
 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 inherited from polyfem::time_integrator::ImplicitTimeIntegrator
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.
 

Detailed Description

Implicit Newmark-beta method.

\[ x^{t+1} = x^t + \Delta t v^t + \frac{\Delta t^2}{2}((1-2\beta)a^t + 2 \beta a^{t+1})\newline v^{t+1} = v^t + (1-\gamma)\Delta ta^t + \gamma \Delta ta^{t+1} \]

See also
https://en.wikipedia.org/wiki/Newmark-beta_method

Definition at line 13 of file ImplicitNewmark.hpp.

Constructor & Destructor Documentation

◆ ImplicitNewmark()

polyfem::time_integrator::ImplicitNewmark::ImplicitNewmark ( )
inline

Definition at line 16 of file ImplicitNewmark.hpp.

Member Function Documentation

◆ acceleration_scaling()

double polyfem::time_integrator::ImplicitNewmark::acceleration_scaling ( ) const
overridevirtual

Compute the acceleration scaling used to scale forces when integrating a second order ODE.

\[ \beta \Delta t^2 \]

Implements polyfem::time_integrator::ImplicitTimeIntegrator.

Definition at line 35 of file ImplicitNewmark.cpp.

References beta(), and polyfem::time_integrator::ImplicitTimeIntegrator::dt().

Here is the call graph for this function:

◆ beta()

double polyfem::time_integrator::ImplicitNewmark::beta ( ) const
inline

\(\beta\) parameter for blending accelerations in the solution update.

Definition at line 73 of file ImplicitNewmark.hpp.

References beta_.

Referenced by acceleration_scaling(), compute_velocity(), dv_dx(), and x_tilde().

Here is the caller graph for this function:

◆ compute_acceleration()

Eigen::VectorXd polyfem::time_integrator::ImplicitNewmark::compute_acceleration ( const Eigen::VectorXd &  v) const
overridevirtual

Compute the current acceleration given the current velocity and using the stored previous velocity(s).

\[ a = \frac{v - v^t - (1 - \gamma) \Delta t a^t)}{\gamma \Delta t} \]

Parameters
vcurrent velocity
Returns
value for \(a\)

Implements polyfem::time_integrator::ImplicitTimeIntegrator.

Definition at line 30 of file ImplicitNewmark.cpp.

References polyfem::time_integrator::ImplicitTimeIntegrator::a_prev(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), gamma(), and polyfem::time_integrator::ImplicitTimeIntegrator::v_prev().

Referenced by update_quantities().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_velocity()

Eigen::VectorXd polyfem::time_integrator::ImplicitNewmark::compute_velocity ( const Eigen::VectorXd &  x) const
overridevirtual

Compute the current velocity given the current solution and using the stored previous solution(s).

\[ a^{t+1} = \frac{x - (x^t + \Delta t v^t + \Delta t^2 (0.5 - \beta) a^t)}{\beta \Delta t^2}\newline v = v^t + \Delta t (1 - \gamma) a^t + \Delta t \gamma a^{t+1} = (\gamma / \beta) / \Delta t (x - x^t) + (1 - \gamma / \beta) v^t + (1 - \gamma / \beta / 2) \Delta t a^t \]

Parameters
xcurrent solution vector
Returns
value for \(v\)

Implements polyfem::time_integrator::ImplicitTimeIntegrator.

Definition at line 24 of file ImplicitNewmark.cpp.

References polyfem::time_integrator::ImplicitTimeIntegrator::a_prev(), beta(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), gamma(), polyfem::time_integrator::ImplicitTimeIntegrator::v_prev(), x, and polyfem::time_integrator::ImplicitTimeIntegrator::x_prev().

Referenced by update_quantities().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ da_dx()

double polyfem::time_integrator::ImplicitNewmark::da_dx ( const unsigned  prev_ti = 0) const

Definition at line 52 of file ImplicitNewmark.cpp.

References da_dx(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), dv_dx(), and gamma().

Referenced by da_dx(), and dv_dx().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ dv_dx()

double polyfem::time_integrator::ImplicitNewmark::dv_dx ( const unsigned  prev_ti = 0) const
overridevirtual

Compute the derivative of the velocity with respect to the solution.

\[ \frac{\partial v}{\partial x} = \frac{\gamma}{\beta\Delta t} \]

\[ \frac{\partial v}{\partial x^t} = ? \]

\[ \frac{\partial v}{\partial x^{t-1}} = ? \]

Parameters
prev_tiindex of the previous solution to use (0 -> current; 1 -> previous; 2 -> second previous; etc.)

Implements polyfem::time_integrator::ImplicitTimeIntegrator.

Definition at line 40 of file ImplicitNewmark.cpp.

References beta(), da_dx(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), dv_dx(), and gamma().

Referenced by da_dx(), and dv_dx().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ gamma()

double polyfem::time_integrator::ImplicitNewmark::gamma ( ) const
inline

\(\gamma\) parameter for blending accelerations in the velocity update.

Definition at line 75 of file ImplicitNewmark.hpp.

References gamma_.

Referenced by compute_acceleration(), compute_velocity(), da_dx(), and dv_dx().

Here is the caller graph for this function:

◆ set_parameters()

void polyfem::time_integrator::ImplicitNewmark::set_parameters ( const json params)
overridevirtual

Set the gamma and beta parameters from a json object.

Parameters
paramsjson containing {"gamma": 0.5, "beta": 0.25}

Reimplemented from polyfem::time_integrator::ImplicitTimeIntegrator.

Definition at line 5 of file ImplicitNewmark.cpp.

References beta_, and gamma_.

◆ update_quantities()

void polyfem::time_integrator::ImplicitNewmark::update_quantities ( const Eigen::VectorXd &  x)
overridevirtual

Update the time integration quantities (i.e., \(x\), \(v\), and \(a\)).

Parameters
xnew solution vector

Implements polyfem::time_integrator::ImplicitTimeIntegrator.

Definition at line 11 of file ImplicitNewmark.cpp.

References compute_acceleration(), compute_velocity(), polyfem::time_integrator::ImplicitTimeIntegrator::set_a_prev(), polyfem::time_integrator::ImplicitTimeIntegrator::set_v_prev(), polyfem::time_integrator::ImplicitTimeIntegrator::set_x_prev(), and x.

Here is the call graph for this function:

◆ x_tilde()

Eigen::VectorXd polyfem::time_integrator::ImplicitNewmark::x_tilde ( ) const
overridevirtual

Compute the predicted solution to be used in the inertia term \((x-\tilde{x})^TM(x-\tilde{x})\).

\[ \tilde{x} = x^t + \Delta t (v^t + (0.5 - \beta) \Delta t a^t) \]

Returns
value for \(\tilde{x}\)

Implements polyfem::time_integrator::ImplicitTimeIntegrator.

Definition at line 19 of file ImplicitNewmark.cpp.

References polyfem::time_integrator::ImplicitTimeIntegrator::a_prev(), beta(), polyfem::time_integrator::ImplicitTimeIntegrator::dt(), polyfem::time_integrator::ImplicitTimeIntegrator::v_prev(), and polyfem::time_integrator::ImplicitTimeIntegrator::x_prev().

Here is the call graph for this function:

Member Data Documentation

◆ beta_

double polyfem::time_integrator::ImplicitNewmark::beta_ = 0.25
protected

\(\beta\) parameter for blending accelerations in the solution update.

Definition at line 79 of file ImplicitNewmark.hpp.

Referenced by beta(), and set_parameters().

◆ gamma_

double polyfem::time_integrator::ImplicitNewmark::gamma_ = 0.5
protected

\(\gamma\) parameter for blending accelerations in the velocity update.

Definition at line 81 of file ImplicitNewmark.hpp.

Referenced by gamma(), and set_parameters().


The documentation for this class was generated from the following files: