PolyFEM
Loading...
Searching...
No Matches
ElasticForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
8
10
11namespace polyfem::solver
12{
14 class ElasticForm : public Form
15 {
16 public:
19 ElasticForm(const int n_bases,
20 const std::vector<basis::ElementBases> &bases,
21 const std::vector<basis::ElementBases> &geom_bases,
22 const assembler::Assembler &assembler,
23 const assembler::AssemblyValsCache &ass_vals_cache,
24 const double t, const double dt,
25 const bool is_volume);
26
27 std::string name() const override { return "elastic"; }
28
29 protected:
33 virtual double value_unweighted(const Eigen::VectorXd &x) const override;
34
38 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
39
43 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
44
48 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
49
50 public:
55 bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
56
60 void update_quantities(const double t, const Eigen::VectorXd &x) override
61 {
62 t_ = t;
63 x_prev_ = x;
64 }
65
71 void force_material_derivative(const double t, const Eigen::MatrixXd &x, const Eigen::MatrixXd &x_prev, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term);
72
79 void force_shape_derivative(const double t, const int n_verts, const Eigen::MatrixXd &x, const Eigen::MatrixXd &x_prev, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term);
80
81 private:
82 const int n_bases_;
83 const std::vector<basis::ElementBases> &bases_;
84 const std::vector<basis::ElementBases> &geom_bases_;
85
88 double t_;
89 const double dt_;
90 const bool is_volume_;
91
93 mutable std::unique_ptr<utils::MatrixCache> mat_cache_;
94
97
98 Eigen::VectorXd x_prev_;
99 };
100} // namespace polyfem::solver
int x
Caches basis evaluation and geometric mapping at every element.
Form of the elasticity potential and forces.
const assembler::Assembler & assembler_
Reference to the assembler.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
void compute_cached_stiffness()
Compute the stiffness matrix (cached)
std::unique_ptr< utils::MatrixCache > mat_cache_
Matrix cache (mutable because it is modified in second_derivative_unweighted)
void force_material_derivative(const double t, const Eigen::MatrixXd &x, const Eigen::MatrixXd &x_prev, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term)
Compute the derivative of the force wrt lame/damping parameters, then multiply the resulting matrix w...
void force_shape_derivative(const double t, const int n_verts, const Eigen::MatrixXd &x, const Eigen::MatrixXd &x_prev, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term)
Compute the derivative of the force wrt vertex positions, then multiply the resulting matrix with adj...
virtual double value_unweighted(const Eigen::VectorXd &x) const override
Compute the elastic potential value.
const std::vector< basis::ElementBases > & geom_bases_
StiffnessMatrix cached_stiffness_
Cached stiffness matrix for linear elasticity.
virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
const assembler::AssemblyValsCache & ass_vals_cache_
Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied with the weigth.
const std::vector< basis::ElementBases > & bases_
bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine if a step from solution x0 to solution x1 is allowed.
std::string name() const override
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time-dependent fields.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22