PolyFEM
Loading...
Searching...
No Matches
LinearElasticity.hpp
Go to the documentation of this file.
1#pragma once
2
7
8// local assembler for linear elasticity
9namespace polyfem::assembler
10{
12 {
13 public:
18
20 // vals stores the evaluation for that element
21 // da contains both the quadrature weight and the change of metric in the integral
22 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>
23 assemble(const LinearAssemblerData &data) const override;
24
25 // compute elastic energy
26 double compute_energy(const NonLinearAssemblerData &data) const override;
27 // neccessary for mixing linear model with non-linear collision response
28 Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override;
29 // compute gradient of elastic energy, as assembler
30 Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override;
31
32 // kernel of the pde, used in kernel problem
33 Eigen::Matrix<AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1> kernel(const int dim, const AutodiffGradPt &r, const AutodiffScalarGrad &) const override;
34
35 // uses autodiff to compute the rhs for a fabbricated solution
36 // uses autogenerated code to compute div(sigma)
37 // pt is the evaluation of the solution at a point
38 VectorNd compute_rhs(const AutodiffHessianPt &pt) const override;
39
40 void compute_stiffness_value(const double t,
42 const Eigen::MatrixXd &local_pts,
43 const Eigen::MatrixXd &displacement,
44 Eigen::MatrixXd &tensor) const override;
45
47 const Eigen::MatrixXd &mat,
48 Eigen::MatrixXd &stress,
49 Eigen::MatrixXd &result) const override;
50
52 Eigen::MatrixXd &stress,
53 Eigen::MatrixXd &result) const override;
54
56 Eigen::MatrixXd &dstress_dmu,
57 Eigen::MatrixXd &dstress_dlambda) const override;
58
59 // inialize material parameter
60 void add_multimaterial(const int index, const json &params, const Units &units) override;
61
62 // class that stores and compute lame parameters per point
63 const LameParameters &lame_params() const { return params_; }
64 void set_params(const LameParameters &params) { params_ = params; }
65
66 void update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus) override
67 {
68 params_.lambda_mat_ = lambdas;
69 params_.mu_mat_ = mus;
70 }
71
72 virtual bool is_linear() const override { return true; }
73
74 std::string name() const override { return "LinearElasticity"; }
75 bool allow_inversion() const override { return true; }
76 std::map<std::string, ParamFunc> parameters() const override;
77
78 void assign_stress_tensor(const OutputData &data,
79 const int all_size,
80 const ElasticityTensorType &type,
81 Eigen::MatrixXd &all,
82 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override;
83
84 private:
85 // class that stores and compute lame parameters per point
87
88 // aux function that computes energy
89 // double compute_energy is the same with T=double
90 // assemble_gradient is the same with T=DScalar1 and return .getGradient()
91 template <typename T>
92 T compute_energy_aux(const NonLinearAssemblerData &data) const;
93 };
94} // namespace polyfem::assembler
ElementAssemblyValues vals
Definition Assembler.cpp:22
stores per element basis values at given quadrature points and geometric mapping
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 ...
void compute_stress_grad_multiply_stress(const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
Eigen::Matrix< AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1 > kernel(const int dim, const AutodiffGradPt &r, const AutodiffScalarGrad &) const override
std::map< std::string, ParamFunc > parameters() const override
void update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus) override
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
void add_multimaterial(const int index, const json &params, const Units &units) override
void compute_dstress_dmu_dlambda(const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) const override
virtual bool is_linear() const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
void set_params(const LameParameters &params)
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
const LameParameters & lame_params() const
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
computes local stiffness matrix is R^{dimĀ²} for bases i,j
std::string name() const override
double compute_energy(const NonLinearAssemblerData &data) const override
T compute_energy_aux(const NonLinearAssemblerData &data) const
void compute_stress_grad_multiply_mat(const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_stiffness_value(const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const override
double assemble_energy(const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const override
void assemble_gradient(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, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, Eigen::MatrixXd &rhs) const override
void assemble_hessian(const bool is_volume, const int n_basis, const bool project_to_psd, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, utils::MatrixCache &mat_cache, StiffnessMatrix &grad) const override
Used for test only.
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffGradPt
Automatic differentiation scalar with first-order derivatives.
Definition autodiff.h:112