PolyFEM
Loading...
Searching...
No Matches
NormalAdhesionForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
6#include <polysolve/nonlinear/PostStepData.hpp>
7
8#include <ipc/broad_phase/broad_phase.hpp>
9#include <ipc/broad_phase/create_broad_phase.hpp>
10#include <ipc/candidates/candidates.hpp>
11#include <ipc/ccd/tight_inclusion_ccd.hpp>
12#include <ipc/collisions/normal/normal_collisions.hpp>
13#include <ipc/collision_mesh.hpp>
14#include <ipc/potentials/normal_adhesion_potential.hpp>
15
16#include <memory>
17
18namespace polyfem::solver
19{
21 class NormalAdhesionForm : public Form
22 {
24
25 public:
36 NormalAdhesionForm(const ipc::CollisionMesh &collision_mesh,
37 const double dhat_p,
38 const double dhat_a,
39 const double Y,
40 const bool is_time_dependent,
41 const bool enable_shape_derivatives,
42 const ipc::BroadPhaseMethod broad_phase_method,
43 const double ccd_tolerance,
44 const int ccd_max_iterations);
45 virtual ~NormalAdhesionForm() = default;
46
47 std::string name() const override { return "normal adhesion"; }
48
51 void init(const Eigen::VectorXd &x) override;
52
53 protected:
57 virtual double value_unweighted(const Eigen::VectorXd &x) const override;
58
62 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
63
67 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
68
72 virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
73
74 public:
78 void update_quantities(const double t, const Eigen::VectorXd &x) override;
79
83 void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
84
86 void line_search_end() override;
87
90 void solution_changed(const Eigen::VectorXd &new_x) override;
91
95 void post_step(const polysolve::nonlinear::PostStepData &data) override;
96
98 Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const;
99
101
104
105 double dhat_a() const { return dhat_a_; }
106 double dhat_p() const { return dhat_p_; }
107 double Y() const { return Y_; }
108 const ipc::NormalCollisions &collision_set() const { return collision_set_; }
109 const ipc::NormalAdhesionPotential &normal_adhesion_potential() const { return normal_adhesion_potential_; }
110
111 protected:
114 void update_collision_set(const Eigen::MatrixXd &displaced_surface);
115
117 const ipc::CollisionMesh &collision_mesh_;
118
120 const double dhat_p_;
121
123 const double dhat_a_;
124
126 const double Y_;
127
129 const double dmin_ = 0;
130
133
136
138 const ipc::BroadPhaseMethod broad_phase_method_;
139 const std::shared_ptr<ipc::BroadPhase> broad_phase_;
141 const ipc::TightInclusionCCD tight_inclusion_ccd_;
142
145
149 ipc::NormalCollisions collision_set_;
151 ipc::Candidates candidates_;
152
153 const ipc::NormalAdhesionPotential normal_adhesion_potential_;
154 };
155} // namespace polyfem::solver
int x
Form representing the contact potential and forces.
const ipc::TightInclusionCCD tight_inclusion_ccd_
Continuous collision detection specification object.
void init(const Eigen::VectorXd &x) override
Initialize the form.
std::string name() const override
const ipc::NormalAdhesionPotential & normal_adhesion_potential() const
bool use_cached_candidates_
If true, use the cached candidate set for the current solution.
const ipc::CollisionMesh & collision_mesh_
Collision mesh.
const std::shared_ptr< ipc::BroadPhase > broad_phase_
ipc::NormalCollisions collision_set_
Cached constraint set for the current solution.
virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
const double dhat_p_
Maximum adhesion strength distance.
const ipc::NormalAdhesionPotential normal_adhesion_potential_
bool save_ccd_debug_meshes
If true, output debug files.
const double Y_
Adhesion strength.
double prev_distance_
Previous minimum distance between all elements.
ipc::Candidates candidates_
Cached candidate set for the current solution.
const bool enable_shape_derivatives_
Enable shape derivatives computation.
const double dhat_a_
Adhesion activation distance.
virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
const double dmin_
Minimum distance between elements.
const bool is_time_dependent_
Is the simulation time dependent?
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time-dependent fields.
virtual double value_unweighted(const Eigen::VectorXd &x) const override
Compute the normal adhesion potential value.
void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
Initialize 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.
const ipc::NormalCollisions & collision_set() const
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::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
void update_collision_set(const Eigen::MatrixXd &displaced_surface)
Update the cached candidate set for the current solution.
void post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24