PolyFEM
Loading...
Searching...
No Matches
FrictionForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4#include "ContactForm.hpp"
5
8
9#include <ipc/ipc.hpp>
10#include <ipc/collision_mesh.hpp>
11#include <ipc/collisions/tangential/tangential_collisions.hpp>
12#include <ipc/potentials/friction_potential.hpp>
13#include <ipc/broad_phase/create_broad_phase.hpp>
14
15#include <memory>
16
17namespace polyfem::solver
18{
20 class FrictionForm : public Form
21 {
23
24 public:
35 const ipc::CollisionMesh &collision_mesh,
36 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator,
37 const double epsv,
38 const double mu,
39 const ipc::BroadPhaseMethod broad_phase_method,
40 const ContactForm &contact_form,
41 const int n_lagging_iters);
42
43 std::string name() const override { return "friction"; }
44
45 protected:
49 double value_unweighted(const Eigen::VectorXd &x) const override;
50
54 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
55
59 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
60
61 public:
64 void init_lagging(const Eigen::VectorXd &x) override { update_lagging(x, 0); }
65
68 void update_lagging(const Eigen::VectorXd &x, const int iter_num) override;
69
72 void update_lagging(const Eigen::VectorXd &x) { update_lagging(x, -1); };
73
75 int max_lagging_iterations() const override { return n_lagging_iters_; }
76
79 bool uses_lagging() const override { return true; }
80
82 Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const;
84 Eigen::MatrixXd compute_surface_velocities(const Eigen::VectorXd &x) const;
86 double dv_dx() const;
87
88 double mu() const { return mu_; }
89 double epsv() const { return epsv_; }
90 const ipc::TangentialCollisions &friction_collision_set() const { return friction_collision_set_; }
91 const ipc::FrictionPotential &friction_potential() const { return friction_potential_; }
92
93 private:
95 const ipc::CollisionMesh &collision_mesh_;
96
98 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator_;
99
100 const double epsv_;
101 const double mu_;
102 const ipc::BroadPhaseMethod broad_phase_method_;
104
105 ipc::TangentialCollisions friction_collision_set_;
106
108
109 const ipc::FrictionPotential friction_potential_;
110 };
111} // namespace polyfem::solver
int x
Form representing the contact potential and forces.
Form of the lagged friction disapative potential and forces.
const ipc::BroadPhaseMethod broad_phase_method_
Broad-phase method used for distance computation and collision detection.
void init_lagging(const Eigen::VectorXd &x) override
Initialize lagged fields.
const int n_lagging_iters_
Number of lagging iterations.
void update_lagging(const Eigen::VectorXd &x, const int iter_num) override
Update lagged fields.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
const ipc::CollisionMesh & collision_mesh_
Reference to the collision mesh.
std::string name() const override
const std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator_
Pointer to the time integrator.
Eigen::MatrixXd compute_surface_velocities(const Eigen::VectorXd &x) const
Compute the surface velocities.
const ipc::FrictionPotential friction_potential_
Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
const ipc::FrictionPotential & friction_potential() const
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
double dv_dx() const
Compute the derivative of the velocities wrt x.
void update_lagging(const Eigen::VectorXd &x)
Update lagged fields.
int max_lagging_iterations() const override
Get the maximum number of lagging iteration allowable.
const ipc::TangentialCollisions & friction_collision_set() const
bool uses_lagging() const override
Does this form require lagging?
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
const ContactForm & contact_form_
necessary to have the barrier stiffnes, maybe clean me
ipc::TangentialCollisions friction_collision_set_
Lagged friction constraint set.
const double epsv_
Smoothing factor between static and dynamic friction.
const double mu_
Global coefficient of friction.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24