PolyFEM
Loading...
Searching...
No Matches
BarrierContactForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "ContactForm.hpp"
4
6#include <polysolve/nonlinear/PostStepData.hpp>
7
8#include <ipc/collisions/normal/normal_collisions.hpp>
9#include <ipc/potentials/barrier_potential.hpp>
10
11namespace polyfem::solver
12{
14 {
16
17 public:
18 BarrierContactForm(const ipc::CollisionMesh &collision_mesh,
19 const double dhat,
20 const double avg_mass,
21 const bool use_area_weighting,
23 const bool use_physical_barrier,
25 const bool is_time_dependent,
26 const bool enable_shape_derivatives,
27 const ipc::BroadPhaseMethod broad_phase_method,
28 const double ccd_tolerance,
29 const int ccd_max_iterations);
30
31 virtual std::string name() const override { return "barrier-contact"; }
32
33 virtual void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override;
34
38 void post_step(const polysolve::nonlinear::PostStepData &data) override;
39
41
43 bool use_area_weighting() const { return collision_set().use_area_weighting(); }
44
46 bool use_improved_max_operator() const { return collision_set().use_improved_max_approximator(); }
47
49 bool use_physical_barrier() const { return barrier_potential_.use_physical_barrier(); }
50
51 const ipc::NormalCollisions &collision_set() const { return collision_set_; }
52 const ipc::BarrierPotential &barrier_potential() const { return barrier_potential_; }
53
54 protected:
58 virtual double value_unweighted(const Eigen::VectorXd &x) const override;
59
63 Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override;
64
68 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
69
73 virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
74
75 void update_collision_set(const Eigen::MatrixXd &displaced_surface) override;
76
78 ipc::NormalCollisions collision_set_;
79
81 const ipc::BarrierPotential barrier_potential_;
82 };
83} // namespace polyfem::solver
int x
virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
bool use_improved_max_operator() const
Get use_improved_max_operator.
virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
const ipc::BarrierPotential & barrier_potential() const
bool use_convergent_formulation() const override
Get use_convergent_formulation.
bool use_area_weighting() const
Get use_area_weighting.
virtual double value_unweighted(const Eigen::VectorXd &x) const override
Compute the contact barrier potential value.
virtual std::string name() const override
Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied per element.
const ipc::NormalCollisions & collision_set() const
virtual void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override
Update the barrier stiffness based on the current elasticity energy.
void update_collision_set(const Eigen::MatrixXd &displaced_surface) override
Update the cached candidate set for the current solution.
const ipc::BarrierPotential barrier_potential_
Contact potential.
ipc::NormalCollisions collision_set_
Cached constraint set for the current solution.
void post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.
bool use_physical_barrier() const
Get use_physical_barrier.
Form representing the contact potential and forces.
bool use_adaptive_barrier_stiffness() const
Get use_adaptive_barrier_stiffness.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24