Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
13namespace polyfem::solver
14{
15 class NormalAdhesionForm;
16
19 {
20 public:
31 const ipc::CollisionMesh &collision_mesh,
32 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator,
33 const double epsa,
34 const double mu,
35 const ipc::BroadPhaseMethod broad_phase_method,
36 const NormalAdhesionForm &normal_adhesion_form,
37 const int n_lagging_iters);
38
39 std::string name() const override { return "tangential adhesion"; }
40
41 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);
42
43 protected:
47 double value_unweighted(const Eigen::VectorXd &x) const override;
48
52 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
53
57 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
58
59 public:
62 void init_lagging(const Eigen::VectorXd &x) override { update_lagging(x, 0); }
63
66 void update_lagging(const Eigen::VectorXd &x, const int iter_num) override;
67
70 void update_lagging(const Eigen::VectorXd &x) { update_lagging(x, -1); };
71
73 int max_lagging_iterations() const override { return n_lagging_iters_; }
74
77 bool uses_lagging() const override { return true; }
78
80 Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const;
82 Eigen::MatrixXd compute_surface_velocities(const Eigen::VectorXd &x) const;
84 double dv_dx() const;
85
86 double mu() const { return mu_; }
87 double epsa() const { return epsa_; }
88 const ipc::TangentialCollisions &tangential_collision_set() const { return tangential_collision_set_; }
89 const ipc::TangentialAdhesionPotential &tangential_adhesion_potential() const { return tangential_adhesion_potential_; }
90
91 private:
93 const ipc::CollisionMesh &collision_mesh_;
94
96 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator_;
97
98 const double epsa_;
99 const double mu_;
100 const ipc::BroadPhaseMethod broad_phase_method_;
102
103 ipc::TangentialCollisions tangential_collision_set_;
104
106
107 const ipc::TangentialAdhesionPotential tangential_adhesion_potential_;
108 };
109} // 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.
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:22