PolyFEM
Loading...
Searching...
No Matches
PeriodicContactForm.hpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <Eigen/Core>
7#include <polysolve/nonlinear/PostStepData.hpp>
8
9namespace polyfem::solver
10{
11 class PeriodicContactForceDerivative;
12
16 {
18
19 public:
23 PeriodicContactForm(const ipc::CollisionMesh &periodic_collision_mesh,
24 const Eigen::VectorXi &tiled_to_single,
25 const double dhat,
26 const double avg_mass,
27 const bool use_area_weighting,
29 const bool use_physical_barrier,
31 const bool is_time_dependent,
32 const bool enable_shape_derivatives,
33 const ipc::BroadPhaseMethod broad_phase_method,
34 const double ccd_tolerance,
35 const int ccd_max_iterations);
36
37 void init(const Eigen::VectorXd &x) override;
38
39 Eigen::VectorXd single_to_tiled(const Eigen::VectorXd &x) const;
40 Eigen::VectorXd tiled_to_single_grad(const Eigen::VectorXd &grad) const;
41
42 protected:
46 double value_unweighted(const Eigen::VectorXd &x) const override;
47
51 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
52
56 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
57
58 public:
62 void update_quantities(const double t, const Eigen::VectorXd &x) override;
63
68 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
69
73 void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
74
77 void solution_changed(const Eigen::VectorXd &new_x) override;
78
82 void post_step(const polysolve::nonlinear::PostStepData &data) override;
83
86 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
87
90 void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override;
91
92 private:
93 void update_projection() const;
94
95 const Eigen::VectorXi tiled_to_single_;
96 const int n_single_dof_;
98 };
99} // namespace polyfem::solver
int x
bool use_improved_max_operator() const
Get use_improved_max_operator.
bool use_area_weighting() const
Get use_area_weighting.
bool use_physical_barrier() const
Get use_physical_barrier.
bool use_adaptive_barrier_stiffness() const
Get use_adaptive_barrier_stiffness.
Form representing the contact potential and forces on a periodic mesh This form has a different input...
void post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.
Eigen::VectorXd single_to_tiled(const Eigen::VectorXd &x) const
void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the contact barrier potential value.
void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
Initialize variables used during the line search.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine the maximum step size allowable between the current and next solution.
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
void init(const Eigen::VectorXd &x) override
Initialize the form.
Eigen::VectorXd tiled_to_single_grad(const Eigen::VectorXd &grad) const
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time-dependent fields.
void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override
Update the barrier stiffness based on the current elasticity energy.
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24