PolyFEM
Loading...
Searching...
No Matches
IncompressibleLinElast.hpp
Go to the documentation of this file.
1#pragma once
2
7
8// local assembler for incompressible model, pressure is separate (see Stokes)
9namespace polyfem::assembler
10{
11 // displacement assembler
13 {
14 public:
16 // res is R^{dimĀ²}
17 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>
18 assemble(const LinearAssemblerData &data) const override;
19
20 void add_multimaterial(const int index, const json &params, const Units &units) override;
21 void set_params(const LameParameters &params) { params_ = params; }
22
23 std::string name() const override { return "IncompressibleLinearElasticity"; }
24 bool allow_inversion() const override { return true; }
25 std::map<std::string, ParamFunc> parameters() const override;
26
27 protected:
28 void assign_stress_tensor(const OutputData &data,
29 const int all_size,
30 const ElasticityTensorType &type,
31 Eigen::MatrixXd &all,
32 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override;
33
34 private:
36 };
37
38 // mixed, displacement and pressure
40 {
41 public:
42 std::string name() const override { return "IncompressibleLinearElasticityMixed"; }
43
44 // res is R^{dim}
45 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1>
46 assemble(const MixedAssemblerData &data) const override;
47
48 inline int rows() const override { return size_; }
49 inline int cols() const override { return 1; }
50 };
51
52 // pressure only part
54 {
55 public:
57
58 // res is R^{1}
59 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>
60 assemble(const LinearAssemblerData &data) const override;
61
62 void add_multimaterial(const int index, const json &params, const Units &units) override;
63 void set_params(const LameParameters &params) { params_ = params; }
64
65 std::string name() const override { return "IncompressibleLinearElasticityPressure"; }
66 std::map<std::string, ParamFunc> parameters() const override { return std::map<std::string, ParamFunc>(); }
67
68 void set_size(const int size) override
69 {
71 size_ = 1;
72 }
73
74 private:
76 int disp_size_ = 0;
77 };
78} // namespace polyfem::assembler
void add_multimaterial(const int index, const json &params, const Units &units) override
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
local assembly function that defines the bilinear form (LHS) computes and returns a single local stif...
std::map< std::string, ParamFunc > parameters() const override
void assign_stress_tensor(const OutputData &data, const int all_size, const ElasticityTensorType &type, Eigen::MatrixXd &all, const std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > assemble(const MixedAssemblerData &data) const override
std::map< std::string, ParamFunc > parameters() const override
void add_multimaterial(const int index, const json &params, const Units &units) override
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
local assembly function that defines the bilinear form (LHS) computes and returns a single local stif...
assemble matrix based on the local assembler local assembler is eg Laplace, LinearElasticity etc
void assemble(const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, StiffnessMatrix &stiffness, const bool is_mass=false) const override
assembles the stiffness matrix for the given basis the bilinear form (local assembler) is encoded by ...
Used for test only.
nlohmann::json json
Definition Common.hpp:9