PolyFEM
Loading...
Searching...
No Matches
SplineParametrizations.cpp
Go to the documentation of this file.
3#include <polyfem/State.hpp>
5#include <igl/bbw.h>
6#include <igl/boundary_conditions.h>
7// #include <igl/normalize_row_sums.h>
8#include <igl/boundary_loop.h>
9#include <igl/exact_geodesic.h>
10#include <igl/bounding_box.h>
11#include <igl/writeOBJ.h>
12
14
15#include <unsupported/Eigen/SparseExtra>
16
17#include <unordered_map>
18
19namespace polyfem::solver
20{
21 namespace
22 {
23
24 Eigen::Matrix<double, 3, 1> affine_transformation(const Eigen::Matrix<double, 3, 1> &control_pt, const Eigen::Matrix<double, 3, 1> &point, const Eigen::Matrix<double, 6, 1> &param)
25 {
26 Eigen::Matrix<double, 3, 1> transformed_point(3);
27
28 const double helper_0 = cos(param(2));
29 const double helper_1 = control_pt(0) - point(0);
30 const double helper_2 = cos(param(1));
31 const double helper_3 = helper_1 * helper_2;
32 const double helper_4 = control_pt(2) - point(2);
33 const double helper_5 = sin(param(0));
34 const double helper_6 = sin(param(2));
35 const double helper_7 = helper_5 * helper_6;
36 const double helper_8 = sin(param(1));
37 const double helper_9 = cos(param(0));
38 const double helper_10 = helper_0 * helper_9;
39 const double helper_11 = control_pt(1) - point(1);
40 const double helper_12 = helper_6 * helper_9;
41 const double helper_13 = helper_0 * helper_5;
42 transformed_point(0) = control_pt(0) - helper_0 * helper_3 - helper_11 * (-helper_12 + helper_13 * helper_8) - helper_4 * (helper_10 * helper_8 + helper_7) + param(3);
43 transformed_point(1) = control_pt(1) - helper_11 * (helper_10 + helper_7 * helper_8) - helper_3 * helper_6 + helper_4 * (-helper_12 * helper_8 + helper_13) + param(4);
44 transformed_point(2) = control_pt(2) + helper_1 * helper_8 - helper_11 * helper_2 * helper_5 - helper_2 * helper_4 * helper_9 + param(5);
45
46 return transformed_point;
47 }
48
49 Eigen::MatrixXd grad_affine_transformation(const Eigen::Matrix<double, 3, 1> &control_pt, const Eigen::VectorXd &point, const Eigen::VectorXd &param)
50 {
51 Eigen::MatrixXd grad(6, 3);
52
53 const double helper_0 = control_pt(1) - point(1);
54 const double helper_1 = sin(param(0));
55 const double helper_2 = sin(param(2));
56 const double helper_3 = helper_1 * helper_2;
57 const double helper_4 = sin(param(1));
58 const double helper_5 = cos(param(0));
59 const double helper_6 = cos(param(2));
60 const double helper_7 = helper_5 * helper_6;
61 const double helper_8 = helper_3 + helper_4 * helper_7;
62 const double helper_9 = control_pt(2) - point(2);
63 const double helper_10 = helper_2 * helper_5;
64 const double helper_11 = helper_1 * helper_6;
65 const double helper_12 = -helper_10 + helper_11 * helper_4;
66 const double helper_13 = control_pt(0) - point(0);
67 const double helper_14 = cos(param(1));
68 const double helper_15 = helper_0 * helper_1;
69 const double helper_16 = helper_5 * helper_9;
70 const double helper_17 = helper_13 * helper_4 - helper_14 * helper_15 - helper_14 * helper_16;
71 const double helper_18 = helper_13 * helper_14;
72 const double helper_19 = helper_3 * helper_4 + helper_7;
73 const double helper_20 = -helper_10 * helper_4 + helper_11;
74 grad(0) = -helper_0 * helper_8 + helper_12 * helper_9;
75 grad(1) = helper_17 * helper_6;
76 grad(2) = helper_0 * helper_19 + helper_18 * helper_2 - helper_20 * helper_9;
77 grad(3) = 1;
78 grad(4) = 0;
79 grad(5) = 0;
80 grad(6) = helper_0 * helper_20 + helper_19 * helper_9;
81 grad(7) = helper_17 * helper_2;
82 grad(8) = -helper_0 * helper_12 - helper_18 * helper_6 - helper_8 * helper_9;
83 grad(9) = 0;
84 grad(10) = 1;
85 grad(11) = 0;
86 grad(12) = helper_14 * (-helper_0 * helper_5 + helper_1 * helper_9);
87 grad(13) = helper_15 * helper_4 + helper_16 * helper_4 + helper_18;
88 grad(14) = 0;
89 grad(15) = 0;
90 grad(16) = 0;
91 grad(17) = 1;
92
93 return grad.transpose();
94 }
95 } // namespace
96 Eigen::VectorXd BSplineParametrization1DTo2D::inverse_eval(const Eigen::VectorXd &y)
97 {
98 spline_ = std::make_shared<BSplineParametrization2D>(initial_control_points_, knots_, utils::unflatten(y, 2));
100 assert(size_ == spline_->vertex_size());
101 if (exclude_ends_)
102 return utils::flatten(initial_control_points_).segment(2, (initial_control_points_.rows() - 2) * 2);
103 else
105 }
106
107 Eigen::VectorXd BSplineParametrization1DTo2D::eval(const Eigen::VectorXd &x) const
108 {
110 log_and_throw_error("Must call inverse eval on this parametrization first!");
111 Eigen::MatrixXd new_control_points;
112 if (exclude_ends_)
113 {
114 new_control_points = initial_control_points_;
115 for (int i = 1; i < new_control_points.rows() - 1; ++i)
116 new_control_points.row(i) = x.segment(2 * i - 2, 2);
117 }
118 else
119 {
120 new_control_points = utils::unflatten(x, 2);
121 }
122 Eigen::MatrixXd new_vertices;
123 spline_->reparametrize(new_control_points, new_vertices);
124 Eigen::VectorXd y = utils::flatten(new_vertices);
125 return y;
126 }
127
128 Eigen::VectorXd BSplineParametrization1DTo2D::apply_jacobian(const Eigen::VectorXd &grad_full, const Eigen::VectorXd &x) const
129 {
130 Eigen::VectorXd grad;
131 spline_->derivative_wrt_params(grad_full, grad);
132 if (exclude_ends_)
133 return grad.segment(2, (initial_control_points_.rows() - 2) * 2);
134 else
135 return grad;
136 }
137
138 Eigen::VectorXd BSplineParametrization2DTo3D::inverse_eval(const Eigen::VectorXd &y)
139 {
140 spline_ = std::make_shared<BSplineParametrization3D>(initial_control_point_grid_, knots_u_, knots_v_, y);
142 return Eigen::VectorXd();
143 }
144
145 Eigen::VectorXd BSplineParametrization2DTo3D::eval(const Eigen::VectorXd &x) const
146 {
148 log_and_throw_error("Must call inverse eval on this parametrization first!");
149 return Eigen::VectorXd();
150 }
151
152 Eigen::VectorXd BSplineParametrization2DTo3D::apply_jacobian(const Eigen::VectorXd &grad_full, const Eigen::VectorXd &x) const
153 {
154 return Eigen::VectorXd();
155 }
156
157 BoundedBiharmonicWeights2Dto3D::BoundedBiharmonicWeights2Dto3D(const int num_control_vertices, const int num_vertices, const State &state, const bool allow_rotations)
158 : num_control_vertices_(num_control_vertices), num_vertices_(num_vertices), allow_rotations_(allow_rotations)
159 {
160 Eigen::MatrixXd V;
161 state.get_vertices(V);
162
163 auto map = state.node_to_primitive();
164
165 int f_size = 0;
166 const auto &mesh = state.mesh;
167 const auto &bases = state.bases;
168 const auto &gbases = state.geom_bases();
169 for (const auto &lb : state.total_local_boundary)
170 {
171 const int e = lb.element_id();
172 for (int i = 0; i < lb.size(); ++i)
173 {
174 const int primitive_global_id = lb.global_primitive_id(i);
175 const int boundary_id = mesh->get_boundary_id(primitive_global_id);
176 const auto nodes = gbases[e].local_nodes_for_primitive(primitive_global_id, *mesh);
177 F_surface_.conservativeResize(++f_size, 3);
178 for (int f = 0; f < nodes.size(); ++f)
179 {
180 F_surface_(f_size - 1, f) = map[gbases[e].bases[nodes(f)].global()[0].index];
181 }
182 }
183 }
184 V_surface_.resizeLike(V);
185 for (int e = 0; e < gbases.size(); e++)
186 {
187 for (const auto &gbs : gbases[e].bases)
188 V_surface_.row(map[gbs.global()[0].index]) = gbs.global()[0].node;
189 }
190 }
191
192 int BoundedBiharmonicWeights2Dto3D::optimal_new_control_point_idx(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::VectorXi &boundary_loop, const std::vector<int> &existing_points) const
193 {
194 std::set<int> fixed_vertices;
195 {
196 for (int i = 0; i < boundary_loop.size(); ++i)
197 fixed_vertices.insert(boundary_loop(i));
198 for (const auto &i : existing_points)
199 fixed_vertices.insert(i);
200 }
201
202 Eigen::VectorXi free_vertices(V.rows() - fixed_vertices.size());
203 int s = 0;
204 for (int i = 0; i < V.rows(); ++i)
205 if (fixed_vertices.find(i) == fixed_vertices.end())
206 free_vertices(s++) = i;
207
208 Eigen::VectorXi VS, FS, FT;
209 VS.resize(fixed_vertices.size());
210 s = 0;
211 for (const auto &j : fixed_vertices)
212 VS(s++) = j;
213
214 Eigen::VectorXd d;
215 igl::exact_geodesic(V, F, VS, FS, free_vertices, FT, d);
216 int opt_idx = -1;
217 double max_min_dist = -1;
218 for (int i = 0; i < d.size(); ++i)
219 {
220 if (opt_idx == -1)
221 {
222 max_min_dist = d(i);
223 opt_idx = free_vertices(i);
224 continue;
225 }
226 else if (d(i) > max_min_dist)
227 {
228 max_min_dist = d(i);
229 opt_idx = free_vertices(i);
230 }
231 }
232
233 return opt_idx;
234 }
235
236 Eigen::VectorXd BoundedBiharmonicWeights2Dto3D::inverse_eval(const Eigen::VectorXd &y)
237 {
238 y_start = y;
239
240 Eigen::MatrixXd V = utils::unflatten(y, 3);
241 Eigen::MatrixXi F;
243
244 Eigen::VectorXi outer_loop;
245 std::vector<std::vector<int>> loops;
246 igl::boundary_loop(F, loops);
247 if (loops.size() == 1)
248 {
249 outer_loop.resize(loops[0].size());
250 for (int i = 0; i < loops[0].size(); ++i)
251 outer_loop(i) = loops[0][i];
252 }
253 else
254 {
255 logger().error("More than 1 boundary loop! Concatenating and continuing as normal.");
256 for (const auto &l : loops)
257 {
258 int size = outer_loop.size();
259 outer_loop.conservativeResize(outer_loop.size() + l.size());
260 for (int i = 0; i < l.size(); ++i)
261 outer_loop(size + i) = l[i];
262 }
263 }
264 Eigen::MatrixXd V_outer_loop = V(outer_loop, Eigen::all);
265
266 Eigen::MatrixXd point_handles(num_control_vertices_ + outer_loop.size(), 3);
268 std::vector<int> control_indices;
269 {
270 std::set<int> possible_control_vertices;
271 for (int i = 0; i < F.rows(); ++i)
272 for (int j = 0; j < F.cols(); ++j)
273 possible_control_vertices.insert(F(i, j));
274 for (int i = 0; i < outer_loop.size(); ++i)
275 possible_control_vertices.erase(outer_loop(i));
276 for (int i = 0; i < num_control_vertices_; ++i)
277 control_indices.push_back(optimal_new_control_point_idx(V, F, outer_loop, control_indices));
278
279 const int recompute_loops = 5;
280 for (int r = 0; r < recompute_loops; ++r)
281 {
282 for (int i = 0; i < num_control_vertices_; ++i)
283 {
284 std::vector<int> indices = control_indices;
285 indices.erase(indices.begin() + i);
286 int new_idx = optimal_new_control_point_idx(V, F, outer_loop, indices);
287 control_indices[i] = new_idx;
288 }
289 }
290 }
291 for (int i = 0; i < num_control_vertices_; ++i)
292 control_points_.row(i) = V.row(control_indices[i]);
293 point_handles.block(0, 0, num_control_vertices_, 3) = control_points_;
294 point_handles.block(num_control_vertices_, 0, outer_loop.size(), 3) = V_outer_loop;
295
296 igl::writeOBJ("bbw_control_points_" + std::to_string(V.rows()) + ".obj", control_points_, Eigen::MatrixXi::Zero(0, 3));
297 igl::writeOBJ("bbw_outer_loop_" + std::to_string(V.rows()) + ".obj", V_outer_loop, Eigen::MatrixXi::Zero(0, 3));
298
299 Eigen::VectorXi b;
300 Eigen::MatrixXd bc;
301 Eigen::VectorXi point_handles_idx(point_handles.rows());
302 for (int i = 0; i < point_handles_idx.size(); ++i)
303 point_handles_idx(i) = i;
304 igl::boundary_conditions(V, F, point_handles, point_handles_idx, Eigen::VectorXi(), Eigen::VectorXi(), Eigen::VectorXi(), b, bc);
305
306 igl::BBWData bbw_data;
307 bbw_data.active_set_params.max_iter = 500;
308 bbw_data.verbosity = 1;
309 bbw_data.partition_unity = false; // Not implemented in libigl
310 Eigen::MatrixXd complete_bbw_weights;
311 bool computation = igl::bbw(V, F, b, bc, bbw_data, complete_bbw_weights);
312 if (!computation)
313 log_and_throw_error("Bounded Bihamonic Weight computation failed!");
314 // Deprecated: igl::normalize_row_sums(complete_bbw_weights, complete_bbw_weights);
315 complete_bbw_weights = (complete_bbw_weights.array().colwise() / complete_bbw_weights.array().rowwise().sum()).eval();
316 bbw_weights_ = complete_bbw_weights.block(0, 0, V.rows(), num_control_vertices_).matrix();
317 boundary_bbw_weights_ = complete_bbw_weights.block(0, num_control_vertices_, V.rows(), V_outer_loop.rows()).matrix();
318
319 igl::writeOBJ("surface_mesh.obj", V, F);
320 Eigen::saveMarket(bbw_weights_.sparseView(0, 1e-8).eval(), "bbw_control_weights.mat");
321 Eigen::saveMarket(boundary_bbw_weights_.sparseView(0, 1e-8).eval(), "bbw_boundary_weights.mat");
322
324
325 return Eigen::VectorXd::Zero(num_control_vertices_ * (allow_rotations_ ? 6 : 3));
326 }
327
328 Eigen::VectorXd BoundedBiharmonicWeights2Dto3D::eval(const Eigen::VectorXd &x) const
329 {
331 log_and_throw_error("Must call inverse eval on this parametrization first!");
332 Eigen::VectorXd y = Eigen::VectorXd::Zero(y_start.size());
334 {
335 for (int j = 0; j < bbw_weights_.cols(); ++j)
336 {
337 Eigen::Matrix<double, 6, 1> affine_params = x.segment(j * 6, 6);
338 for (int i = 0; i < bbw_weights_.rows(); ++i)
339 y.segment(i * 3, 3) += bbw_weights_(i, j) * affine_transformation(control_points_.row(j), y_start.segment(i * 3, 3), affine_params);
340 }
341 }
342 else
343 {
344 for (int j = 0; j < bbw_weights_.cols(); ++j)
345 for (int i = 0; i < bbw_weights_.rows(); ++i)
346 y.segment(i * 3, 3) += bbw_weights_(i, j) * (y_start.segment(i * 3, 3) + x.segment(j * 3, 3));
347 }
348
349 for (int j = 0; j < boundary_bbw_weights_.cols(); ++j)
350 for (int i = 0; i < boundary_bbw_weights_.rows(); ++i)
351 y.segment(i * 3, 3) += boundary_bbw_weights_(i, j) * y_start.segment(i * 3, 3);
352
353 return y;
354 }
355
356 Eigen::VectorXd BoundedBiharmonicWeights2Dto3D::apply_jacobian(const Eigen::VectorXd &grad_full, const Eigen::VectorXd &x) const
357 {
358 Eigen::VectorXd grad = Eigen::VectorXd::Zero(x.size());
360 {
361 for (int j = 0; j < bbw_weights_.cols(); ++j)
362 for (int i = 0; i < bbw_weights_.rows(); ++i)
363 grad.segment(j * 6, 6) += bbw_weights_(i, j) * grad_affine_transformation(control_points_.row(j), y_start.segment(i * 3, 3), x.segment(j * 6, 6)).transpose() * grad_full.segment(i * 3, 3);
364 }
365 else
366 {
367 for (int j = 0; j < bbw_weights_.cols(); ++j)
368 for (int i = 0; i < bbw_weights_.rows(); ++i)
369 grad.segment(j * 3, 3) += bbw_weights_(i, j) * grad_full.segment(i * 3, 3);
370 }
371 return grad;
372 }
373
374 void BoundedBiharmonicWeights2Dto3D::compute_faces_for_partial_vertices(const Eigen::MatrixXd &V, Eigen::MatrixXi &F) const
375 {
376 // The following implementation is maybe a bit wasteful, but is independent of state or surface selections
377 std::unordered_map<int, int> full_to_reduced_indices;
378
379 Eigen::MatrixXd BV;
380 BV.setZero(3, 2);
381 for (int i = 0; i < V.rows(); ++i)
382 {
383 if (i == 0)
384 {
385 BV.col(0) = V.row(i);
386 BV.col(1) = V.row(i);
387 }
388 else
389 {
390 for (int j = 0; j < 3; ++j)
391 {
392 BV(j, 0) = std::min(BV(j, 0), V(i, j));
393 BV(j, 1) = std::max(BV(j, 1), V(i, j));
394 }
395 }
396 }
397
398 Eigen::VectorXd bbox_width = (BV.col(1) - BV.col(0));
399 for (int j = 0; j < 3; ++j)
400 if (bbox_width(j) < 1e-12)
401 bbox_width(j) = 1e-3;
402
403 // Pad the bbox to make it conservative
404 BV.col(0) -= 0.05 * (BV.col(1) - BV.col(0));
405 BV.col(1) += 0.05 * (BV.col(1) - BV.col(0));
406
407 auto in_bbox = [&BV](const Eigen::VectorXd &x) {
408 bool in = true;
409 in &= (x(0) >= BV(0, 0)) && (x(0) <= BV(0, 1));
410 in &= (x(1) >= BV(1, 0)) && (x(1) <= BV(1, 1));
411 in &= (x(2) >= BV(2, 0)) && (x(0) <= BV(2, 1));
412 return in;
413 };
414
415 for (int i = 0; i < V_surface_.rows(); ++i)
416 for (int j = 0; j < V.rows(); ++j)
417 {
418 // if (in_bbox(V_surface_.row(i)))
419 if ((V_surface_.row(i) - V.row(j)).norm() < 1e-12)
420 full_to_reduced_indices[i] = j;
421 }
422
423 F.resize(0, 3);
424 for (int i = 0; i < F_surface_.rows(); ++i)
425 {
426 bool contains_face = true;
427 for (int j = 0; j < 3; ++j)
428 contains_face &= (full_to_reduced_indices.count(F_surface_(i, j)) == 1);
429
430 if (contains_face)
431 {
432 F.conservativeResize(F.rows() + 1, 3);
433 for (int j = 0; j < 3; ++j)
434 F(F.rows() - 1, j) = full_to_reduced_indices.at(F_surface_(i, j));
435 }
436 }
437 }
438} // namespace polyfem::solver
int V
int y
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:69
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:466
std::vector< int > node_to_primitive() const
Definition State.cpp:235
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
std::vector< mesh::LocalBoundary > total_local_boundary
mapping from elements to nodes for all mesh
Definition State.hpp:431
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad_full, const Eigen::VectorXd &x) const override
std::shared_ptr< BSplineParametrization2D > spline_
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad_full, const Eigen::VectorXd &x) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
std::shared_ptr< BSplineParametrization3D > spline_
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad_full, const Eigen::VectorXd &x) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
void compute_faces_for_partial_vertices(const Eigen::MatrixXd &V, Eigen::MatrixXi &F) const
BoundedBiharmonicWeights2Dto3D(const int num_control_vertices, const int num_vertices, const Eigen::MatrixXd &V_surface, const Eigen::MatrixXi &F_surface)
int optimal_new_control_point_idx(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::VectorXi &boundary_loop, const std::vector< int > &existing_points) const
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 & logger()
Retrieves the current logger.
Definition Logger.cpp:42
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71