Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
NormalAdhesionForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4#include "ContactForm.hpp"
5
6#include <polyfem/Common.hpp>
8
9#include <ipc/collisions/normal/normal_collisions.hpp>
10#include <ipc/collision_mesh.hpp>
11#include <ipc/broad_phase/broad_phase.hpp>
12#include <ipc/potentials/normal_adhesion_potential.hpp>
13
14
15namespace polyfem::solver
16{
18 class NormalAdhesionForm : public Form
19 {
20 public:
31 NormalAdhesionForm(const ipc::CollisionMesh &collision_mesh,
32 const double dhat_p,
33 const double dhat_a,
34 const double Y,
35 const bool is_time_dependent,
36 const bool enable_shape_derivatives,
37 const ipc::BroadPhaseMethod broad_phase_method,
38 const double ccd_tolerance,
39 const int ccd_max_iterations);
40 virtual ~NormalAdhesionForm() = default;
41
42 std::string name() const override { return "normal adhesion"; }
43
46 void init(const Eigen::VectorXd &x) override;
47
48 virtual void force_shape_derivative(const ipc::NormalCollisions &collision_set, const Eigen::MatrixXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term);
49
50 protected:
54 virtual double value_unweighted(const Eigen::VectorXd &x) const override;
55
59 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
60
64 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
65
69 virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
70
71 public:
75 void update_quantities(const double t, const Eigen::VectorXd &x) override;
76
80 void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
81
83 void line_search_end() override;
84
87 void solution_changed(const Eigen::VectorXd &new_x) override;
88
92 void post_step(const polysolve::nonlinear::PostStepData &data) override;
93
95 Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const;
96
98
101
102 double dhat_a() const { return dhat_a_; }
103 double dhat_p() const { return dhat_p_; }
104 double Y() const { return Y_; }
105 const ipc::NormalCollisions &collision_set() const { return collision_set_; }
106 const ipc::NormalAdhesionPotential &normal_adhesion_potential() const { return normal_adhesion_potential_; }
107
108 protected:
111 void update_collision_set(const Eigen::MatrixXd &displaced_surface);
112
114 const ipc::CollisionMesh &collision_mesh_;
115
117 const double dhat_p_;
118
120 const double dhat_a_;
121
123 const double Y_;
124
126 const double dmin_ = 0;
127
130
133
135 const ipc::BroadPhaseMethod broad_phase_method_;
137 const ipc::TightInclusionCCD tight_inclusion_ccd_;
138
141
145 ipc::NormalCollisions collision_set_;
147 ipc::Candidates candidates_;
148
149 const ipc::NormalAdhesionPotential normal_adhesion_potential_;
150 };
151} // 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.
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.
virtual void force_shape_derivative(const ipc::NormalCollisions &collision_set, const Eigen::MatrixXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term)
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:22