PolyFEM
Loading...
Searching...
No Matches
AMIPSEnergy.cpp
Go to the documentation of this file.
1#include "AMIPSEnergy.hpp"
2
4
9
10namespace polyfem::assembler
11{
12 void AMIPSEnergy::add_multimaterial(const int index, const json &params, const Units &units, const std::string &root_path)
13 {
14 assert(size() == 2 || size() == 3);
15
16 if (params.contains("use_rest_pose"))
17 {
18 use_rest_pose_ = params["use_rest_pose"].get<bool>();
19 }
20
21 if (energy_weights_.size() <= index)
22 energy_weights_.resize(index + 1, 1.0);
23
24 if (params.contains("weight"))
25 energy_weights_[index] = params["weight"].get<double>();
26 else
27 energy_weights_[index] = 1.0;
28 }
29
30 double AMIPSEnergy::get_energy_weight(const int el_id) const
31 {
32 if (energy_weights_.empty())
33 return 1.0;
34 if (energy_weights_.size() == 1)
35 return energy_weights_[0];
36 if (el_id >= 0 && el_id < (int)energy_weights_.size())
37 return energy_weights_[el_id];
38 return 1.0;
39 }
40
41 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> AMIPSEnergy::gradient(
42 const RowVectorNd &p,
43 const double t,
44 const int el_id,
45 const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> &F) const
46 {
47 const double det = polyfem::utils::determinant(F);
48 if (det <= 0)
49 {
50 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> grad(size(), size());
51 grad.setConstant(std::nan(""));
52 return grad;
53 }
54
55 const double weight = get_energy_weight(el_id);
56
58 {
59 if (size() == 2)
60 return weight * autogen::AMIPS2drest_gradient(p, t, el_id, F);
61 else
62 return weight * autogen::AMIPS3drest_gradient(p, t, el_id, F);
63 }
64 else
65 {
66 if (size() == 2)
67 return weight * autogen::AMIPS2d_gradient(p, t, el_id, F);
68 else
69 return weight * autogen::AMIPS3d_gradient(p, t, el_id, F);
70 }
71 }
72
73 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9> AMIPSEnergy::hessian(
74 const RowVectorNd &p,
75 const double t,
76 const int el_id,
77 const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> &F) const
78 {
79 const double det = polyfem::utils::determinant(F);
80 if (det <= 0)
81 {
82 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9> hessian(size() * size(), size() * size());
83 hessian.setConstant(std::nan(""));
84 return hessian;
85 }
86
87 const double weight = get_energy_weight(el_id);
88
90 {
91 if (size() == 2)
92 return weight * autogen::AMIPS2drest_hessian(p, t, el_id, F);
93 else
94 return weight * autogen::AMIPS3drest_hessian(p, t, el_id, F);
95 }
96 else
97 {
98 if (size() == 2)
99 return weight * autogen::AMIPS2d_hessian(p, t, el_id, F);
100 else
101 return weight * autogen::AMIPS3d_hessian(p, t, el_id, F);
102 }
103 }
104
105} // namespace polyfem::assembler
double get_energy_weight(const int el_id) const
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > gradient(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F) const override
std::vector< double > energy_weights_
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9 > hessian(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F) const override
void add_multimaterial(const int index, const json &params, const Units &units, const std::string &root_path) override
Used for test only.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > AMIPS2drest_gradient(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > AMIPS3drest_gradient(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > AMIPS2d_gradient(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Definition AMIPS2d.hpp:9
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9 > AMIPS3d_hessian(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Definition AMIPS3d.hpp:58
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9 > AMIPS3drest_hessian(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9 > AMIPS2drest_hessian(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 9, 9 > AMIPS2d_hessian(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Definition AMIPS2d.hpp:34
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > AMIPS3d_gradient(const RowVectorNd &p, const double t, const int el_id, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > &F)
Definition AMIPS3d.hpp:9
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