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/potential.hpp>
12
13// map BroadPhaseMethod values to JSON as strings
14namespace ipc
15{
16 // map ipc::BroadPhaseMethod values to JSON as strings
17 NLOHMANN_JSON_SERIALIZE_ENUM(
18 ipc::BroadPhaseMethod,
19 {{ipc::BroadPhaseMethod::HASH_GRID, "hash_grid"}, // also default
20 {ipc::BroadPhaseMethod::HASH_GRID, "HG"},
21 {ipc::BroadPhaseMethod::BRUTE_FORCE, "brute_force"},
22 {ipc::BroadPhaseMethod::BRUTE_FORCE, "BF"},
23 {ipc::BroadPhaseMethod::SPATIAL_HASH, "spatial_hash"},
24 {ipc::BroadPhaseMethod::SPATIAL_HASH, "SH"},
25 {ipc::BroadPhaseMethod::BVH, "bvh"},
26 {ipc::BroadPhaseMethod::BVH, "BVH"},
27 {ipc::BroadPhaseMethod::SWEEP_AND_TINIEST_QUEUE, "sweep_and_tiniest_queue"},
28 {ipc::BroadPhaseMethod::SWEEP_AND_TINIEST_QUEUE, "STQ"}})
29}
30
31namespace polyfem::solver
32{
34 class ContactForm : public Form
35 {
36 public:
46 ContactForm(const ipc::CollisionMesh &collision_mesh,
47 const double dhat,
48 const double avg_mass,
50 const bool is_time_dependent,
51 const bool enable_shape_derivatives,
52 const ipc::BroadPhaseMethod broad_phase_method,
53 const double ccd_tolerance,
54 const int ccd_max_iterations);
55 virtual ~ContactForm() = default;
56
57 virtual std::string name() const override { return "contact"; }
58
61 virtual void init(const Eigen::VectorXd &x) override;
62
63 public:
67 virtual void update_quantities(const double t, const Eigen::VectorXd &x) override;
68
73 virtual double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
74
78 virtual void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
79
81 void line_search_end() override;
82
85 virtual void solution_changed(const Eigen::VectorXd &new_x) override;
86
89 virtual bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
90
93 virtual void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) = 0;
94
96 Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const;
97
99 double barrier_stiffness() const { return barrier_stiffness_; }
105 virtual bool use_convergent_formulation() const { return false; }
106
108
109 double weight() const override { return weight_ * barrier_stiffness_; }
110
113
114 double dhat() const { return dhat_; }
115
116 std::shared_ptr<ipc::BroadPhase> get_broad_phase() const { return broad_phase_; }
117
118 protected:
121 virtual void update_collision_set(const Eigen::MatrixXd &displaced_surface) = 0;
122
123 virtual double barrier_support_size() const { return dhat_; }
124
126 const ipc::CollisionMesh &collision_mesh_;
127
129 const double dhat_;
130
132 const double dmin_ = 0;
133
140
142 const double avg_mass_;
143
146
149
151 const ipc::BroadPhaseMethod broad_phase_method_;
152 const std::shared_ptr<ipc::BroadPhase> broad_phase_;
154 const ipc::TightInclusionCCD tight_inclusion_ccd_;
155
158
162 ipc::Candidates candidates_;
163 };
164} // namespace polyfem::solver
int x
Form representing the contact potential and forces.
virtual void update_collision_set(const Eigen::MatrixXd &displaced_surface)=0
Update the cached candidate set for the current solution.
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.
virtual void init(const Eigen::VectorXd &x) override
Initialize the form.
std::shared_ptr< ipc::BroadPhase > get_broad_phase() const
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.
double barrier_stiffness() const
Get the current barrier stiffness.
void set_barrier_stiffness(const double barrier_stiffness)
Get the current barrier stiffness.
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.
double barrier_stiffness_
Barrier stiffness.
double max_barrier_stiffness_
Maximum barrier stiffness to use when using adaptive barrier stiffness.
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.
double weight() const override
Get the form's multiplicative constant weight.
const bool enable_shape_derivatives_
Enable shape derivatives computation.
virtual std::string name() const override
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.
bool use_adaptive_barrier_stiffness() const
Get use_adaptive_barrier_stiffness.
ipc::Candidates candidates_
Cached candidate set for the current solution.
virtual double barrier_support_size() const
virtual bool use_convergent_formulation() const
Get use_convergent_formulation.
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.
virtual void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy)=0
Update the barrier stiffness based on the current elasticity energy.
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