PolyFEM
Loading...
Searching...
No Matches
NeoHookeanElasticity.hpp
Go to the documentation of this file.
1#pragma once
2
7
8// non linear NeoHookean material model
9namespace polyfem::assembler
10{
12 {
13 public:
15
19
20 // energy, gradient, and hessian used in newton method
21 double compute_energy(const NonLinearAssemblerData &data) const override;
22 Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override;
23 Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override;
24
25 // rhs for fabbricated solution, compute with automatic sympy code
26 VectorNd compute_rhs(const AutodiffHessianPt &pt) const override;
27
28 void compute_stiffness_value(const double t,
30 const Eigen::MatrixXd &local_pts,
31 const Eigen::MatrixXd &displacement,
32 Eigen::MatrixXd &tensor) const override;
33
35 const Eigen::MatrixXd &mat,
36 Eigen::MatrixXd &stress,
37 Eigen::MatrixXd &result) const override;
38
40 Eigen::MatrixXd &stress,
41 Eigen::MatrixXd &result) const override;
42
44 const Eigen::MatrixXd &vect,
45 Eigen::MatrixXd &stress,
46 Eigen::MatrixXd &result) const override;
47
49 Eigen::MatrixXd &dstress_dmu,
50 Eigen::MatrixXd &dstress_dlambda) const override;
51
52 // sets material params
53 void add_multimaterial(const int index, const json &params, const Units &units) override;
54
55 void set_params(const LameParameters &params) { params_ = params; }
56 const LameParameters &lame_params() const { return params_; }
57
58 void update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus) override
59 {
60 params_.lambda_mat_ = lambdas;
61 params_.mu_mat_ = mus;
62 }
63
64 std::string name() const override { return "NeoHookean"; }
65 bool allow_inversion() const override { return false; }
66 std::map<std::string, ParamFunc> parameters() const override;
67
68 void assign_stress_tensor(const OutputData &data,
69 const int all_size,
70 const ElasticityTensorType &type,
71 Eigen::MatrixXd &all,
72 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override;
73
74 private:
76
77 // utility function that computes energy, the template is used for double, DScalar1, and DScalar2 in energy, gradient and hessian
78 template <typename T>
79 T compute_energy_aux(const NonLinearAssemblerData &data) const;
80 template <int n_basis, int dim>
81 void compute_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const;
82 template <int n_basis, int dim>
83 void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const;
84 };
85} // namespace polyfem::assembler
ElementAssemblyValues vals
Definition Assembler.cpp:22
stores per element basis values at given quadrature points and geometric mapping
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
double compute_energy(const NonLinearAssemblerData &data) const override
std::map< std::string, ParamFunc > parameters() const override
void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
void compute_stress_grad_multiply_stress(const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const
void compute_dstress_dmu_dlambda(const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) const override
void compute_stress_grad_multiply_vect(const OptAssemblerData &data, const Eigen::MatrixXd &vect, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_stress_grad_multiply_mat(const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) 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
void update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus) override
void add_multimaterial(const int index, const json &params, const Units &units) override
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
T compute_energy_aux(const NonLinearAssemblerData &data) const
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
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