PolyFEM
Loading...
Searching...
No Matches
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 broad_phase_(ipc::create_broad_phase(broad_phase_method)),
23 n_lagging_iters_(n_lagging_iters < 0 ? std::numeric_limits<int>::max() : n_lagging_iters),
24 normal_adhesion_form_(normal_adhesion_form),
25 tangential_adhesion_potential_(epsa)
26 {
27 assert(epsa_ > 0);
28 }
29
30 Eigen::MatrixXd TangentialAdhesionForm::compute_displaced_surface(const Eigen::VectorXd &x) const
31 {
33 }
34
35 Eigen::MatrixXd TangentialAdhesionForm::compute_surface_velocities(const Eigen::VectorXd &x) const
36 {
37 // In the case of a static problem, the velocity is the displacement
38 const Eigen::VectorXd v = time_integrator_ != nullptr ? time_integrator_->compute_velocity(x) : x;
39 return collision_mesh_.map_displacements(utils::unflatten(v, collision_mesh_.dim()));
40 }
41
43 {
44 return time_integrator_ != nullptr ? time_integrator_->dv_dx() : 1;
45 }
46
51
52 void TangentialAdhesionForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
53 {
54 const Eigen::VectorXd grad_tangential_adhesion = tangential_adhesion_potential_.gradient(
56 gradv = collision_mesh_.to_full_dof(grad_tangential_adhesion);
57 }
58
59 void TangentialAdhesionForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
60 {
61 POLYFEM_SCOPED_TIMER("tangential adhesion hessian");
62
63 ipc::PSDProjectionMethod psd_projection_method;
64
66 {
67 psd_projection_method = ipc::PSDProjectionMethod::CLAMP;
68 }
69 else
70 {
71 psd_projection_method = ipc::PSDProjectionMethod::NONE;
72 }
73
74 hessian = dv_dx() * tangential_adhesion_potential_.hessian( //
76
77 hessian = collision_mesh_.to_full_dof(hessian);
78 }
79
80 void TangentialAdhesionForm::update_lagging(const Eigen::VectorXd &x, const int iter_num)
81 {
82 const Eigen::MatrixXd displaced_surface = compute_displaced_surface(x);
83
84 ipc::NormalCollisions collision_set;
85
86 collision_set.build(
87 collision_mesh_, displaced_surface, normal_adhesion_form_.dhat_a(), /*dmin=*/0, broad_phase_);
88
90 collision_mesh_, displaced_surface, collision_set,
92 1., mu_);
93 }
94} // 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.
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 std::shared_ptr< ipc::BroadPhase > broad_phase_
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:24