PolyFEM
Loading...
Searching...
No Matches
FixedCorotational.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
26 const Eigen::MatrixXd &mat,
27 Eigen::MatrixXd &stress,
28 Eigen::MatrixXd &result) const override;
29
31 Eigen::MatrixXd &stress,
32 Eigen::MatrixXd &result) const override;
33
35 Eigen::MatrixXd &dstress_dmu,
36 Eigen::MatrixXd &dstress_dlambda) const override;
37
38 // sets material params
39 void add_multimaterial(const int index, const json &params, const Units &units) override;
40
41 void set_params(const LameParameters &params) { params_ = params; }
42 const LameParameters &lame_params() const { return params_; }
43
44 void update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus) override
45 {
46 params_.lambda_mat_ = lambdas;
47 params_.mu_mat_ = mus;
48 }
49
50 std::string name() const override { return "FixedCorotational"; }
51 std::map<std::string, ParamFunc> parameters() const override;
52
53 void assign_stress_tensor(const OutputData &data,
54 const int all_size,
55 const ElasticityTensorType &type,
56 Eigen::MatrixXd &all,
57 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override;
58
59 private:
61
62 // utility function that computes energy, the template is used for double, DScalar1, and DScalar2 in energy, gradient and hessian
63 template <int dim>
64 double compute_energy_aux(const NonLinearAssemblerData &data) const;
65 template <int n_basis, int dim>
66 void compute_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const;
67 template <int n_basis, int dim>
68 void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const;
69
70 template <int dim>
71 static double compute_energy_from_singular_values(const Eigen::Vector<double, dim> &sigmas, const double lambda, const double mu);
72 template <int dim>
73 static Eigen::Vector<double, dim> compute_stress_from_singular_values(const Eigen::Vector<double, dim> &sigmas, const double lambda, const double mu);
74 template <int dim>
75 static Eigen::Matrix<double, dim, dim> compute_stiffness_from_singular_values(const Eigen::Vector<double, dim> &sigmas, const double lambda, const double mu);
76
77 template <int dim>
78 static double compute_energy_from_def_grad(const Eigen::Matrix<double, dim, dim> &F, const double lambda, const double mu);
79 template <int dim>
80 static Eigen::Matrix<double, dim, dim> compute_stress_from_def_grad(const Eigen::Matrix<double, dim, dim> &F, const double lambda, const double mu);
81 template <int dim>
82 static Eigen::Matrix<double, dim*dim, dim*dim> compute_stiffness_from_def_grad(const Eigen::Matrix<double, dim, dim> &F, const double lambda, const double mu);
83 };
84} // namespace polyfem::assembler
const LameParameters & lame_params() const
static double compute_energy_from_def_grad(const Eigen::Matrix< double, dim, dim > &F, const double lambda, const double mu)
static Eigen::Matrix< double, dim, dim > compute_stress_from_def_grad(const Eigen::Matrix< double, dim, dim > &F, const double lambda, const double mu)
void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
static Eigen::Vector< double, dim > compute_stress_from_singular_values(const Eigen::Vector< double, dim > &sigmas, const double lambda, const double mu)
static double compute_energy_from_singular_values(const Eigen::Vector< double, dim > &sigmas, const double lambda, const double mu)
void add_multimaterial(const int index, const json &params, const Units &units) override
static Eigen::Matrix< double, dim *dim, dim *dim > compute_stiffness_from_def_grad(const Eigen::Matrix< double, dim, dim > &F, const double lambda, const double mu)
double compute_energy_aux(const NonLinearAssemblerData &data) const
double compute_energy(const NonLinearAssemblerData &data) const override
static Eigen::Matrix< double, dim, dim > compute_stiffness_from_singular_values(const Eigen::Vector< double, dim > &sigmas, const double lambda, const double mu)
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
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 update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus) override
void compute_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
void compute_stress_grad_multiply_stress(const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_dstress_dmu_dlambda(const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) 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.
nlohmann::json json
Definition Common.hpp:9