PolyFEM
Loading...
Searching...
No Matches
PeriodicContactForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "ContactForm.hpp"
4
5namespace polyfem
6{
7 class State;
8
9 namespace solver {
10 class PeriodicMeshToMesh;
11 }
12}
13
14namespace polyfem::solver
15{
19 {
20 public:
24 PeriodicContactForm(const ipc::CollisionMesh &periodic_collision_mesh,
25 const Eigen::VectorXi &tiled_to_single,
26 const double dhat,
27 const double avg_mass,
30 const bool is_time_dependent,
31 const bool enable_shape_derivatives,
32 const ipc::BroadPhaseMethod broad_phase_method,
33 const double ccd_tolerance,
34 const int ccd_max_iterations);
35
36 void init(const Eigen::VectorXd &x) override;
37
38 void force_periodic_shape_derivative(const State& state, const PeriodicMeshToMesh &periodic_mesh_map, const Eigen::VectorXd &periodic_mesh_representation, const ipc::Collisions &contact_set, const Eigen::VectorXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term);
39
40 Eigen::VectorXd single_to_tiled(const Eigen::VectorXd &x) const;
41 Eigen::VectorXd tiled_to_single_grad(const Eigen::VectorXd &grad) const;
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:
63 void update_quantities(const double t, const Eigen::VectorXd &x) override;
64
69 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
70
74 void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
75
78 void solution_changed(const Eigen::VectorXd &new_x) override;
79
83 void post_step(const polysolve::nonlinear::PostStepData &data) override;
84
87 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
88
91 void update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy) override;
92
93 private:
94 void update_projection() const;
95
96 const Eigen::VectorXi tiled_to_single_;
97 const int n_single_dof_;
99 };
100}
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
Form representing the contact potential and forces.
bool use_adaptive_barrier_stiffness() const
Get use_adaptive_barrier_stiffness.
bool use_convergent_formulation() const
Get use_convergent_formulation.
Form representing the contact potential and forces on a periodic mesh This form has a different input...
void force_periodic_shape_derivative(const State &state, const PeriodicMeshToMesh &periodic_mesh_map, const Eigen::VectorXd &periodic_mesh_representation, const ipc::Collisions &contact_set, const Eigen::VectorXd &solution, const Eigen::VectorXd &adjoint_sol, Eigen::VectorXd &term)
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:22