PolyFEM
Loading...
Searching...
No Matches
TargetForms.cpp
Go to the documentation of this file.
1#include "TargetForms.hpp"
5
6#include <polyfem/State.hpp>
10
12
13using namespace polyfem::utils;
14
15namespace polyfem::solver
16{
17 namespace
18 {
19 class LocalThreadScalarStorage
20 {
21 public:
22 double val;
24
25 LocalThreadScalarStorage()
26 {
27 val = 0;
28 }
29 };
30 } // namespace
31
33 {
35 if (target_state_)
36 {
37 assert(target_state_->diff_cached.size() > 0);
38
39 auto j_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
40 int e_ref = params.elem;
41 if (auto search = e_to_ref_e_.find(params.elem); search != e_to_ref_e_.end())
42 e_ref = search->second;
43
44 Eigen::MatrixXd pts_ref;
45 target_state_->geom_bases()[e_ref].eval_geom_mapping(local_pts, pts_ref);
46
47 Eigen::MatrixXd u_ref, grad_u_ref;
48 const Eigen::VectorXd &sol_ref = target_state_->diff_cached.u(target_state_->problem->is_time_dependent() ? params.step : 0);
49 io::Evaluator::interpolate_at_local_vals(*(target_state_->mesh), target_state_->problem->is_scalar(), target_state_->bases, target_state_->geom_bases(), e_ref, local_pts, sol_ref, u_ref, grad_u_ref);
50
51 val = (u_ref + pts_ref - u - pts).rowwise().squaredNorm();
52 };
53
54 auto djdu_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
55 int e_ref = params.elem;
56 if (auto search = e_to_ref_e_.find(params.elem); search != e_to_ref_e_.end())
57 e_ref = search->second;
58
59 Eigen::MatrixXd pts_ref;
60 target_state_->geom_bases()[e_ref].eval_geom_mapping(local_pts, pts_ref);
61
62 Eigen::MatrixXd u_ref, grad_u_ref;
63 const Eigen::VectorXd &sol_ref = target_state_->diff_cached.u(target_state_->problem->is_time_dependent() ? params.step : 0);
64 io::Evaluator::interpolate_at_local_vals(*(target_state_->mesh), target_state_->problem->is_scalar(), target_state_->bases, target_state_->geom_bases(), e_ref, local_pts, sol_ref, u_ref, grad_u_ref);
65
66 val = 2 * (u + pts - u_ref - pts_ref);
67 };
68
69 j.set_j(j_func);
70 j.set_dj_du(djdu_func);
71 j.set_dj_dx(djdu_func); // only used for shape derivative
72 }
73 else if (have_target_func)
74 {
75 auto j_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
76 val.setZero(u.rows(), 1);
77
78 const Eigen::MatrixXd X = u + pts;
79 for (int q = 0; q < u.rows(); q++)
80 val(q) = target_func(X(q, 0), X(q, 1), X.cols() == 2 ? 0 : X(q, 2), 0, params.elem);
81 };
82
83 auto djdu_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
84 val.setZero(u.rows(), u.cols());
85
86 const Eigen::MatrixXd X = u + pts;
87 for (int q = 0; q < u.rows(); q++)
88 for (int d = 0; d < val.cols(); d++)
89 val(q, d) = target_func_grad[d](X(q, 0), X(q, 1), X.cols() == 2 ? 0 : X(q, 2), 0, params.elem);
90 };
91
92 j.set_j(j_func);
93 j.set_dj_du(djdu_func);
94 j.set_dj_dx(djdu_func); // only used for shape derivative
95 }
96 else // error wrt. a constant displacement
97 {
98 if (target_disp.size() == state_.mesh->dimension())
99 {
100 auto j_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
101 val.setZero(u.rows(), 1);
102
103 for (int q = 0; q < u.rows(); q++)
104 {
105 Eigen::VectorXd err = u.row(q) - this->target_disp.transpose();
106 for (int d = 0; d < active_dimension_mask.size(); d++)
107 if (!active_dimension_mask[d])
108 err(d) = 0;
109 val(q) = err.squaredNorm();
110 }
111 };
112 auto djdu_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
113 val.setZero(u.rows(), u.cols());
114
115 for (int q = 0; q < u.rows(); q++)
116 {
117 Eigen::VectorXd err = u.row(q) - this->target_disp.transpose();
118 for (int d = 0; d < active_dimension_mask.size(); d++)
119 if (!active_dimension_mask[d])
120 err(d) = 0;
121 val.row(q) = 2 * err;
122 }
123 };
124
125 j.set_j(j_func);
126 j.set_dj_du(djdu_func);
127 }
128 else
129 log_and_throw_adjoint_error("[{}] Only constant target displacement is supported!", name());
130 }
131
132 return j;
133 }
134
135 void TargetForm::set_reference(const std::shared_ptr<const State> &target_state, const std::set<int> &reference_cached_body_ids)
136 {
137 target_state_ = target_state;
138
139 std::map<int, std::vector<int>> ref_interested_body_id_to_e;
140 int ref_count = 0;
141 for (int e = 0; e < target_state_->bases.size(); ++e)
142 {
143 int body_id = target_state_->mesh->get_body_id(e);
144 if (reference_cached_body_ids.size() > 0 && reference_cached_body_ids.count(body_id) == 0)
145 continue;
146 if (ref_interested_body_id_to_e.find(body_id) != ref_interested_body_id_to_e.end())
147 ref_interested_body_id_to_e[body_id].push_back(e);
148 else
149 ref_interested_body_id_to_e[body_id] = {e};
150 ref_count++;
151 }
152
153 std::map<int, std::vector<int>> interested_body_id_to_e;
154 int count = 0;
155 for (int e = 0; e < state_.bases.size(); ++e)
156 {
157 int body_id = state_.mesh->get_body_id(e);
158 if (reference_cached_body_ids.size() > 0 && reference_cached_body_ids.count(body_id) == 0)
159 continue;
160 if (interested_body_id_to_e.find(body_id) != interested_body_id_to_e.end())
161 interested_body_id_to_e[body_id].push_back(e);
162 else
163 interested_body_id_to_e[body_id] = {e};
164 count++;
165 }
166
167 if (count != ref_count)
168 adjoint_logger().error("[{}] Number of interested elements in the reference and optimization examples do not match! {} {}", name(), count, ref_count);
169 else
170 adjoint_logger().trace("[{}] Found {} matching elements.", name(), count);
171
172 for (const auto &kv : interested_body_id_to_e)
173 {
174 for (int i = 0; i < kv.second.size(); ++i)
175 {
176 e_to_ref_e_[kv.second[i]] = ref_interested_body_id_to_e[kv.first][i];
177 }
178 }
179 }
180
181 void TargetForm::set_reference(const json &func, const json &grad_func)
182 {
183 target_func.init(func);
184 for (size_t k = 0; k < grad_func.size(); k++)
185 target_func_grad[k].init(grad_func[k]);
186 have_target_func = true;
187 }
188
189 void SDFTargetForm::solution_changed_step(const int time_step, const Eigen::VectorXd &x)
190 {
191 const auto &bases = state_.bases;
192 const auto &gbases = state_.geom_bases();
193 const int actual_dim = state_.problem->is_scalar() ? 1 : dim;
194
195 auto storage = utils::create_thread_storage(LocalThreadScalarStorage());
196 utils::maybe_parallel_for(state_.total_local_boundary.size(), [&](int start, int end, int thread_id) {
197 LocalThreadScalarStorage &local_storage = utils::get_local_thread_storage(storage, thread_id);
198
199 Eigen::MatrixXd uv, samples, gtmp;
200 Eigen::MatrixXd points, normal;
201 Eigen::VectorXd weights;
202
203 Eigen::MatrixXd u, grad_u;
204
205 for (int lb_id = start; lb_id < end; ++lb_id)
206 {
207 const auto &lb = state_.total_local_boundary[lb_id];
208 const int e = lb.element_id();
209
210 for (int i = 0; i < lb.size(); i++)
211 {
212 const int global_primitive_id = lb.global_primitive_id(i);
213 if (ids_.size() != 0 && ids_.find(state_.mesh->get_boundary_id(global_primitive_id)) == ids_.end())
214 continue;
215
216 utils::BoundarySampler::boundary_quadrature(lb, state_.n_boundary_samples(), *state_.mesh, i, false, uv, points, normal, weights);
217
218 assembler::ElementAssemblyValues &vals = local_storage.vals;
219 vals.compute(e, state_.mesh->is_volume(), points, bases[e], gbases[e]);
220 io::Evaluator::interpolate_at_local_vals(e, dim, actual_dim, vals, state_.diff_cached.u(time_step), u, grad_u);
221
222 normal = normal * vals.jac_it[0]; // assuming linear geometry
223
224 for (int q = 0; q < u.rows(); q++)
225 interpolation_fn->cache_grid([this](const Eigen::MatrixXd &point, double &distance) { compute_distance(point, distance); }, vals.val.row(q) + u.row(q));
226 }
227 }
228 });
229 }
230
231 void SDFTargetForm::set_bspline_target(const Eigen::MatrixXd &control_points, const Eigen::VectorXd &knots, const double delta)
232 {
233 dim = control_points.cols();
234 delta_ = delta;
235 if ((dim != 2) || (state_.mesh->dimension() != 2))
236 log_and_throw_error("SDFTargetForm specified for 2d.");
237
238 samples = 100;
239
240 nanospline::BSpline<double, 2, 3> curve;
241 curve.set_control_points(control_points);
242 curve.set_knots(knots);
243
244 t_or_uv_sampling = Eigen::VectorXd::LinSpaced(samples, 0, 1);
245 point_sampling.setZero(samples, 2);
246 for (int i = 0; i < t_or_uv_sampling.size(); ++i)
247 point_sampling.row(i) = curve.evaluate(t_or_uv_sampling(i));
248
249 {
250 Eigen::MatrixXi edges(samples - 1, 2);
251 edges.col(0) = Eigen::VectorXi::LinSpaced(samples - 1, 0, samples - 2);
252 edges.col(1) = Eigen::VectorXi::LinSpaced(samples - 1, 1, samples - 1);
253 io::OBJWriter::write(fmt::format("spline_target_{:d}.obj", rand() % 100), point_sampling, edges);
254 }
255
256 interpolation_fn = std::make_unique<LazyCubicInterpolator>(dim, delta_);
257 }
258
259 void SDFTargetForm::set_bspline_target(const Eigen::MatrixXd &control_points, const Eigen::VectorXd &knots_u, const Eigen::VectorXd &knots_v, const double delta)
260 {
261
262 dim = control_points.cols();
263 delta_ = delta;
264 if ((dim != 3) || (state_.mesh->dimension() != 3))
265 log_and_throw_error("SDFTargetForm specified for 3d.");
266
267 samples = 100;
268
269 nanospline::BSplinePatch<double, 3, 3, 3> patch;
270 patch.set_control_grid(control_points);
271 patch.set_knots_u(knots_u);
272 patch.set_knots_v(knots_v);
273 patch.initialize();
274
275 t_or_uv_sampling.resize(samples * samples, 2);
276 for (int i = 0; i < samples; ++i)
277 {
278 t_or_uv_sampling.block(i * samples, 0, samples, 1) = Eigen::VectorXd::LinSpaced(samples, 0, 1);
279 t_or_uv_sampling.block(i * samples, 1, samples, 1) = (double)i / (samples - 1) * Eigen::VectorXd::Ones(samples);
280 }
281 point_sampling.setZero(samples * samples, 3);
282 for (int i = 0; i < t_or_uv_sampling.rows(); ++i)
283 {
284 point_sampling.row(i) = patch.evaluate(t_or_uv_sampling(i, 0), t_or_uv_sampling(i, 1));
285 }
286
287 {
288 Eigen::MatrixXi F(2 * ((samples - 1) * (samples - 1)), 3);
289 int f = 0;
290 for (int i = 0; i < samples - 1; ++i)
291 for (int j = 0; j < samples - 1; ++j)
292 {
293 Eigen::MatrixXi F_local(2, 3);
294 F_local << (i * samples + j), ((i + 1) * samples + j), (i * samples + j + 1),
295 (i * samples + j + 1), ((i + 1) * samples + j), ((i + 1) * samples + j + 1);
296 F.block(f, 0, 2, 3) = F_local;
297 f += 2;
298 }
299 io::OBJWriter::write(fmt::format("spline_target_{:d}.obj", rand() % 100), point_sampling, F);
300 }
301
302 interpolation_fn = std::make_unique<LazyCubicInterpolator>(dim, delta_);
303 }
304
305 void SDFTargetForm::compute_distance(const Eigen::MatrixXd &point, double &distance) const
306 {
307 distance = DBL_MAX;
308 Eigen::MatrixXd p = point.transpose();
309
310 if (dim == 2)
311 for (int i = 0; i < t_or_uv_sampling.size() - 1; ++i)
312 {
313 const double l = (point_sampling.row(i + 1) - point_sampling.row(i)).squaredNorm();
314 double distance_to_perpendicular = ((p - point_sampling.row(i)) * (point_sampling.row(i + 1) - point_sampling.row(i)).transpose())(0) / l;
315 const double t = std::max(0., std::min(1., distance_to_perpendicular));
316 const auto project = point_sampling.row(i) * (1 - t) + point_sampling.row(i + 1) * t;
317 const double project_distance = (p - project).norm();
318 if (project_distance < distance)
319 distance = project_distance;
320 }
321 else if (dim == 3)
322 {
323 for (int i = 0; i < samples - 1; ++i)
324 for (int j = 0; j < samples - 1; ++j)
325 {
326 int loc = samples * i + j;
327 const double l1 = (point_sampling.row(loc + 1) - point_sampling.row(loc)).squaredNorm();
328 double distance_to_perpendicular = ((p - point_sampling.row(loc)) * (point_sampling.row(loc + 1) - point_sampling.row(loc)).transpose())(0) / l1;
329 const double u = std::max(0., std::min(1., distance_to_perpendicular));
330
331 const double l2 = (point_sampling.row(loc + samples) - point_sampling.row(loc)).squaredNorm();
332 distance_to_perpendicular = ((p - point_sampling.row(loc)) * (point_sampling.row(loc + samples) - point_sampling.row(loc)).transpose())(0) / l2;
333 const double v = std::max(0., std::min(1., distance_to_perpendicular));
334
335 Eigen::MatrixXd project = point_sampling.row(loc) * (1 - u) + point_sampling.row(loc + 1) * u;
336 project += v * (point_sampling.row(loc + samples) - point_sampling.row(loc));
337 const double project_distance = (p - project).norm();
338 if (project_distance < distance)
339 distance = project_distance;
340 }
341 }
342 }
343
345 {
347 auto j_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
348 val.setZero(u.rows(), 1);
349
350 for (int q = 0; q < u.rows(); q++)
351 {
352 double distance;
353 Eigen::MatrixXd unused_grad;
354 interpolation_fn->evaluate(u.row(q) + pts.row(q), distance, unused_grad);
355 val(q) = pow(distance, 2);
356 }
357 };
358
359 auto djdu_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
360 val.setZero(u.rows(), u.cols());
361
362 for (int q = 0; q < u.rows(); q++)
363 {
364 double distance;
365 Eigen::MatrixXd grad;
366 interpolation_fn->evaluate(u.row(q) + pts.row(q), distance, grad);
367 val.row(q) = 2 * distance * grad.transpose();
368 }
369 };
370
371 j.set_j(j_func);
372 j.set_dj_du(djdu_func);
373 j.set_dj_dx(djdu_func);
374
375 return j;
376 }
377
378 void MeshTargetForm::set_surface_mesh_target(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const double delta)
379 {
380 dim = V.cols();
381 delta_ = delta;
382 if ((dim != 3) || (state_.mesh->dimension() != 3))
383 log_and_throw_error("MeshTargetForm is only available for 3d scenes.");
384
385 tree_.init(V, F);
386 V_ = V;
387 F_ = F;
388
389 interpolation_fn = std::make_unique<LazyCubicInterpolator>(dim, delta_);
390 }
391
392 void MeshTargetForm::solution_changed_step(const int time_step, const Eigen::VectorXd &x)
393 {
394 const auto &bases = state_.bases;
395 const auto &gbases = state_.geom_bases();
396 const int actual_dim = state_.problem->is_scalar() ? 1 : dim;
397
398 auto storage = utils::create_thread_storage(LocalThreadScalarStorage());
399 utils::maybe_parallel_for(state_.total_local_boundary.size(), [&](int start, int end, int thread_id) {
400 LocalThreadScalarStorage &local_storage = utils::get_local_thread_storage(storage, thread_id);
401
402 Eigen::MatrixXd uv, samples, gtmp;
403 Eigen::MatrixXd points, normal;
404 Eigen::VectorXd weights;
405
406 Eigen::MatrixXd u, grad_u;
407
408 for (int lb_id = start; lb_id < end; ++lb_id)
409 {
410 const auto &lb = state_.total_local_boundary[lb_id];
411 const int e = lb.element_id();
412
413 for (int i = 0; i < lb.size(); i++)
414 {
415 const int global_primitive_id = lb.global_primitive_id(i);
416 if (ids_.size() != 0 && ids_.find(state_.mesh->get_boundary_id(global_primitive_id)) == ids_.end())
417 continue;
418
419 utils::BoundarySampler::boundary_quadrature(lb, state_.n_boundary_samples(), *state_.mesh, i, false, uv, points, normal, weights);
420
421 assembler::ElementAssemblyValues &vals = local_storage.vals;
422 vals.compute(e, state_.mesh->is_volume(), points, bases[e], gbases[e]);
423 io::Evaluator::interpolate_at_local_vals(e, dim, actual_dim, vals, state_.diff_cached.u(time_step), u, grad_u);
424
425 normal = normal * vals.jac_it[0]; // assuming linear geometry
426
427 for (int q = 0; q < u.rows(); q++)
428 interpolation_fn->cache_grid([this](const Eigen::MatrixXd &point, double &distance) {
429 int idx;
430 Eigen::Matrix<double, 1, 3> closest;
431 distance = pow(tree_.squared_distance(V_, F_, point.col(0), idx, closest), 0.5);
432 },
433 vals.val.row(q) + u.row(q));
434 }
435 }
436 });
437 }
438
440 {
442 auto j_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
443 val.setZero(u.rows(), 1);
444
445 for (int q = 0; q < u.rows(); q++)
446 {
447 double distance;
448 Eigen::MatrixXd unused_grad;
449 interpolation_fn->evaluate(u.row(q) + pts.row(q), distance, unused_grad);
450 val(q) = pow(distance, 2);
451 }
452 };
453
454 auto djdu_func = [this](const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::VectorXd &lambda, const Eigen::VectorXd &mu, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, const IntegrableFunctional::ParameterType &params, Eigen::MatrixXd &val) {
455 val.setZero(u.rows(), u.cols());
456
457 for (int q = 0; q < u.rows(); q++)
458 {
459 double distance;
460 Eigen::MatrixXd grad;
461 interpolation_fn->evaluate(u.row(q) + pts.row(q), distance, grad);
462 val.row(q) = 2 * distance * grad.transpose();
463 }
464 };
465
466 j.set_j(j_func);
467 j.set_dj_du(djdu_func);
468 j.set_dj_dx(djdu_func);
469
470 return j;
471 }
472
473 NodeTargetForm::NodeTargetForm(const State &state, const VariableToSimulationGroup &variable_to_simulations, const json &args) : StaticForm(variable_to_simulations), state_(state)
474 {
475 const int dim = state.mesh->dimension();
476 const std::string target_data_path = args["target_data_path"];
477 if (!std::filesystem::is_regular_file(target_data_path))
478 {
479 throw std::runtime_error("Marker path invalid!");
480 }
481 // N x (dim * 2), each row is [rest_x, rest_y, rest_z, deform_x, deform_y, deform_z]
482 Eigen::MatrixXd data;
483 io::read_matrix(target_data_path, data);
484
485 // markers to nodes
486 target_vertex_positions.setZero(data.rows(), dim);
487 active_nodes.reserve(data.rows());
488 for (int s = 0; s < data.rows(); s++)
489 {
490 target_vertex_positions.row(s) = data.block(s, dim, 1, dim);
491
492 const RowVectorNd node = data.block(s, 0, 1, dim);
493 bool not_found = true;
494 double min_dist = std::numeric_limits<double>::max();
495 for (int v = 0; v < state_.mesh_nodes->n_nodes(); v++)
496 {
497 min_dist = std::min(min_dist, (state_.mesh_nodes->node_position(v) - node).norm());
498 if ((state_.mesh_nodes->node_position(v) - node).norm() < args["tolerance"])
499 {
500 active_nodes.push_back(v);
501 not_found = false;
502 break;
503 }
504 }
505 if (not_found)
506 log_and_throw_adjoint_error("Failed to find corresponding node for {}! Minimum distance {}", node, min_dist);
507 }
508 }
509
510 NodeTargetForm::NodeTargetForm(const State &state, const VariableToSimulationGroup &variable_to_simulations, const std::vector<int> &active_nodes_, const Eigen::MatrixXd &target_vertex_positions_) : StaticForm(variable_to_simulations), state_(state), target_vertex_positions(target_vertex_positions_), active_nodes(active_nodes_)
511 {
512 // log_and_throw_adjoint_error("[{}] Constructor not implemented!", name());
513 }
514
515 Eigen::VectorXd NodeTargetForm::compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const
516 {
517 Eigen::VectorXd rhs;
518 rhs.setZero(state.diff_cached.u(0).size());
519
520 const int dim = state_.mesh->dimension();
521
522 if (&state == &state_)
523 {
524 int i = 0;
525 const Eigen::VectorXd disp = state_.diff_cached.u(time_step);
526 for (int v : active_nodes)
527 {
528 const RowVectorNd cur_pos = state_.mesh_nodes->node_position(v) + disp.segment(v * dim, dim).transpose();
529
530 rhs.segment(v * dim, dim) = 2 * (cur_pos - target_vertex_positions.row(i++));
531 }
532 }
533
534 return rhs * weight();
535 }
536
537 double NodeTargetForm::value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const
538 {
539 const int dim = state_.mesh->dimension();
540 double val = 0;
541 int i = 0;
542 const Eigen::VectorXd disp = state_.diff_cached.u(time_step);
543 for (int v : active_nodes)
544 {
545 const RowVectorNd cur_pos = state_.mesh_nodes->node_position(v) + disp.segment(v * dim, dim).transpose();
546 val += (cur_pos - target_vertex_positions.row(i++)).squaredNorm();
547 }
548 return val;
549 }
550
551 void NodeTargetForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
552 {
553 gradv.setZero(x.size());
555 log_and_throw_adjoint_error("[{}] Doesn't support derivatives wrt. shape!", name());
556 return Eigen::VectorXd::Zero(0).eval();
557 });
558 }
559
560 BarycenterTargetForm::BarycenterTargetForm(const VariableToSimulationGroup &variable_to_simulations, const json &args, const std::shared_ptr<State> &state1, const std::shared_ptr<State> &state2) : StaticForm(variable_to_simulations)
561 {
562 dim = state1->mesh->dimension();
563 json tmp_args = args;
564 for (int d = 0; d < dim; d++)
565 {
566 tmp_args["dim"] = d;
567 center1.push_back(std::make_unique<PositionForm>(variable_to_simulations, *state1, tmp_args));
568 center2.push_back(std::make_unique<PositionForm>(variable_to_simulations, *state2, tmp_args));
569 }
570 }
571
572 Eigen::VectorXd BarycenterTargetForm::compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const
573 {
574 Eigen::VectorXd term;
575 term.setZero(state.ndof());
576 for (int d = 0; d < dim; d++)
577 {
578 double value = center1[d]->value_unweighted_step(time_step, x) - center2[d]->value_unweighted_step(time_step, x);
579 term += (2 * value) * (center1[d]->compute_adjoint_rhs_step(time_step, x, state) - center2[d]->compute_adjoint_rhs_step(time_step, x, state));
580 }
581 return term * weight();
582 }
583 void BarycenterTargetForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
584 {
585 gradv.setZero(x.size());
586 Eigen::VectorXd tmp1, tmp2;
587 for (int d = 0; d < dim; d++)
588 {
589 double value = center1[d]->value_unweighted_step(time_step, x) - center2[d]->value_unweighted_step(time_step, x);
590 center1[d]->compute_partial_gradient_step(time_step, x, tmp1);
591 center2[d]->compute_partial_gradient_step(time_step, x, tmp2);
592 gradv += (2 * value) * (tmp1 - tmp2);
593 }
594 gradv *= weight();
595 }
596 double BarycenterTargetForm::value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const
597 {
598 double dist = 0;
599 for (int d = 0; d < dim; d++)
600 dist += std::pow(center1[d]->value_unweighted_step(time_step, x) - center2[d]->value_unweighted_step(time_step, x), 2);
601
602 return dist;
603 }
604
605 MinTargetDistForm::MinTargetDistForm(const VariableToSimulationGroup &variable_to_simulations, const std::vector<int> &steps, const Eigen::VectorXd &target, const json &args, const std::shared_ptr<State> &state)
606 : AdjointForm(variable_to_simulations), steps_(steps), target_(target)
607 {
608 dim = state->mesh->dimension();
609 json tmp_args = args;
610 for (int d = 0; d < dim; d++)
611 {
612 tmp_args["dim"] = d;
613 objs.push_back(std::make_unique<PositionForm>(variable_to_simulations, *state, tmp_args));
614 }
615 objs.push_back(std::make_unique<VolumeForm>(variable_to_simulations, *state, args));
616 }
617 Eigen::MatrixXd MinTargetDistForm::compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const
618 {
619 Eigen::VectorXd values(steps_.size());
620 std::vector<Eigen::MatrixXd> grads(steps_.size(), Eigen::MatrixXd::Zero(state.ndof(), objs.size()));
621 Eigen::MatrixXd g2(steps_.size(), objs.size());
622 int i = 0;
623 for (int s : steps_)
624 {
625 Eigen::VectorXd input(objs.size());
626 Eigen::VectorXd tmp;
627 for (int d = 0; d < objs.size(); d++)
628 {
629 input(d) = objs[d]->value_unweighted_step(s, x);
630 grads[i].col(d) = objs[d]->compute_adjoint_rhs_step(s, x, state);
631 }
632 values[i] = eval2(input);
633 g2.row(i++) = eval2_grad(input);
634 }
635
636 Eigen::VectorXd g1 = eval1_grad(values);
637 Eigen::MatrixXd terms = Eigen::MatrixXd::Zero(state.ndof(), state.diff_cached.size());
638 i = 0;
639 for (int s : steps_)
640 {
641 terms.col(s) += g1(i) * grads[i] * g2.row(i).transpose();
642 i++;
643 }
644
645 return terms * weight();
646 }
647 void MinTargetDistForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
648 {
649 Eigen::VectorXd values(steps_.size());
650 std::vector<Eigen::MatrixXd> grads(steps_.size(), Eigen::MatrixXd::Zero(x.size(), objs.size()));
651 Eigen::MatrixXd g2(steps_.size(), objs.size());
652 int i = 0;
653 for (int s : steps_)
654 {
655 Eigen::VectorXd input(objs.size());
656 Eigen::VectorXd tmp;
657 for (int d = 0; d < objs.size(); d++)
658 {
659 input(d) = objs[d]->value_unweighted_step(s, x);
660 objs[d]->compute_partial_gradient_step(s, x, tmp);
661 grads[i].col(d) = tmp;
662 }
663 values[i] = eval2(input);
664 g2.row(i++) = eval2_grad(input);
665 }
666
667 Eigen::VectorXd g1 = eval1_grad(values);
668 gradv.setZero(x.size());
669 i = 0;
670 for (int s : steps_)
671 {
672 gradv += g1(i) * grads[i] * g2.row(i).transpose();
673 i++;
674 }
675 gradv *= weight();
676 }
677 double MinTargetDistForm::value_unweighted(const Eigen::VectorXd &x) const
678 {
679 Eigen::VectorXd values(steps_.size());
680 int i = 0;
681 for (int s : steps_)
682 {
683 Eigen::VectorXd input(objs.size());
684 for (int d = 0; d < objs.size(); d++)
685 input(d) = objs[d]->value_unweighted_step(s, x);
686 values[i++] = eval2(input);
687 }
688
689 return eval1(values);
690 }
691} // namespace polyfem::solver
int V
double val
Definition Assembler.cpp:86
ElementAssemblyValues vals
Definition Assembler.cpp:22
int x
void set_dj_du(const functionalType &dj_du)
void evaluate(const Eigen::MatrixXd &elastic_params, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &u, const Eigen::MatrixXd &grad_u, const Eigen::MatrixXd &reference_normals, const assembler::ElementAssemblyValues &vals, ParameterType &params, Eigen::MatrixXd &val) const
void set_j(const functionalType &j)
void set_dj_dx(const functionalType &dj_dx)
main class that contains the polyfem solver and all its state
Definition State.hpp:79
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::shared_ptr< polyfem::mesh::MeshNodes > mesh_nodes
Mapping from input nodes to FE nodes.
Definition State.hpp:193
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
solver::DiffCache diff_cached
Definition State.hpp:649
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
int ndof() const
Definition State.hpp:653
stores per element basis values at given quadrature points and geometric mapping
static void interpolate_at_local_vals(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const int el_index, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &fun, Eigen::MatrixXd &result, Eigen::MatrixXd &result_grad)
interpolate solution and gradient at element (calls interpolate_at_local_vals with sol)
static bool write(const std::string &path, const Eigen::MatrixXd &v, const Eigen::MatrixXi &e, const Eigen::MatrixXi &f)
Definition OBJWriter.cpp:18
double value(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied with the weigth.
const VariableToSimulationGroup variable_to_simulations_
std::vector< std::unique_ptr< PositionForm > > center1
std::vector< std::unique_ptr< PositionForm > > center2
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override
BarycenterTargetForm(const VariableToSimulationGroup &variable_to_simulations, const json &args, const std::shared_ptr< State > &state1, const std::shared_ptr< State > &state2)
Eigen::VectorXd u(int step) const
virtual void init(const Eigen::VectorXd &x)
Initialize the form.
Definition Form.hpp:19
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:128
void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) override
IntegrableFunctional get_integral_functional() const override
igl::AABB< Eigen::MatrixXd, 3 > tree_
void set_surface_mesh_target(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const double delta)
std::unique_ptr< LazyCubicInterpolator > interpolation_fn
MinTargetDistForm(const VariableToSimulationGroup &variable_to_simulations, const std::vector< int > &steps, const Eigen::VectorXd &target, const json &args, const std::shared_ptr< State > &state)
static double eval1(Eigen::VectorXd &x)
double eval2(Eigen::VectorXd &x) const
std::vector< std::unique_ptr< StaticForm > > objs
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.
const std::vector< int > steps_
static Eigen::VectorXd eval1_grad(Eigen::VectorXd &x)
Eigen::VectorXd eval2_grad(Eigen::VectorXd &x) const
Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const override
std::string name() const override
NodeTargetForm(const State &state, const VariableToSimulationGroup &variable_to_simulations, const json &args)
double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Eigen::MatrixXd target_vertex_positions
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override
void compute_distance(const Eigen::MatrixXd &point, double &distance) const
void set_bspline_target(const Eigen::MatrixXd &control_points, const Eigen::VectorXd &knots, const double delta)
std::unique_ptr< LazyCubicInterpolator > interpolation_fn
void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) override
IntegrableFunctional get_integral_functional() const override
std::shared_ptr< const State > target_state_
void set_reference(const std::shared_ptr< const State > &target_state, const std::set< int > &reference_cached_body_ids)
IntegrableFunctional get_integral_functional() const override
std::array< utils::ExpressionValue, 3 > target_func_grad
utils::ExpressionValue target_func
std::map< int, int > e_to_ref_e_
std::vector< bool > active_dimension_mask
virtual std::string name() const override
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.
bool read_matrix(const std::string &path, Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &mat)
Reads a matrix from a file. Determines the file format based on the path's extension.
Definition MatrixIO.cpp:18
auto create_thread_storage(const LocalStorage &initial_local_storage)
void maybe_parallel_for(int size, const std::function< void(int, int, int)> &partial_for)
spdlog::logger & adjoint_logger()
Retrieves the current logger for adjoint.
Definition Logger.cpp:28
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
Parameters for the functional evaluation.