23 class QuadraticBarrier :
public ipc::Barrier
26 QuadraticBarrier(
const double weight = 1) :
weight_(weight) {}
28 double operator()(
const double d,
const double dhat)
const override
33 return weight_ * (d - dhat) * (d - dhat);
35 double first_derivative(
const double d,
const double dhat)
const override
40 return 2 *
weight_ * (d - dhat);
42 double second_derivative(
const double d,
const double dhat)
const override
49 double units(
const double dhat)
const override
61 :
AdjointForm(variable_to_simulation), state_(std::move(state)), dhat_(dhat), dmin_(dmin), barrier_potential_(dhat)
66 [
this](
const std::string &p) { return this->state_->resolve_input_path(p); },
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);
105 if ((
V1 -
V0).lpNorm<Eigen::Infinity>() == 0.0)
108 const ipc::TightInclusionCCD tight_inclusion_ccd(1e-6, 1e6);
109 bool is_valid = ipc::is_step_collision_free(
115 tight_inclusion_ccd);
125 const ipc::TightInclusionCCD tight_inclusion_ccd(1e-6, 1e6);
126 double max_step = ipc::compute_collision_free_stepsize(
132 tight_inclusion_ccd);
141 static Eigen::MatrixXd cached_displaced_surface;
142 if (cached_displaced_surface.size() == displaced_surface.size() && cached_displaced_surface == displaced_surface)
147 cached_displaced_surface = displaced_surface;
152 Eigen::VectorXd X =
X_init;
158 std::shared_ptr<const State> state,
159 const std::vector<int> &boundary_ids,
161 const bool use_log_barrier,
164 boundary_ids_(boundary_ids)
171 if (!use_log_barrier)
177 Eigen::MatrixXd node_positions;
178 Eigen::MatrixXi boundary_edges, boundary_triangles;
179 std::vector<Eigen::Triplet<double>> displacement_map_entries;
181 node_positions, boundary_edges, boundary_triangles, displacement_map_entries);
183 std::vector<bool> is_on_surface;
184 is_on_surface.resize(node_positions.rows(),
false);
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)
192 const int e = lb.element_id();
200 vals.compute(e,
state_->mesh->is_volume(), points, gbs, gbs);
202 for (
int i = 0; i < lb.size(); ++i)
204 const int primitive_global_id = lb.global_primitive_id(i);
206 const int boundary_id =
state_->mesh->get_boundary_id(primitive_global_id);
211 for (
long n = 0; n < nodes.size(); ++n)
214 is_on_surface[v.
global[0].index] =
true;
215 assert(v.
global[0].index < node_positions.rows());
221 Eigen::SparseMatrix<double> displacement_map;
222 if (!displacement_map_entries.empty())
224 displacement_map.resize(node_positions.rows(),
state_->n_geom_bases);
225 displacement_map.setFromTriplets(displacement_map_entries.begin(), displacement_map_entries.end());
229 Eigen::MatrixXi boundary_edges_alt(0, 2), boundary_triangles_alt(0, 3);
231 for (
int i = 0; i < boundary_edges.rows(); ++i)
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)];
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);
243 if (
state_->mesh->is_volume())
245 for (
int i = 0; i < boundary_triangles.rows(); ++i)
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)];
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);
258 boundary_triangles_alt.resize(0, 0);
262 std::vector<bool>(is_on_surface.size(),
false),
265 boundary_triangles_alt,
272 if (!is_on_surface[dof_idx_i])
277 if (!is_on_surface[dof_idx_j])
280 bool collision_allowed =
true;
283 collision_allowed =
false;
296 :
AdjointForm(variable_to_simulation), state_(std::move(state)), diff_cache_(std::move(diff_cache)), dhat_(dhat), barrier_potential_(dhat)
304 [
this](
const std::string &p) { return this->state_->resolve_input_path(p); },
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);
373 static Eigen::MatrixXd cached_displaced_surface;
374 if (cached_displaced_surface.size() == displaced_surface.size() && cached_displaced_surface == displaced_surface)
379 cached_displaced_surface = displaced_surface;
384 Eigen::VectorXd X =
X_init;
391 std::shared_ptr<const State> state,
392 std::shared_ptr<const DiffCache> diff_cache,
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),
400 auto tmp_ids = args[
"surface_selection"].get<std::vector<int>>();
417 std::vector<int> is_obstacle(
state_->n_bases);
418 for (
int e = 0; e <
state_->bases.size(); e++)
420 const auto &b =
state_->bases[e];
421 if (
state_->mesh->get_body_id(e) == 1)
422 for (
const auto &bs : b.bases)
424 for (
const auto &g : bs.global())
426 is_obstacle[g.index] =
true;
431 collision_mesh_.can_collide = [
this, is_obstacle](
size_t vi,
size_t vj) {
438 ipc::SmoothCollisions collisions;
440 collisions.build(
collision_mesh_, displaced_surface, smooth_contact->get_params(), smooth_contact->using_adaptive_dhat(), smooth_contact->get_broad_phase());
446 assert(
state_->solve_data.contact_form !=
nullptr);
454 Eigen::VectorXd coeff(forces.size());
457 return (coeff.array() * forces.array()).matrix().squaredNorm() / 2;
462 assert(
state_->solve_data.contact_form !=
nullptr);
472 Eigen::VectorXd coeff(forces.size());
475 return weight() * (hessian * (coeff.array() * forces.array()).matrix());
480 assert(
state_->solve_data.contact_form !=
nullptr);
493 Eigen::VectorXd coeff(forces.size());
495 coeff(Eigen::seq(1, coeff.size(), collision_mesh_.dim())).array() = 1;
496 Eigen::VectorXd grads = (hessian * (coeff.array() * forces.array()).matrix());
498 grads = diff_cache_->basis_nodes_to_gbasis_nodes() * grads;
499 return AdjointTools::map_node_to_primitive_order(*state_, grads);
ElementAssemblyValues vals
Storage for additional data required by differntial code.
main class that contains the polyfem solver and all its state
void build_collision_mesh()
extracts the boundary mesh for collision, called in build_basis
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
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)
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)
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.
void log_and_throw_adjoint_error(const std::string &msg)
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix