PolyFEM
Loading...
Searching...
No Matches
AMIPSEnergy.hpp
Go to the documentation of this file.
1#pragma once
2
7
8#include <polyfem/Common.hpp>
9
10#include <Eigen/Dense>
11
12#include <functional>
13#include <iostream>
14#include <vector>
15
16namespace polyfem::assembler
17{
19 {
20 public:
22
26
27 // energy, gradient, and hessian used in newton method
28 double compute_energy(const NonLinearAssemblerData &data) const override;
29 Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override;
30 Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override;
31
32 // sets material params
33 void add_multimaterial(const int index, const json &params, const Units &units) override {}
34
35 std::string name() const override { return "AMIPS"; }
36 std::map<std::string, ParamFunc> parameters() const override { return std::map<std::string, ParamFunc>(); }
37
38 void assign_stress_tensor(const OutputData &data,
39 const int all_size,
40 const ElasticityTensorType &type,
41 Eigen::MatrixXd &all,
42 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override;
43
44 bool allow_inversion() const override { return true; }
45
46 private:
47 // utility function that computes energy, the template is used for double, DScalar1, and DScalar2 in energy, gradient and hessian
48 template <typename T>
49 T compute_energy_aux(const NonLinearAssemblerData &data) const;
50 template <int n_basis, int dim>
51 void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const;
52 template <int n_basis, int dim>
53 void compute_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const;
54 };
55
56 class AMIPSEnergyAutodiff : public GenericElastic<AMIPSEnergyAutodiff>
57 {
58 public:
60
61 // sets material params
62 void add_multimaterial(const int index, const json &params, const Units &units) override;
63
64 std::string name() const override { return "AMIPSAutodiff"; }
65 std::map<std::string, ParamFunc> parameters() const override;
66
67 template <typename T>
69 const RowVectorNd &p,
70 const double t,
71 const int el_id,
72 const DefGradMatrix<T> &def_grad) const
73 {
74 using std::pow;
75 using MatrixNT = Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3>;
76
77 MatrixNT F = MatrixNT::Zero(size(), size());
78 for (int i = 0; i < size(); ++i)
79 {
80 for (int j = 0; j < size(); ++j)
81 {
82 for (int k = 0; k < size(); ++k)
83 {
84 F(i, j) += def_grad(i, k) * canonical_transformation_[el_id](k, j);
85 }
86 }
87 }
88
90 if (J <= 0)
91 J = T(std::nan(""));
92 return (F.transpose() * F).trace() / pow(J, 2. / size());
93 }
94
95 private:
96 std::vector<Eigen::MatrixXd> canonical_transformation_;
97 };
98} // namespace polyfem::assembler
std::map< std::string, ParamFunc > parameters() const override
std::vector< Eigen::MatrixXd > canonical_transformation_
void add_multimaterial(const int index, const json &params, const Units &units) override
T elastic_energy(const RowVectorNd &p, const double t, const int el_id, const DefGradMatrix< T > &def_grad) const
std::string name() const override
double compute_energy(const NonLinearAssemblerData &data) const override
std::string name() const override
T compute_energy_aux(const NonLinearAssemblerData &data) const
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_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const
void add_multimaterial(const int index, const json &params, const Units &units) override
bool allow_inversion() const override
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
std::map< std::string, ParamFunc > parameters() 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< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > DefGradMatrix
T determinant(const Eigen::Matrix< T, rows, cols, option, maxRow, maxCol > &mat)
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13