PolyFEM
Loading...
Searching...
No Matches
ContactForm.cpp
Go to the documentation of this file.
1#include "ContactForm.hpp"
2
10
12
13#include <ipc/barrier/adaptive_stiffness.hpp>
14#include <ipc/utils/world_bbox_diagonal_length.hpp>
15
16#include <igl/writePLY.h>
17
18namespace polyfem::solver
19{
20 ContactForm::ContactForm(const ipc::CollisionMesh &collision_mesh,
21 const double dhat,
22 const double avg_mass,
23 const bool use_adaptive_barrier_stiffness,
24 const bool is_time_dependent,
25 const bool enable_shape_derivatives,
26 const ipc::BroadPhaseMethod broad_phase_method,
27 const double ccd_tolerance,
28 const int ccd_max_iterations)
29 : collision_mesh_(collision_mesh),
30 dhat_(dhat),
31 use_adaptive_barrier_stiffness_(use_adaptive_barrier_stiffness),
32 barrier_stiffness_(1.0),
33 avg_mass_(avg_mass),
34 is_time_dependent_(is_time_dependent),
35 enable_shape_derivatives_(enable_shape_derivatives),
36 broad_phase_method_(broad_phase_method),
37 broad_phase_(ipc::create_broad_phase(broad_phase_method)),
38 tight_inclusion_ccd_(ccd_tolerance, ccd_max_iterations)
39 {
40 assert(dhat_ > 0);
41 assert(ccd_tolerance > 0);
42
43 prev_distance_ = -1;
44 }
45
46 void ContactForm::init(const Eigen::VectorXd &x)
47 {
49 }
50
51 void ContactForm::update_quantities(const double t, const Eigen::VectorXd &x)
52 {
54 }
55
56 Eigen::MatrixXd ContactForm::compute_displaced_surface(const Eigen::VectorXd &x) const
57 {
58 return collision_mesh_.displace_vertices(utils::unflatten(x, collision_mesh_.dim()));
59 }
60
61 void ContactForm::solution_changed(const Eigen::VectorXd &new_x)
62 {
64 }
65
66 double ContactForm::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
67 {
68 // Extract surface only
69 const Eigen::MatrixXd V0 = compute_displaced_surface(x0);
70 const Eigen::MatrixXd V1 = compute_displaced_surface(x1);
71
73 {
74 const Eigen::MatrixXi E = collision_mesh_.dim() == 2 ? Eigen::MatrixXi() : collision_mesh_.edges();
75 const Eigen::MatrixXi &F = collision_mesh_.faces();
76 igl::writePLY(resolve_output_path("debug_ccd_0.ply"), V0, F, E);
77 igl::writePLY(resolve_output_path("debug_ccd_1.ply"), V1, F, E);
78 }
79
80 double max_step;
81 if (use_cached_candidates_ && broad_phase_method_ != ipc::BroadPhaseMethod::SWEEP_AND_TINIEST_QUEUE)
82 max_step = candidates_.compute_collision_free_stepsize(
84 else
85 max_step = ipc::compute_collision_free_stepsize(
87
88 if (save_ccd_debug_meshes && ipc::has_intersections(collision_mesh_, (V1 - V0) * max_step + V0, broad_phase_.get()))
89 {
90 log_and_throw_error("Taking max_step results in intersections (max_step={})", max_step);
91 }
92
93#ifndef NDEBUG
94 // This will check for static intersections as a failsafe. Not needed if we use our conservative CCD.
95 Eigen::MatrixXd V_toi = (V1 - V0) * max_step + V0;
96
97 while (ipc::has_intersections(collision_mesh_, V_toi, broad_phase_.get()))
98 {
99 logger().error("Taking max_step results in intersections (max_step={:g})", max_step);
100 max_step /= 2.0;
101
102 const double Linf = (V_toi - V0).lpNorm<Eigen::Infinity>();
103 if (max_step <= 0 || Linf == 0)
104 log_and_throw_error("Unable to find an intersection free step size (max_step={:g} L∞={:g})", max_step, Linf);
105
106 V_toi = (V1 - V0) * max_step + V0;
107 }
108#endif
109
110 return max_step;
111 }
112
113 void ContactForm::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
114 {
115 candidates_.build(
119 /*inflation_radius=*/barrier_support_size() / 2,
120 broad_phase_.get());
121
123 }
124
126 {
127 candidates_.clear();
129 }
130
131 bool ContactForm::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
132 {
133 const auto displaced0 = compute_displaced_surface(x0);
134 const auto displaced1 = compute_displaced_surface(x1);
135
136 // Skip CCD if the displacement is zero.
137 if ((displaced1 - displaced0).lpNorm<Eigen::Infinity>() == 0.0)
138 {
139 // Assumes initially intersection-free
140 return true;
141 }
142
143 bool is_valid;
145 is_valid = candidates_.is_step_collision_free(
146 collision_mesh_, displaced0, displaced1, dmin_,
148 else
149 is_valid = ipc::is_step_collision_free(
150 collision_mesh_, displaced0, displaced1, dmin_, broad_phase_.get(),
152
153 return is_valid;
154 }
155} // namespace polyfem::solver
int x
virtual void update_collision_set(const Eigen::MatrixXd &displaced_surface)=0
Update the cached candidate set for the current solution.
const double dmin_
Minimum distance between elements.
virtual void init(const Eigen::VectorXd &x) override
Initialize the form.
bool use_cached_candidates_
If true, use the cached candidate set for the current solution.
Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
virtual 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.
virtual void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
double prev_distance_
Previous minimum distance between all elements.
const ipc::BroadPhaseMethod broad_phase_method_
Broad phase method to use for distance and CCD evaluations.
void line_search_end() override
Clear variables used during the line search.
const ipc::CollisionMesh & collision_mesh_
Collision mesh.
const ipc::TightInclusionCCD tight_inclusion_ccd_
Continuous collision detection specification object.
virtual void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
Initialize variables used during the line search.
ipc::Candidates candidates_
Cached candidate set for the current solution.
virtual double barrier_support_size() const
ContactForm(const ipc::CollisionMesh &collision_mesh, const double dhat, const double avg_mass, const bool use_adaptive_barrier_stiffness, const bool is_time_dependent, const bool enable_shape_derivatives, const ipc::BroadPhaseMethod broad_phase_method, const double ccd_tolerance, const int ccd_max_iterations)
Construct a new Contact Form object.
virtual void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time-dependent fields.
bool save_ccd_debug_meshes
If true, output debug files.
const double dhat_
Barrier activation distance.
const std::shared_ptr< ipc::BroadPhase > broad_phase_
virtual bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
std::string resolve_output_path(const std::string &path) const
Definition Form.hpp:155
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.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73