PolyFEM
Loading...
Searching...
No Matches
TangentialAdhesionForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
7
8#include <ipc/ipc.hpp>
9#include <ipc/collision_mesh.hpp>
10#include <ipc/collisions/tangential/tangential_collisions.hpp>
11#include <ipc/potentials/tangential_adhesion_potential.hpp>
12#include <ipc/broad_phase/create_broad_phase.hpp>
13
14namespace polyfem::solver
15{
16 class NormalAdhesionForm;
17
20 {
21 public:
32 const ipc::CollisionMesh &collision_mesh,
33 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator,
34 const double epsa,
35 const double mu,
36 const ipc::BroadPhaseMethod broad_phase_method,
37 const NormalAdhesionForm &normal_adhesion_form,
38 const int n_lagging_iters);
39
40 std::string name() const override { return "tangential adhesion"; }
41
42 void force_shape_derivative(const Eigen::MatrixXd &prev_solution, const Eigen::MatrixXd &solution, const Eigen::MatrixXd &adjoint, const ipc::TangentialCollisions &tangential_constraints_set, Eigen::VectorXd &term);
43
44 protected:
48 double value_unweighted(const Eigen::VectorXd &x) const override;
49
53 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
54
58 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
59
60 public:
63 void init_lagging(const Eigen::VectorXd &x) override { update_lagging(x, 0); }
64
67 void update_lagging(const Eigen::VectorXd &x, const int iter_num) override;
68
71 void update_lagging(const Eigen::VectorXd &x) { update_lagging(x, -1); };
72
74 int max_lagging_iterations() const override { return n_lagging_iters_; }
75
78 bool uses_lagging() const override { return true; }
79
81 Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const;
83 Eigen::MatrixXd compute_surface_velocities(const Eigen::VectorXd &x) const;
85 double dv_dx() const;
86
87 double mu() const { return mu_; }
88 double epsa() const { return epsa_; }
89 const ipc::TangentialCollisions &tangential_collision_set() const { return tangential_collision_set_; }
90 const ipc::TangentialAdhesionPotential &tangential_adhesion_potential() const { return tangential_adhesion_potential_; }
91
92 private:
94 const ipc::CollisionMesh &collision_mesh_;
95
97 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator_;
98
99 const double epsa_;
100 const double mu_;
101 const ipc::BroadPhaseMethod broad_phase_method_;
102 const std::shared_ptr<ipc::BroadPhase> broad_phase_;
104
105 ipc::TangentialCollisions tangential_collision_set_;
106
108
109 const ipc::TangentialAdhesionPotential tangential_adhesion_potential_;
110 };
111} // namespace polyfem::solver
int x
Form representing the contact potential and forces.
Form of the lagged tangential adhesion disapative potential and forces.
const ipc::BroadPhaseMethod broad_phase_method_
Broad-phase method used for distance computation and collision detection.
bool uses_lagging() const override
Does this form require lagging?
const ipc::TangentialAdhesionPotential & tangential_adhesion_potential() const
void force_shape_derivative(const Eigen::MatrixXd &prev_solution, const Eigen::MatrixXd &solution, const Eigen::MatrixXd &adjoint, const ipc::TangentialCollisions &tangential_constraints_set, Eigen::VectorXd &term)
const ipc::TangentialCollisions & tangential_collision_set() const
ipc::TangentialCollisions tangential_collision_set_
Lagged tangential constraint set.
double dv_dx() const
Compute the derivative of the velocities wrt x.
const int n_lagging_iters_
Number of lagging iterations.
Eigen::MatrixXd compute_surface_velocities(const Eigen::VectorXd &x) const
Compute the surface velocities.
void init_lagging(const Eigen::VectorXd &x) override
Initialize lagged fields.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
const std::shared_ptr< ipc::BroadPhase > broad_phase_
int max_lagging_iterations() const override
Get the maximum number of lagging iteration allowable.
const ipc::CollisionMesh & collision_mesh_
Reference to the collision mesh.
void update_lagging(const Eigen::VectorXd &x, const int iter_num) override
Update lagged fields.
const double epsa_
Smoothing factor for turning on/off tangential adhesion.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
const NormalAdhesionForm & normal_adhesion_form_
necessary to have the barrier stiffnes, maybe clean me
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
void update_lagging(const Eigen::VectorXd &x)
Update lagged fields.
Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
const double mu_
Global coefficient of tangential adhesion.
const std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator_
Pointer to the time integrator.
const ipc::TangentialAdhesionPotential tangential_adhesion_potential_
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24