Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ContactForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
5#include <polyfem/Common.hpp>
7
8#include <ipc/collisions/normal/normal_collisions.hpp>
9#include <ipc/collision_mesh.hpp>
10#include <ipc/broad_phase/broad_phase.hpp>
11#include <ipc/potentials/barrier_potential.hpp>
12
13// map BroadPhaseMethod values to JSON as strings
14namespace ipc
15{
16 NLOHMANN_JSON_SERIALIZE_ENUM(
17 ipc::BroadPhaseMethod,
18 {{ipc::BroadPhaseMethod::HASH_GRID, "hash_grid"}, // also default
19 {ipc::BroadPhaseMethod::HASH_GRID, "HG"},
20 {ipc::BroadPhaseMethod::BRUTE_FORCE, "brute_force"},
21 {ipc::BroadPhaseMethod::BRUTE_FORCE, "BF"},
22 {ipc::BroadPhaseMethod::SPATIAL_HASH, "spatial_hash"},
23 {ipc::BroadPhaseMethod::SPATIAL_HASH, "SH"},
24 {ipc::BroadPhaseMethod::BVH, "bvh"},
25 {ipc::BroadPhaseMethod::BVH, "BVH"},
26 {ipc::BroadPhaseMethod::SWEEP_AND_PRUNE, "sweep_and_prune"},
27 {ipc::BroadPhaseMethod::SWEEP_AND_PRUNE, "SAP"},
28 {ipc::BroadPhaseMethod::SWEEP_AND_TINIEST_QUEUE, "sweep_and_tiniest_queue"},
29 {ipc::BroadPhaseMethod::SWEEP_AND_TINIEST_QUEUE, "STQ"}})
30} // namespace ipc
31
32namespace polyfem::solver
33{
35 class ContactForm : public Form
36 {
37 public:
47 ContactForm(const ipc::CollisionMesh &collision_mesh,
48 const double dhat,
49 const double avg_mass,
50 const bool use_area_weighting,
52 const bool use_physical_barrier,
54 const bool is_time_dependent,
55 const bool enable_shape_derivatives,
56 const ipc::BroadPhaseMethod broad_phase_method,
57 const double ccd_tolerance,
58 const int ccd_max_iterations);
59 virtual ~ContactForm() = default;
60
61 std::string name() const override { return "contact"; }
62
65 void init(const Eigen::VectorXd &x) override;
66
67 virtual void force_shape_derivative(const ipc::NormalCollisions &collision_set, const Eigen::MatrixXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term);
68
69 protected:
73 virtual double value_unweighted(const Eigen::VectorXd &x) const override;
74
78 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
79
83 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
84
88 virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
89
90 public:
94 void update_quantities(const double t, const Eigen::VectorXd &x) override;
95
100 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
101
105 void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
106
108 void line_search_end() override;
109
112 void solution_changed(const Eigen::VectorXd &new_x) override;
113
117 void post_step(const polysolve::nonlinear::PostStepData &data) override;
118
121 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
122
125 virtual void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy);
126
128 Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const;
129
131 double barrier_stiffness() const { return barrier_stiffness_; }
138
140 bool use_area_weighting() const {return collision_set_.use_area_weighting();}
141
143 bool use_improved_max_operator() const {return collision_set_.use_improved_max_approximator();}
144
146 bool use_physical_barrier() const {return barrier_potential_.use_physical_barrier();}
147
149
150 double weight() const override { return weight_ * barrier_stiffness_; }
151
154
155 double dhat() const { return dhat_; }
156 const ipc::NormalCollisions &collision_set() const { return collision_set_; }
157 const ipc::BarrierPotential &barrier_potential() const { return barrier_potential_; }
158
159 protected:
162 void update_collision_set(const Eigen::MatrixXd &displaced_surface);
163
165 const ipc::CollisionMesh &collision_mesh_;
166
168 const double dhat_;
169
171 const double dmin_ = 0;
172
179
181 const double avg_mass_;
182
185
188
190 const ipc::BroadPhaseMethod broad_phase_method_;
192 const ipc::TightInclusionCCD tight_inclusion_ccd_;
193
196
200 ipc::NormalCollisions collision_set_;
202 ipc::Candidates candidates_;
203
204 const ipc::BarrierPotential barrier_potential_;
205 };
206} // namespace polyfem::solver
int x
Form representing the contact potential and forces.
const double avg_mass_
Average mass of the mesh (used for adaptive barrier stiffness)
virtual ~ContactForm()=default
const bool use_adaptive_barrier_stiffness_
If true, use an adaptive barrier stiffness.
const double dmin_
Minimum distance between elements.
const ipc::BarrierPotential & barrier_potential() const
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.
std::string name() const override
bool use_physical_barrier() const
Get use_physical_barrier.
Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
double barrier_stiffness() const
Get the current barrier stiffness.
void set_barrier_stiffness(const double barrier_stiffness)
Get the current barrier stiffness.
void post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.
virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
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 force_shape_derivative(const ipc::NormalCollisions &collision_set, const Eigen::MatrixXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term)
double barrier_stiffness_
Barrier stiffness.
double max_barrier_stiffness_
Maximum barrier stiffness to use when using adaptive barrier stiffness.
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.
bool use_area_weighting() const
Get use_area_weighting.
double weight() const override
Get the form's multiplicative constant weight.
const bool enable_shape_derivatives_
Enable shape derivatives computation.
const ipc::BarrierPotential barrier_potential_
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.
Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied per element.
virtual double value_unweighted(const Eigen::VectorXd &x) const override
Compute the contact barrier potential value.
const ipc::CollisionMesh & collision_mesh_
Collision mesh.
const ipc::TightInclusionCCD tight_inclusion_ccd_
Continuous collision detection specification object.
void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
Initialize variables used during the line search.
bool use_adaptive_barrier_stiffness() const
Get use_adaptive_barrier_stiffness.
ipc::Candidates candidates_
Cached candidate set for the current solution.
bool use_improved_max_operator() const
Get use_improved_max_operator.
virtual void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy)
Update the barrier stiffness based on the current elasticity energy.
void update_collision_set(const Eigen::MatrixXd &displaced_surface)
Update the cached candidate set for the current solution.
virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
const ipc::NormalCollisions & collision_set() const
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.
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
ipc::NormalCollisions collision_set_
Cached constraint set for the current solution.
bool use_convergent_formulation() const
Get use_convergent_formulation.
const bool is_time_dependent_
Is the simulation time dependent?
double weight_
weight of the form (e.g., AL penalty weight or Δt²)
Definition Form.hpp:149
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22