PolyFEM
Loading...
Searching...
No Matches
BarrierForms.cpp
Go to the documentation of this file.
2
3#include <polyfem/State.hpp>
4#include <polyfem/Common.hpp>
9
10#include <Eigen/Core>
11
12#include <cassert>
13#include <cstddef>
14#include <memory>
15#include <string>
16#include <utility>
17#include <vector>
18
19namespace polyfem::solver
20{
21 namespace
22 {
23 class QuadraticBarrier : public ipc::Barrier
24 {
25 public:
26 QuadraticBarrier(const double weight = 1) : weight_(weight) {}
27
28 double operator()(const double d, const double dhat) const override
29 {
30 if (d > dhat)
31 return 0;
32 else
33 return weight_ * (d - dhat) * (d - dhat);
34 }
35 double first_derivative(const double d, const double dhat) const override
36 {
37 if (d > dhat)
38 return 0;
39 else
40 return 2 * weight_ * (d - dhat);
41 }
42 double second_derivative(const double d, const double dhat) const override
43 {
44 if (d > dhat)
45 return 0;
46 else
47 return 2 * weight_;
48 }
49 double units(const double dhat) const override
50 {
51 return dhat * dhat;
52 }
53
54 private:
55 const double weight_;
56 };
57
58 } // namespace
59
60 CollisionBarrierForm::CollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, std::shared_ptr<const State> state, const double dhat, const double dmin)
61 : AdjointForm(variable_to_simulation), state_(std::move(state)), dhat_(dhat), dmin_(dmin), barrier_potential_(dhat)
62 {
64 *state_->mesh, state_->n_geom_bases, state_->geom_bases(), state_->geom_bases(),
65 state_->total_local_boundary, state_->obstacle, state_->args,
66 [this](const std::string &p) { return this->state_->resolve_input_path(p); },
67 state_->in_node_to_node, collision_mesh_);
68
69 Eigen::MatrixXd V;
70 state_->get_vertices(V);
72
73 broad_phase_method_ = ipc::BroadPhaseMethod::HASH_GRID;
74 }
75
76 double CollisionBarrierForm::value_unweighted(const Eigen::VectorXd &x) const
77 {
78 const Eigen::MatrixXd displaced_surface = collision_mesh_.vertices(utils::unflatten(get_updated_mesh_nodes(x), state_->mesh->dimension()));
79 return barrier_potential_(collision_set, collision_mesh_, displaced_surface);
80 }
81
82 void CollisionBarrierForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
83 {
85 const Eigen::MatrixXd displaced_surface = collision_mesh_.vertices(utils::unflatten(get_updated_mesh_nodes(x), state_->mesh->dimension()));
86 const Eigen::VectorXd grad = collision_mesh_.to_full_dof(barrier_potential_.gradient(collision_set, collision_mesh_, displaced_surface));
87 return AdjointTools::map_node_to_primitive_order(*state_, grad);
88 });
89 }
90
91 void CollisionBarrierForm::solution_changed(const Eigen::VectorXd &x)
92 {
94
95 const Eigen::MatrixXd displaced_surface = collision_mesh_.vertices(utils::unflatten(get_updated_mesh_nodes(x), state_->mesh->dimension()));
96 build_collision_set(displaced_surface);
97 }
98
99 bool CollisionBarrierForm::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
100 {
101 const Eigen::MatrixXd V0 = utils::unflatten(get_updated_mesh_nodes(x0), state_->mesh->dimension());
102 const Eigen::MatrixXd V1 = utils::unflatten(get_updated_mesh_nodes(x1), state_->mesh->dimension());
103
104 // Skip CCD if the displacement is zero.
105 if ((V1 - V0).lpNorm<Eigen::Infinity>() == 0.0)
106 return true;
107
108 const ipc::TightInclusionCCD tight_inclusion_ccd(1e-6, 1e6);
109 bool is_valid = ipc::is_step_collision_free(
111 collision_mesh_.vertices(V0),
112 collision_mesh_.vertices(V1),
113 dmin_,
114 ipc::create_broad_phase(broad_phase_method_),
115 tight_inclusion_ccd);
116
117 return is_valid;
118 }
119
120 double CollisionBarrierForm::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
121 {
122 const Eigen::MatrixXd V0 = utils::unflatten(get_updated_mesh_nodes(x0), state_->mesh->dimension());
123 const Eigen::MatrixXd V1 = utils::unflatten(get_updated_mesh_nodes(x1), state_->mesh->dimension());
124
125 const ipc::TightInclusionCCD tight_inclusion_ccd(1e-6, 1e6);
126 double max_step = ipc::compute_collision_free_stepsize(
128 collision_mesh_.vertices(V0),
129 collision_mesh_.vertices(V1),
130 dmin_,
131 ipc::create_broad_phase(broad_phase_method_),
132 tight_inclusion_ccd);
133
134 adjoint_logger().info("Objective {}: max step size is {}.", name(), max_step);
135
136 return max_step;
137 }
138
139 void CollisionBarrierForm::build_collision_set(const Eigen::MatrixXd &displaced_surface)
140 {
141 static Eigen::MatrixXd cached_displaced_surface;
142 if (cached_displaced_surface.size() == displaced_surface.size() && cached_displaced_surface == displaced_surface)
143 return;
144
145 collision_set.build(collision_mesh_, displaced_surface, dhat_, dmin_, ipc::create_broad_phase(broad_phase_method_));
146
147 cached_displaced_surface = displaced_surface;
148 }
149
150 Eigen::VectorXd CollisionBarrierForm::get_updated_mesh_nodes(const Eigen::VectorXd &x) const
151 {
152 Eigen::VectorXd X = X_init;
155 }
156
158 std::shared_ptr<const State> state,
159 const std::vector<int> &boundary_ids,
160 const double dhat,
161 const bool use_log_barrier,
162 const double dmin)
163 : CollisionBarrierForm(variable_to_simulations, std::move(state), dhat, dmin),
164 boundary_ids_(boundary_ids)
165 {
166 for (const auto &id : boundary_ids_)
167 boundary_ids_to_dof_[id] = std::set<int>();
168
170
171 if (!use_log_barrier)
172 barrier_potential_.set_barrier(std::make_shared<QuadraticBarrier>());
173 }
174
176 {
177 Eigen::MatrixXd node_positions;
178 Eigen::MatrixXi boundary_edges, boundary_triangles;
179 std::vector<Eigen::Triplet<double>> displacement_map_entries;
180 io::OutGeometryData::extract_boundary_mesh(*state_->mesh, state_->n_geom_bases, state_->geom_bases(), state_->total_local_boundary,
181 node_positions, boundary_edges, boundary_triangles, displacement_map_entries);
182
183 std::vector<bool> is_on_surface;
184 is_on_surface.resize(node_positions.rows(), false);
185
187 Eigen::MatrixXd points, uv, normals;
188 Eigen::VectorXd weights;
189 Eigen::VectorXi global_primitive_ids;
190 for (const auto &lb : state_->total_local_boundary)
191 {
192 const int e = lb.element_id();
193 bool has_samples = utils::BoundarySampler::boundary_quadrature(lb, state_->n_boundary_samples(), *state_->mesh, false, uv, points, normals, weights, global_primitive_ids);
194
195 if (!has_samples)
196 continue;
197
198 const basis::ElementBases &gbs = state_->geom_bases()[e];
199
200 vals.compute(e, state_->mesh->is_volume(), points, gbs, gbs);
201
202 for (int i = 0; i < lb.size(); ++i)
203 {
204 const int primitive_global_id = lb.global_primitive_id(i);
205 const auto nodes = gbs.local_nodes_for_primitive(primitive_global_id, *state_->mesh);
206 const int boundary_id = state_->mesh->get_boundary_id(primitive_global_id);
207
208 if (!std::count(boundary_ids_.begin(), boundary_ids_.end(), boundary_id))
209 continue;
210
211 for (long n = 0; n < nodes.size(); ++n)
212 {
213 const assembler::AssemblyValues &v = vals.basis_values[nodes(n)];
214 is_on_surface[v.global[0].index] = true;
215 assert(v.global[0].index < node_positions.rows());
216 boundary_ids_to_dof_[boundary_id].insert(v.global[0].index);
217 }
218 }
219 }
220
221 Eigen::SparseMatrix<double> displacement_map;
222 if (!displacement_map_entries.empty())
223 {
224 displacement_map.resize(node_positions.rows(), state_->n_geom_bases);
225 displacement_map.setFromTriplets(displacement_map_entries.begin(), displacement_map_entries.end());
226 }
227
228 // Fix boundary edges and boundary triangles to exclude vertices not on triangles
229 Eigen::MatrixXi boundary_edges_alt(0, 2), boundary_triangles_alt(0, 3);
230 {
231 for (int i = 0; i < boundary_edges.rows(); ++i)
232 {
233 bool on_surface = true;
234 for (int j = 0; j < boundary_edges.cols(); ++j)
235 on_surface &= is_on_surface[boundary_edges(i, j)];
236 if (on_surface)
237 {
238 boundary_edges_alt.conservativeResize(boundary_edges_alt.rows() + 1, 2);
239 boundary_edges_alt.row(boundary_edges_alt.rows() - 1) = boundary_edges.row(i);
240 }
241 }
242
243 if (state_->mesh->is_volume())
244 {
245 for (int i = 0; i < boundary_triangles.rows(); ++i)
246 {
247 bool on_surface = true;
248 for (int j = 0; j < boundary_triangles.cols(); ++j)
249 on_surface &= is_on_surface[boundary_triangles(i, j)];
250 if (on_surface)
251 {
252 boundary_triangles_alt.conservativeResize(boundary_triangles_alt.rows() + 1, 3);
253 boundary_triangles_alt.row(boundary_triangles_alt.rows() - 1) = boundary_triangles.row(i);
254 }
255 }
256 }
257 else
258 boundary_triangles_alt.resize(0, 0);
259 }
260
261 collision_mesh_ = ipc::CollisionMesh(is_on_surface,
262 std::vector<bool>(is_on_surface.size(), false),
263 node_positions,
264 boundary_edges_alt,
265 boundary_triangles_alt,
266 displacement_map);
267
268 can_collide_cache_.resize(collision_mesh_.num_vertices(), collision_mesh_.num_vertices());
269 for (int i = 0; i < can_collide_cache_.rows(); ++i)
270 {
271 int dof_idx_i = collision_mesh_.to_full_vertex_id(i);
272 if (!is_on_surface[dof_idx_i])
273 continue;
274 for (int j = 0; j < can_collide_cache_.cols(); ++j)
275 {
276 int dof_idx_j = collision_mesh_.to_full_vertex_id(j);
277 if (!is_on_surface[dof_idx_j])
278 continue;
279
280 bool collision_allowed = true;
281 for (const auto &id : boundary_ids_)
282 if (boundary_ids_to_dof_[id].count(dof_idx_i) && boundary_ids_to_dof_[id].count(dof_idx_j))
283 collision_allowed = false;
284 can_collide_cache_(i, j) = collision_allowed;
285 }
286 }
287
288 collision_mesh_.can_collide = [&](size_t vi, size_t vj) {
289 return (bool)can_collide_cache_(vi, vj);
290 };
291
292 collision_mesh_.init_area_jacobians();
293 }
294
295 DeformedCollisionBarrierForm::DeformedCollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, std::shared_ptr<const State> state, std::shared_ptr<const DiffCache> diff_cache, const double dhat)
296 : AdjointForm(variable_to_simulation), state_(std::move(state)), diff_cache_(std::move(diff_cache)), dhat_(dhat), barrier_potential_(dhat)
297 {
298 if (state_->n_bases != state_->n_geom_bases)
299 log_and_throw_adjoint_error("[{}] Should use linear FE basis!", name());
300
302 *state_->mesh, state_->n_geom_bases, state_->geom_bases(), state_->geom_bases(),
303 state_->total_local_boundary, state_->obstacle, state_->args,
304 [this](const std::string &p) { return this->state_->resolve_input_path(p); },
305 state_->in_node_to_node, collision_mesh_);
306
307 Eigen::MatrixXd V;
308 state_->get_vertices(V);
310
311 broad_phase_method_ = ipc::BroadPhaseMethod::HASH_GRID;
312 }
313
314 double DeformedCollisionBarrierForm::value_unweighted(const Eigen::VectorXd &x) const
315 {
316 const Eigen::MatrixXd displaced_surface = collision_mesh_.vertices(utils::unflatten(get_updated_mesh_nodes(x), state_->mesh->dimension()));
317
318 return barrier_potential_(collision_set, collision_mesh_, displaced_surface);
319 }
320
321 void DeformedCollisionBarrierForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
322 {
324 const Eigen::MatrixXd displaced_surface = collision_mesh_.vertices(utils::unflatten(get_updated_mesh_nodes(x), state_->mesh->dimension()));
325 const Eigen::VectorXd grad = collision_mesh_.to_full_dof(barrier_potential_.gradient(collision_set, collision_mesh_, displaced_surface));
326 return AdjointTools::map_node_to_primitive_order(*state_, grad);
327 });
328 }
329
331 {
333
334 const Eigen::MatrixXd displaced_surface = collision_mesh_.vertices(utils::unflatten(get_updated_mesh_nodes(x), state_->mesh->dimension()));
335 build_collision_set(displaced_surface);
336 }
337
338 bool DeformedCollisionBarrierForm::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
339 {
340 // const Eigen::MatrixXd V0 = utils::unflatten(get_updated_mesh_nodes(x0), state_->mesh->dimension());
341 // const Eigen::MatrixXd V1 = utils::unflatten(get_updated_mesh_nodes(x1), state_->mesh->dimension());
342
343 // // Skip CCD if the displacement is zero.
344 // if ((V1 - V0).lpNorm<Eigen::Infinity>() == 0.0)
345 // return true;
346
347 // bool is_valid = ipc::is_step_collision_free(
348 // collision_mesh_,
349 // collision_mesh_.vertices(V0),
350 // collision_mesh_.vertices(V1),
351 // broad_phase_method_,
352 // 1e-6, 1e6);
353
354 return true; // is_valid;
355 }
356
357 double DeformedCollisionBarrierForm::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
358 {
359 // const Eigen::MatrixXd V0 = utils::unflatten(get_updated_mesh_nodes(x0), state_->mesh->dimension());
360 // const Eigen::MatrixXd V1 = utils::unflatten(get_updated_mesh_nodes(x1), state_->mesh->dimension());
361
362 // double max_step = ipc::compute_collision_free_stepsize(
363 // collision_mesh_,
364 // collision_mesh_.vertices(V0),
365 // collision_mesh_.vertices(V1),
366 // broad_phase_method_, 1e-6, 1e6);
367
368 return 1; // max_step;
369 }
370
371 void DeformedCollisionBarrierForm::build_collision_set(const Eigen::MatrixXd &displaced_surface)
372 {
373 static Eigen::MatrixXd cached_displaced_surface;
374 if (cached_displaced_surface.size() == displaced_surface.size() && cached_displaced_surface == displaced_surface)
375 return;
376
377 collision_set.build(collision_mesh_, displaced_surface, dhat_, 0, ipc::create_broad_phase(broad_phase_method_));
378
379 cached_displaced_surface = displaced_surface;
380 }
381
382 Eigen::VectorXd DeformedCollisionBarrierForm::get_updated_mesh_nodes(const Eigen::VectorXd &x) const
383 {
384 Eigen::VectorXd X = X_init;
387 }
388
390 const VariableToSimulationGroup &variable_to_simulations,
391 std::shared_ptr<const State> state,
392 std::shared_ptr<const DiffCache> diff_cache,
393 const json &args)
394 : StaticForm(variable_to_simulations),
395 state_(std::move(state)),
396 diff_cache_(std::move(diff_cache)),
397 params_(state_->args["contact"]["dhat"], state_->args["contact"]["alpha_t"], 0, state_->args["contact"]["alpha_n"], 0, state_->mesh->is_volume() ? 2 : 1),
398 potential_(params_)
399 {
400 auto tmp_ids = args["surface_selection"].get<std::vector<int>>();
401 boundary_ids_ = std::set(tmp_ids.begin(), tmp_ids.end());
402
404 }
405
407 {
408 // Deep copy and change the can_collide() function
409 collision_mesh_ = state_->collision_mesh;
410
411 // const int num_fe_nodes = state_->n_bases - state_->obstacle.n_vertices();
412
413 // collision_mesh_.can_collide = [this, num_fe_nodes](size_t vi, size_t vj) {
414 // return this->collision_mesh_.to_full_vertex_id(vi) >= num_fe_nodes || this->collision_mesh_.to_full_vertex_id(vj) >= num_fe_nodes;
415 // };
416
417 std::vector<int> is_obstacle(state_->n_bases);
418 for (int e = 0; e < state_->bases.size(); e++)
419 {
420 const auto &b = state_->bases[e];
421 if (state_->mesh->get_body_id(e) == 1)
422 for (const auto &bs : b.bases)
423 {
424 for (const auto &g : bs.global())
425 {
426 is_obstacle[g.index] = true;
427 }
428 }
429 }
430
431 collision_mesh_.can_collide = [this, is_obstacle](size_t vi, size_t vj) {
432 return is_obstacle[this->collision_mesh_.to_full_vertex_id(vi)] || is_obstacle[this->collision_mesh_.to_full_vertex_id(vj)];
433 };
434 }
435
436 ipc::SmoothCollisions SmoothContactForceForm::get_smooth_collision_set(const Eigen::MatrixXd &displaced_surface)
437 {
438 ipc::SmoothCollisions collisions;
439 const auto smooth_contact = dynamic_cast<const SmoothContactForm *>(state_->solve_data.contact_form.get());
440 collisions.build(collision_mesh_, displaced_surface, smooth_contact->get_params(), smooth_contact->using_adaptive_dhat(), smooth_contact->get_broad_phase());
441 return collisions;
442 }
443
444 double SmoothContactForceForm::value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const
445 {
446 assert(state_->solve_data.contact_form != nullptr);
447
448 const Eigen::MatrixXd displaced_surface = collision_mesh_.displace_vertices(utils::unflatten(diff_cache_->u(time_step), collision_mesh_.dim()));
449
450 Eigen::VectorXd forces = collision_mesh_.to_full_dof(potential_.gradient(collisions_, collision_mesh_, displaced_surface));
451
452 // return forces.squaredNorm();
453
454 Eigen::VectorXd coeff(forces.size());
455 coeff.setZero();
456 coeff(Eigen::seq(1, coeff.size(), collision_mesh_.dim())).array() = 1;
457 return (coeff.array() * forces.array()).matrix().squaredNorm() / 2;
458 }
459
460 Eigen::VectorXd SmoothContactForceForm::compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const
461 {
462 assert(state_->solve_data.contact_form != nullptr);
463
464 const Eigen::MatrixXd displaced_surface = collision_mesh_.displace_vertices(utils::unflatten(diff_cache_->u(time_step), collision_mesh_.dim()));
465
466 Eigen::VectorXd forces = potential_.gradient(collisions_, collision_mesh_, displaced_surface);
467 forces = collision_mesh_.to_full_dof(forces);
468
469 StiffnessMatrix hessian = potential_.hessian(collisions_, collision_mesh_, displaced_surface, ipc::PSDProjectionMethod::NONE);
470 hessian = collision_mesh_.to_full_dof(hessian);
471
472 Eigen::VectorXd coeff(forces.size());
473 coeff.setZero();
474 coeff(Eigen::seq(1, coeff.size(), collision_mesh_.dim())).array() = 1;
475 return weight() * (hessian * (coeff.array() * forces.array()).matrix());
476 }
477
478 void SmoothContactForceForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
479 {
480 assert(state_->solve_data.contact_form != nullptr);
481
482 const Eigen::MatrixXd displaced_surface = collision_mesh_.displace_vertices(utils::unflatten(diff_cache_->u(time_step), collision_mesh_.dim()));
483
484 Eigen::VectorXd forces = potential_.gradient(collisions_, collision_mesh_, displaced_surface);
485 forces = collision_mesh_.to_full_dof(forces);
486
487 StiffnessMatrix hessian = potential_.hessian(collisions_, collision_mesh_, displaced_surface, ipc::PSDProjectionMethod::NONE);
488 hessian = collision_mesh_.to_full_dof(hessian);
489
490 gradv = weight() * variable_to_simulations_.apply_parametrization_jacobian(ParameterType::Shape, state_.get(), x, [this, &x, &forces, &hessian]() {
491 // Eigen::VectorXd grads = 2 * hessian.transpose() * forces;
492
493 Eigen::VectorXd coeff(forces.size());
494 coeff.setZero();
495 coeff(Eigen::seq(1, coeff.size(), collision_mesh_.dim())).array() = 1;
496 Eigen::VectorXd grads = (hessian * (coeff.array() * forces.array()).matrix());
497
498 grads = diff_cache_->basis_nodes_to_gbasis_nodes() * grads;
499 return AdjointTools::map_node_to_primitive_order(*state_, grads);
500 });
501 }
502
503 void SmoothContactForceForm::solution_changed_step(const int time_step, const Eigen::VectorXd &x)
504 {
506 const Eigen::MatrixXd displaced_surface = collision_mesh_.displace_vertices(utils::unflatten(diff_cache_->u(time_step), collision_mesh_.dim()));
507 collisions_ = get_smooth_collision_set(displaced_surface);
508 }
509} // namespace polyfem::solver
int V
ElementAssemblyValues vals
Definition Assembler.cpp:22
const double weight_
int x
Storage for additional data required by differntial code.
Definition DiffCache.hpp:21
main class that contains the polyfem solver and all its state
Definition State.hpp:113
void build_collision_mesh()
extracts the boundary mesh for collision, called in build_basis
Definition State.cpp:1326
stores per local bases evaluations
std::vector< basis::Local2Global > global
stores per element basis values at given quadrature points and geometric mapping
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
Eigen::VectorXi local_nodes_for_primitive(const int local_index, const mesh::Mesh &mesh) const
static void extract_boundary_mesh(const mesh::Mesh &mesh, const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< mesh::LocalBoundary > &total_local_boundary, Eigen::MatrixXd &node_positions, Eigen::MatrixXi &boundary_edges, Eigen::MatrixXi &boundary_triangles, std::vector< Eigen::Triplet< double > > &displacement_map_entries)
extracts the boundary mesh
Definition OutData.cpp:149
virtual void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
const VariableToSimulationGroup variable_to_simulations_
ipc::BarrierPotential barrier_potential_
CollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, std::shared_ptr< const State > state, const double dhat, const double dmin=0)
std::string name() const override
void build_collision_set(const Eigen::MatrixXd &displaced_surface)
ipc::BroadPhaseMethod broad_phase_method_
std::shared_ptr< const State > state_
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine the maximum step size allowable between the current and next solution.
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
void solution_changed(const Eigen::VectorXd &x) override
Update cached fields upon a change in the solution.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
std::shared_ptr< const DiffCache > diff_cache_
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine the maximum step size allowable between the current and next solution.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
DeformedCollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, std::shared_ptr< const State > state, std::shared_ptr< const DiffCache > diff_cache, const double dhat)
const ipc::BarrierPotential barrier_potential_
void build_collision_set(const Eigen::MatrixXd &displaced_surface)
Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const
void solution_changed(const Eigen::VectorXd &x) override
Update cached fields upon a change in the solution.
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:128
std::map< int, std::set< int > > boundary_ids_to_dof_
LayerThicknessForm(const VariableToSimulationGroup &variable_to_simulations, std::shared_ptr< const State > state, const std::vector< int > &boundary_ids, const double dhat, const bool use_log_barrier=false, const double dmin=0)
SmoothContactForceForm(const VariableToSimulationGroup &variable_to_simulations, std::shared_ptr< const State > state, std::shared_ptr< const DiffCache > diff_cache, const json &args)
ipc::SmoothCollisions get_smooth_collision_set(const Eigen::MatrixXd &displaced_surface)
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override
void solution_changed_step(const int time_step, const Eigen::VectorXd &x) override
std::shared_ptr< const State > state_
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
std::shared_ptr< const DiffCache > diff_cache_
ipc::SmoothContactPotential potential_
A collection of VariableToSimulation.
void compute_state_variable(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, Eigen::VectorXd &state_variable) const
Evaluate the variable to simulations and overwrite the state_variable based on x.
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.
static bool boundary_quadrature(const mesh::LocalBoundary &local_boundary, const QuadratureOrders &order, const mesh::Mesh &mesh, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::MatrixXd &normals, Eigen::VectorXd &weights, Eigen::VectorXi &global_primitive_ids)
Eigen::VectorXd map_primitive_to_node_order(const State &state, const Eigen::VectorXd &primitives)
std::tuple< bool, int, Tree > is_valid(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u, const double threshold)
Definition Jacobian.cpp:297
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
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24