PolyFEM
Loading...
Searching...
No Matches
ElasticForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
8
12
13#include <memory>
14#include <vector>
15
16namespace polyfem::solver
17{
19 {
22 };
26 {ElementInversionCheck::Conservative, "Conservative"}})
27
29 class ElasticForm : public Form
30 {
31 friend class ElasticForceDerivative;
32
33 public:
36 ElasticForm(const int n_bases,
37 std::vector<basis::ElementBases> &bases,
38 const std::vector<basis::ElementBases> &geom_bases,
39 const assembler::Assembler &assembler,
40 assembler::AssemblyValsCache &ass_vals_cache,
41 const double t, const double dt,
42 const bool is_volume,
43 const double jacobian_threshold = 0.,
45
46 std::string name() const override { return "elastic"; }
47
48 protected:
52 virtual double value_unweighted(const Eigen::VectorXd &x) const override;
53
57 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
58
62 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
63
67 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
68
69 public:
74 bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
75
78 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
79
83 void update_quantities(const double t, const Eigen::VectorXd &x) override
84 {
85 t_ = t;
86 x_prev_ = x;
87 }
88
93 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
94
97 void solution_changed(const Eigen::VectorXd &new_x) override;
98
100 void finish() override;
101
102 private:
103 const int n_bases_;
104 std::vector<basis::ElementBases> &bases_;
105 const std::vector<basis::ElementBases> &geom_bases_;
106
107 const assembler::Assembler &assembler_;
108 assembler::AssemblyValsCache &ass_vals_cache_;
109 double t_;
110 const double jacobian_threshold_;
111 const ElementInversionCheck check_inversion_;
112 const double dt_;
113 const bool is_volume_;
114
115 StiffnessMatrix cached_stiffness_;
116 mutable std::unique_ptr<utils::MatrixCache> mat_cache_;
117
119 void compute_cached_stiffness();
120
121 Eigen::VectorXd x_prev_;
122
123 mutable std::vector<utils::Tree> quadrature_hierarchy_;
124 int quadrature_order_;
125 };
126} // 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:24