PolyFEM
Loading...
Searching...
No Matches
MultiModel.hpp
Go to the documentation of this file.
1#pragma once
2
12
13namespace polyfem::assembler
14{
16 {
17 public:
21
22 // compute elastic energy
23 double compute_energy(const NonLinearAssemblerData &data) const override;
24 // neccessary for mixing linear model with non-linear collision response
25 Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override;
26 // compute gradient of elastic energy, as assembler
27 Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override;
28
29 // uses autodiff to compute the rhs for a fabbricated solution
30 // uses autogenerated code to compute div(sigma)
31 // pt is the evaluation of the solution at a point
32 VectorNd compute_rhs(const AutodiffHessianPt &pt) const override;
33 void set_size(const int size) override;
34
35 // inialize material parameter
36 void add_multimaterial(const int index, const json &params, const Units &units) override;
37
38 // initialized multi models
39 inline void init_multimodels(const std::vector<std::string> &mats) { multi_material_models_ = mats; }
40
41 std::string name() const override { return "MultiModels"; }
42 std::map<std::string, ParamFunc> parameters() const override;
43
44 protected:
45 void assign_stress_tensor(const OutputData &data,
46 const int all_size,
47 const ElasticityTensorType &type,
48 Eigen::MatrixXd &all,
49 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override;
50
51 private:
52 std::vector<std::string> multi_material_models_;
53
63 };
64} // namespace polyfem::assembler
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
std::vector< std::string > multi_material_models_
UnconstrainedOgdenElasticity unconstrained_ogden_elasticity_
LinearElasticity linear_elasticity_
FixedCorotational fixed_corotational_
MooneyRivlinElasticity mooney_rivlin_elasticity_
void set_size(const int size) override
void add_multimaterial(const int index, const json &params, const Units &units) override
SaintVenantElasticity saint_venant_
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
double compute_energy(const NonLinearAssemblerData &data) const override
MooneyRivlin3ParamElasticity mooney_rivlin_3_param_elasticity_
std::map< std::string, ParamFunc > parameters() const override
IncompressibleOgdenElasticity incompressible_ogden_elasticity_
HookeLinearElasticity hooke_
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
NeoHookeanElasticity neo_hookean_
std::string name() const override
void init_multimodels(const std::vector< std::string > &mats)
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