PolyFEM
Loading...
Searching...
No Matches
RayleighDampingForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
8
9namespace polyfem::solver
10{
13 {
14 public:
17 const Form &form_to_damp,
18 const time_integrator::ImplicitTimeIntegrator &time_integrator,
19 const bool use_stiffness_as_ratio,
20 const double stiffness,
21 const int n_lagging_iters);
22
23 static std::shared_ptr<RayleighDampingForm> create(
24 const json &params,
25 const std::unordered_map<std::string, std::shared_ptr<Form>> &forms,
26 const time_integrator::ImplicitTimeIntegrator &time_integrator);
27
28 std::string name() const override { return "rayleigh-damping"; }
29
30 protected:
34 double value_unweighted(const Eigen::VectorXd &x) const override;
35
39 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
40
44 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
45
46 public:
49 void init_lagging(const Eigen::VectorXd &x) override;
50
54 void update_lagging(const Eigen::VectorXd &x, const int iter_num) override;
55
57 int max_lagging_iterations() const override { return n_lagging_iters_; }
58
61 bool uses_lagging() const override { return true; }
62
64 double stiffness() const;
65
66 private:
70 const double stiffness_;
71 const int n_lagging_iters_;
72
74 };
75} // namespace polyfem::solver
int x
Tikonov regularization form between x and x_lagged.
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
static std::shared_ptr< RayleighDampingForm > create(const json &params, const std::unordered_map< std::string, std::shared_ptr< Form > > &forms, const time_integrator::ImplicitTimeIntegrator &time_integrator)
void update_lagging(const Eigen::VectorXd &x, const int iter_num) override
Update lagged fields.
bool uses_lagging() const override
Does this form require lagging?
void init_lagging(const Eigen::VectorXd &x) override
Initialize lagged fields.
const int n_lagging_iters_
Number of iterations to lag for.
const bool use_stiffness_as_ratio_
Whether to use the stiffness ratio or the stiffness value.
int max_lagging_iterations() const override
Get the maximum number of lagging iteration allowable.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
const double stiffness_
Damping stiffness coefficient.
const time_integrator::ImplicitTimeIntegrator & time_integrator_
Reference to the time integrator.
StiffnessMatrix lagged_stiffness_matrix_
The lagged stiffness matrix.
const Form & form_to_damp_
Reference to the form we are damping.
double stiffness() const
Get the stiffness of the form.
Implicit time integrator of a second order ODE (equivently a system of coupled first order ODEs).
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22