Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
TangentialAdhesionForm.cpp
Go to the documentation of this file.
3
6
7namespace polyfem::solver
8{
10 const ipc::CollisionMesh &collision_mesh,
11 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator,
12 const double epsa,
13 const double mu,
14 const ipc::BroadPhaseMethod broad_phase_method,
15 const NormalAdhesionForm &normal_adhesion_form,
16 const int n_lagging_iters)
17 : collision_mesh_(collision_mesh),
18 time_integrator_(time_integrator),
19 epsa_(epsa),
20 mu_(mu),
21 broad_phase_method_(broad_phase_method),
22 n_lagging_iters_(n_lagging_iters < 0 ? std::numeric_limits<int>::max() : n_lagging_iters),
23 normal_adhesion_form_(normal_adhesion_form),
24 tangential_adhesion_potential_(epsa)
25 {
26 assert(epsa_ > 0);
27 }
28
30 const Eigen::MatrixXd &prev_solution,
31 const Eigen::MatrixXd &solution,
32 const Eigen::MatrixXd &adjoint,
33 const ipc::TangentialCollisions &tangential_constraints_set,
34 Eigen::VectorXd &term)
35 {
36 Eigen::MatrixXd U = collision_mesh_.vertices(utils::unflatten(solution, collision_mesh_.dim()));
37 Eigen::MatrixXd U_prev = collision_mesh_.vertices(utils::unflatten(prev_solution, collision_mesh_.dim()));
38
39 // TODO: use the time integration to compute the velocity
40 const Eigen::MatrixXd velocities = (U - U_prev) / time_integrator_->dt();
41
43 tangential_constraints_set,
44 collision_mesh_, collision_mesh_.rest_positions(),
45 /*lagged_displacements=*/U_prev, velocities,
47 1,
48 ipc::TangentialPotential::DiffWRT::REST_POSITIONS);
49
50 term = collision_mesh_.to_full_dof(hess).transpose() * adjoint;
51 }
52
53 Eigen::MatrixXd TangentialAdhesionForm::compute_displaced_surface(const Eigen::VectorXd &x) const
54 {
56 }
57
58 Eigen::MatrixXd TangentialAdhesionForm::compute_surface_velocities(const Eigen::VectorXd &x) const
59 {
60 // In the case of a static problem, the velocity is the displacement
61 const Eigen::VectorXd v = time_integrator_ != nullptr ? time_integrator_->compute_velocity(x) : x;
62 return collision_mesh_.map_displacements(utils::unflatten(v, collision_mesh_.dim()));
63 }
64
66 {
67 return time_integrator_ != nullptr ? time_integrator_->dv_dx() : 1;
68 }
69
74
75 void TangentialAdhesionForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
76 {
77 const Eigen::VectorXd grad_tangential_adhesion = tangential_adhesion_potential_.gradient(
79 gradv = collision_mesh_.to_full_dof(grad_tangential_adhesion);
80 }
81
82 void TangentialAdhesionForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
83 {
84 POLYFEM_SCOPED_TIMER("tangential adhesion hessian");
85
86 ipc::PSDProjectionMethod psd_projection_method;
87
88 if (project_to_psd_) {
89 psd_projection_method = ipc::PSDProjectionMethod::CLAMP;
90 } else {
91 psd_projection_method = ipc::PSDProjectionMethod::NONE;
92 }
93
94 hessian = dv_dx() * tangential_adhesion_potential_.hessian( //
96
97 hessian = collision_mesh_.to_full_dof(hessian);
98 }
99
100 void TangentialAdhesionForm::update_lagging(const Eigen::VectorXd &x, const int iter_num)
101 {
102 const Eigen::MatrixXd displaced_surface = compute_displaced_surface(x);
103
104 ipc::NormalCollisions collision_set;
105
106 collision_set.build(
107 collision_mesh_, displaced_surface, normal_adhesion_form_.dhat_a(), /*dmin=*/0, broad_phase_method_);
108
110 collision_mesh_, displaced_surface, collision_set,
112 1., mu_);
113 }
114} // namespace polyfem::solver
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
bool project_to_psd_
If true, the form's second derivative is projected to be positive semidefinite.
Definition Form.hpp:147
Form representing the contact potential and forces.
const ipc::NormalAdhesionPotential & normal_adhesion_potential() const
Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
const ipc::BroadPhaseMethod broad_phase_method_
Broad-phase method used for distance computation and collision detection.
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)
ipc::TangentialCollisions tangential_collision_set_
Lagged tangential constraint set.
double dv_dx() const
Compute the derivative of the velocities wrt x.
Eigen::MatrixXd compute_surface_velocities(const Eigen::VectorXd &x) const
Compute the surface velocities.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
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.
TangentialAdhesionForm(const ipc::CollisionMesh &collision_mesh, const std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator, const double epsa, const double mu, const ipc::BroadPhaseMethod broad_phase_method, const NormalAdhesionForm &normal_adhesion_form, const int n_lagging_iters)
Construct a new Tangential Adhesion Form object.
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::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22