PolyFEM
Loading...
Searching...
No Matches
ElasticForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
8
11
12namespace polyfem::solver
13{
18 {ElementInversionCheck::Conservative, "Conservative"}})
19
21 class ElasticForm : public Form
22 {
23 public:
26 ElasticForm(const int n_bases,
27 std::vector<basis::ElementBases> &bases,
28 const std::vector<basis::ElementBases> &geom_bases,
29 const assembler::Assembler &assembler,
30 assembler::AssemblyValsCache &ass_vals_cache,
31 const double t, const double dt,
32 const bool is_volume,
33 const double jacobian_threshold = 0.,
35
36 std::string name() const override { return "elastic"; }
37
38 protected:
42 virtual double value_unweighted(const Eigen::VectorXd &x) const override;
43
47 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
48
52 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
53
57 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
58
59 public:
64 bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
65
68 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
69
73 void update_quantities(const double t, const Eigen::VectorXd &x) override
74 {
75 t_ = t;
76 x_prev_ = x;
77 }
78
83 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
84
87 void solution_changed(const Eigen::VectorXd &new_x) override;
88
94 void force_material_derivative(const double t, const Eigen::MatrixXd &x, const Eigen::MatrixXd &x_prev, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term);
95
102 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);
103
105 void finish() override;
106
107 private:
108 const int n_bases_;
109 std::vector<basis::ElementBases> &bases_;
110 const std::vector<basis::ElementBases> &geom_bases_;
111
112 const assembler::Assembler &assembler_;
113 assembler::AssemblyValsCache &ass_vals_cache_;
114 double t_;
115 const double jacobian_threshold_;
116 const ElementInversionCheck check_inversion_;
117 const double dt_;
118 const bool is_volume_;
119
120 StiffnessMatrix cached_stiffness_;
121 mutable std::unique_ptr<utils::MatrixCache> mat_cache_;
122
124 void compute_cached_stiffness();
125
126 Eigen::VectorXd x_prev_;
127
128 mutable std::vector<utils::Tree> quadrature_hierarchy_;
129 int quadrature_order_;
130 };
131} // namespace polyfem::solver
int x
NLOHMANN_JSON_SERIALIZE_ENUM(CollisionProxyTessellation, {{CollisionProxyTessellation::REGULAR, "regular"}, {CollisionProxyTessellation::IRREGULAR, "irregular"}})
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22