PolyFEM
Loading...
Searching...
No Matches
FrictionForm.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/friction/friction_collisions.hpp>
11#include <ipc/potentials/friction_potential.hpp>
12
13namespace polyfem::solver
14{
15 class ContactForm;
16
18 class FrictionForm : public Form
19 {
20 public:
31 const ipc::CollisionMesh &collision_mesh,
32 const std::shared_ptr<time_integrator::ImplicitTimeIntegrator> time_integrator,
33 const double epsv,
34 const double mu,
35 const ipc::BroadPhaseMethod broad_phase_method,
36 const ContactForm &contact_form,
37 const int n_lagging_iters);
38
39 std::string name() const override { return "friction"; }
40
41 void force_shape_derivative(const Eigen::MatrixXd &prev_solution, const Eigen::MatrixXd &solution, const Eigen::MatrixXd &adjoint, const ipc::FrictionCollisions &friction_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 epsv() const { return epsv_; }
88 const ipc::FrictionCollisions &friction_collision_set() const { return friction_collision_set_; }
89 const ipc::FrictionPotential &friction_potential() const { return friction_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 epsv_;
99 const double mu_;
100 const ipc::BroadPhaseMethod broad_phase_method_;
102
103 ipc::FrictionCollisions friction_collision_set_;
104
106
107 const ipc::FrictionPotential friction_potential_;
108 };
109} // 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.
ipc::FrictionCollisions friction_collision_set_
Lagged friction constraint set.
int max_lagging_iterations() const override
Get the maximum number of lagging iteration allowable.
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
void force_shape_derivative(const Eigen::MatrixXd &prev_solution, const Eigen::MatrixXd &solution, const Eigen::MatrixXd &adjoint, const ipc::FrictionCollisions &friction_constraints_set, Eigen::VectorXd &term)
const ipc::FrictionCollisions & friction_collision_set() const
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:22