PolyFEM
Loading...
Searching...
No Matches
AMIPSForm.cpp
Go to the documentation of this file.
1#include "AMIPSForm.hpp"
2
3#include <polyfem/State.hpp>
4
7
14
15namespace polyfem::solver
16{
17 namespace
18 {
19 void scaled_jacobian(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, Eigen::VectorXd &quality)
20 {
21 const int dim = F.cols() - 1;
22
23 quality.setZero(F.rows());
24 if (dim == 2)
25 {
26 for (int i = 0; i < F.rows(); i++)
27 {
28 Eigen::RowVector3d e0;
29 e0(2) = 0;
30 e0.head(2) = V.row(F(i, 2)) - V.row(F(i, 1));
31 Eigen::RowVector3d e1;
32 e1(2) = 0;
33 e1.head(2) = V.row(F(i, 0)) - V.row(F(i, 2));
34 Eigen::RowVector3d e2;
35 e2(2) = 0;
36 e2.head(2) = V.row(F(i, 1)) - V.row(F(i, 0));
37
38 double l0 = e0.norm();
39 double l1 = e1.norm();
40 double l2 = e2.norm();
41
42 double A = 0.5 * (e0.cross(e1)).norm();
43 double Lmax = std::max(l0 * l1, std::max(l1 * l2, l0 * l2));
44
45 quality(i) = 2 * A * (2 / sqrt(3)) / Lmax;
46 }
47 }
48 else
49 {
50 for (int i = 0; i < F.rows(); i++)
51 {
52 Eigen::RowVector3d e0 = V.row(F(i, 1)) - V.row(F(i, 0));
53 Eigen::RowVector3d e1 = V.row(F(i, 2)) - V.row(F(i, 1));
54 Eigen::RowVector3d e2 = V.row(F(i, 0)) - V.row(F(i, 2));
55 Eigen::RowVector3d e3 = V.row(F(i, 3)) - V.row(F(i, 0));
56 Eigen::RowVector3d e4 = V.row(F(i, 3)) - V.row(F(i, 1));
57 Eigen::RowVector3d e5 = V.row(F(i, 3)) - V.row(F(i, 2));
58
59 double l0 = e0.norm();
60 double l1 = e1.norm();
61 double l2 = e2.norm();
62 double l3 = e3.norm();
63 double l4 = e4.norm();
64 double l5 = e5.norm();
65
66 double J = std::abs((e0.cross(e3)).dot(e2));
67
68 double a1 = l0 * l2 * l3;
69 double a2 = l0 * l1 * l4;
70 double a3 = l1 * l2 * l5;
71 double a4 = l3 * l4 * l5;
72
73 double a = std::max({a1, a2, a3, a4, J});
74 quality(i) = J * sqrt(2) / a;
75 }
76 }
77 }
78 } // namespace
79
80 double MinJacobianForm::value_unweighted(const Eigen::VectorXd &x) const
81 {
82 const bool is_volume = state_.mesh->is_volume();
83 double min_jacs = std::numeric_limits<double>::max();
84 for (size_t e = 0; e < state_.geom_bases().size(); ++e)
85 {
86 if (state_.mesh->is_polytope(e))
87 continue;
88
89 const auto &gbasis = state_.geom_bases()[e];
90 const int n_local_bases = int(gbasis.bases.size());
91
93 gbasis.compute_quadrature(quad);
94
95 std::vector<assembler::AssemblyValues> tmp;
96
97 Eigen::MatrixXd dx = Eigen::MatrixXd::Zero(quad.points.rows(), quad.points.cols());
98 Eigen::MatrixXd dy = Eigen::MatrixXd::Zero(quad.points.rows(), quad.points.cols());
99 Eigen::MatrixXd dz;
100 if (is_volume)
101 dz = Eigen::MatrixXd::Zero(quad.points.rows(), quad.points.cols());
102
103 gbasis.evaluate_grads(quad.points, tmp);
104
105 for (int j = 0; j < n_local_bases; ++j)
106 {
107 const basis::Basis &b = gbasis.bases[j];
108
109 for (std::size_t ii = 0; ii < b.global().size(); ++ii)
110 {
111 dx += tmp[j].grad.col(0) * b.global()[ii].node * b.global()[ii].val;
112 dy += tmp[j].grad.col(1) * b.global()[ii].node * b.global()[ii].val;
113 if (is_volume)
114 dz += tmp[j].grad.col(2) * b.global()[ii].node * b.global()[ii].val;
115 }
116 }
117
118 for (long i = 0; i < dx.rows(); ++i)
119 {
120 if (is_volume)
121 {
122 Eigen::Matrix3d tmp;
123 tmp << dx.row(i), dy.row(i), dz.row(i);
124 min_jacs = std::min(min_jacs, tmp.determinant());
125 }
126 else
127 {
128 Eigen::Matrix2d tmp;
129 tmp << dx.row(i), dy.row(i);
130 min_jacs = std::min(min_jacs, tmp.determinant());
131 }
132 }
133 }
134
135 return min_jacs;
136 }
137
138 void MinJacobianForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
139 {
140 log_and_throw_adjoint_error("{} is not differentiable!", name());
141 }
142
143 AMIPSForm::AMIPSForm(const VariableToSimulationGroup &variable_to_simulation, const State &state)
144 : AdjointForm(variable_to_simulation),
145 state_(state)
146 {
148 amips_energy_->set_size(state.mesh->dimension());
149
150 json transform_params = {};
151 transform_params["canonical_transformation"] = json::array();
152 if (!state.mesh->is_volume())
153 {
154 Eigen::MatrixXd regular_tri(3, 3);
155 regular_tri << 0, 0, 1,
156 1, 0, 1,
157 1. / 2., std::sqrt(3) / 2., 1;
158 regular_tri.transposeInPlace();
159 Eigen::MatrixXd regular_tri_inv = regular_tri.inverse();
160
161 const auto &mesh2d = *dynamic_cast<mesh::Mesh2D *>(state.mesh.get());
162 for (int e = 0; e < state.mesh->n_elements(); e++)
163 {
164 Eigen::MatrixXd transform;
165 mesh2d.compute_face_jacobian(e, regular_tri_inv, transform);
166 transform_params["canonical_transformation"].push_back(json({
167 {
168 transform(0, 0),
169 transform(0, 1),
170 },
171 {
172 transform(1, 0),
173 transform(1, 1),
174 },
175 }));
176 }
177 }
178 else
179 {
180 Eigen::MatrixXd regular_tet(4, 4);
181 regular_tet << 0, 0, 0, 1,
182 1, 0, 0, 1,
183 1. / 2., std::sqrt(3) / 2., 0, 1,
184 1. / 2., 1. / 2. / std::sqrt(3), std::sqrt(3) / 2., 1;
185 regular_tet.transposeInPlace();
186 Eigen::MatrixXd regular_tet_inv = regular_tet.inverse();
187
188 const auto &mesh3d = *dynamic_cast<mesh::Mesh3D *>(state.mesh.get());
189 for (int e = 0; e < state.mesh->n_elements(); e++)
190 {
191 Eigen::MatrixXd transform;
192 mesh3d.compute_cell_jacobian(e, regular_tet_inv, transform);
193 transform_params["canonical_transformation"].push_back(json({
194 {
195 transform(0, 0),
196 transform(0, 1),
197 transform(0, 2),
198 },
199 {
200 transform(1, 0),
201 transform(1, 1),
202 transform(1, 2),
203 },
204 {
205 transform(2, 0),
206 transform(2, 1),
207 transform(2, 2),
208 },
209 }));
210 }
211 }
212 transform_params["solve_displacement"] = true;
213 amips_energy_->add_multimaterial(0, transform_params, state.units);
214
215 Eigen::MatrixXd V;
220 }
221
222 double AMIPSForm::value_unweighted(const Eigen::VectorXd &x) const
223 {
224 Eigen::VectorXd X = get_updated_mesh_nodes(x);
225
226 return amips_energy_->assemble_energy(state_.mesh->is_volume(), init_geom_bases_, init_geom_bases_, init_ass_vals_cache_, 0, 0, AdjointTools::map_primitive_to_node_order(state_, X - X_rest), Eigen::VectorXd());
227 }
228
229 void AMIPSForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
230 {
232 const Eigen::VectorXd X = get_updated_mesh_nodes(x);
233 Eigen::MatrixXd grad;
234 amips_energy_->assemble_gradient(state_.mesh->is_volume(), state_.n_geom_bases, init_geom_bases_, init_geom_bases_, init_ass_vals_cache_, 0, 0, AdjointTools::map_primitive_to_node_order(state_, X - X_rest), Eigen::VectorXd(), grad); // grad wrt. gbases
236 });
237 }
238
239 bool AMIPSForm::is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
240 {
241 Eigen::VectorXd X = get_updated_mesh_nodes(x1);
242 Eigen::MatrixXd V1 = utils::unflatten(X, state_.mesh->dimension());
243 bool flipped = utils::is_flipped(V1, F);
244
245 if (flipped)
246 adjoint_logger().trace("[{}] Step flips elements.", name());
247
248 return !flipped;
249 }
250} // namespace polyfem::solver
int V
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
void get_vertices(Eigen::MatrixXd &vertices) const
Definition StateDiff.cpp:72
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:223
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:471
void get_elements(Eigen::MatrixXi &elements) const
Definition StateDiff.cpp:80
int n_geom_bases
number of geometric bases
Definition State.hpp:182
static std::shared_ptr< Assembler > make_assembler(const std::string &formulation)
Represents one basis function and its gradient.
Definition Basis.hpp:44
void compute_face_jacobian(const int el_id, const Eigen::MatrixXd &reference_map, Eigen::MatrixXd &jacobian) const
Definition Mesh2D.cpp:101
void compute_cell_jacobian(const int el_id, const Eigen::MatrixXd &reference_map, Eigen::MatrixXd &jacobian) const
Definition Mesh3D.cpp:440
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
virtual std::string name() const override
Definition AMIPSForm.hpp:45
bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine if a step from solution x0 to solution x1 is allowed.
Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const
Definition AMIPSForm.hpp:52
std::vector< polyfem::basis::ElementBases > init_geom_bases_
Definition AMIPSForm.hpp:63
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
std::shared_ptr< assembler::Assembler > amips_energy_
Definition AMIPSForm.hpp:65
assembler::AssemblyValsCache init_ass_vals_cache_
Definition AMIPSForm.hpp:64
AMIPSForm(const VariableToSimulationGroup &variable_to_simulation, const State &state)
const VariableToSimulationGroup variable_to_simulations_
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:128
virtual std::string name() const override
Definition AMIPSForm.hpp:31
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
Definition AMIPSForm.cpp:80
A collection of VariableToSimulation.
virtual Eigen::VectorXd apply_parametrization_jacobian(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, const std::function< Eigen::VectorXd()> &grad) const
Maps the partial gradient wrt.
int norm
Definition p_bases.py:264
bool scaled_jacobian(Mesh3DStorage &hmi, Mesh_Quality &mq)
Eigen::VectorXd map_primitive_to_node_order(const State &state, const Eigen::VectorXd &primitives)
Eigen::VectorXd map_node_to_primitive_order(const State &state, const Eigen::VectorXd &nodes)
bool is_flipped(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
Determine if any simplex is inverted or collapses.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
spdlog::logger & adjoint_logger()
Retrieves the current logger for adjoint.
Definition Logger.cpp:30
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:79