PolyFEM
Loading...
Searching...
No Matches
PeriodicContactForm.cpp
Go to the documentation of this file.
2
4#include <polysolve/nonlinear/PostStepData.hpp>
5
6#include <Eigen/Core>
7#include <Eigen/SparseCore>
8
9#include <cassert>
10#include <vector>
11
12namespace polyfem::solver
13{
14 PeriodicContactForm::PeriodicContactForm(const ipc::CollisionMesh &periodic_collision_mesh,
15 const Eigen::VectorXi &tiled_to_single,
16 const double dhat,
17 const double avg_mass,
18 const bool use_area_weighting,
19 const bool use_improved_max_operator,
20 const bool use_physical_barrier,
21 const bool use_adaptive_barrier_stiffness,
22 const bool is_time_dependent,
23 const bool enable_shape_derivatives,
24 const ipc::BroadPhaseMethod broad_phase_method,
25 const double ccd_tolerance,
26 const int ccd_max_iterations) : BarrierContactForm(periodic_collision_mesh, dhat, avg_mass, use_area_weighting, use_improved_max_operator, use_physical_barrier, use_adaptive_barrier_stiffness, is_time_dependent, enable_shape_derivatives, broad_phase_method, ccd_tolerance, ccd_max_iterations), tiled_to_single_(tiled_to_single), n_single_dof_(tiled_to_single_.maxCoeff() + 1)
27 {
28 assert(tiled_to_single_.size() == collision_mesh_.full_num_vertices());
29
31
32 // const Eigen::MatrixXd displaced = collision_mesh_.displace_vertices(
33 // Eigen::MatrixXd::Zero(collision_mesh_.full_num_vertices(), collision_mesh_.dim()));
34
35 // io::OBJWriter::write(
36 // "tiled-rest.obj", displaced,
37 // collision_mesh_.edges(), collision_mesh_.faces());
38 }
39
41 {
42 const int dim = collision_mesh_.dim();
43 const auto &boundary_vertices = collision_mesh_.rest_positions();
44
45 std::vector<Eigen::Triplet<double>> entries;
46 for (int i = 0; i < collision_mesh_.num_vertices(); i++)
47 {
48 const int i_full = collision_mesh_.to_full_vertex_id(i);
49 for (int d = 0; d < dim; d++)
50 entries.emplace_back(tiled_to_single_(i_full) * dim + d, i_full * dim + d, 1);
51
52 for (int p = 0; p < dim; p++)
53 for (int q = 0; q < dim; q++)
54 entries.emplace_back(n_single_dof_ * dim + p * dim + q, i_full * dim + p, boundary_vertices(i, q));
55 }
56
57 proj.resize(n_single_dof_ * dim + dim * dim, tiled_to_single_.size() * dim);
58 proj.setZero();
59 proj.setFromTriplets(entries.begin(), entries.end());
60 }
61
62 Eigen::VectorXd PeriodicContactForm::single_to_tiled(const Eigen::VectorXd &x) const
63 {
64 const int dim = collision_mesh_.dim();
65 const auto &boundary_vertices = collision_mesh_.rest_positions();
66 assert(x.size() == n_single_dof_ * dim + dim * dim);
67
68 Eigen::VectorXd tiled_x;
69 tiled_x.setZero(tiled_to_single_.size() * dim);
70 for (int i = 0; i < collision_mesh_.num_vertices(); i++)
71 {
72 const int i_full = collision_mesh_.to_full_vertex_id(i);
73 for (int d = 0; d < dim; d++)
74 tiled_x(i_full * dim + d) += x(tiled_to_single_(i_full) * dim + d);
75
76 for (int p = 0; p < dim; p++)
77 for (int q = 0; q < dim; q++)
78 tiled_x(i_full * dim + p) += x(n_single_dof_ * dim + p * dim + q) * boundary_vertices(i, q);
79 }
80
81 return tiled_x;
82 }
83 Eigen::VectorXd PeriodicContactForm::tiled_to_single_grad(const Eigen::VectorXd &grad) const
84 {
85 const int dim = collision_mesh_.dim();
86 const auto &boundary_vertices = collision_mesh_.rest_positions();
87 assert(grad.size() == tiled_to_single_.size() * dim);
88
89 Eigen::VectorXd reduced_grad;
90 reduced_grad.setZero(n_single_dof_ * dim + dim * dim);
91 for (int i = 0; i < collision_mesh_.num_vertices(); i++)
92 {
93 const int i_full = collision_mesh_.to_full_vertex_id(i);
94 for (int d = 0; d < dim; d++)
95 reduced_grad(tiled_to_single_(i_full) * dim + d) += grad(i_full * dim + d);
96
97 for (int p = 0; p < dim; p++)
98 for (int q = 0; q < dim; q++)
99 reduced_grad(n_single_dof_ * dim + p * dim + q) += boundary_vertices(i, q) * grad(i_full * dim + p);
100 }
101
102 return reduced_grad;
103 }
104
105 void PeriodicContactForm::init(const Eigen::VectorXd &x)
106 {
108 }
109
110 double PeriodicContactForm::value_unweighted(const Eigen::VectorXd &x) const
111 {
113 }
114
115 void PeriodicContactForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
116 {
118 gradv = tiled_to_single_grad(gradv);
119 }
120
121 void PeriodicContactForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
122 {
123 StiffnessMatrix hessian_full;
125
127 hessian = proj * hessian_full * proj.transpose();
128
129 // const Eigen::MatrixXd displaced = collision_mesh_.displace_vertices(utils::unflatten(single_to_tiled(x), collision_mesh_.dim()));
130
131 // io::OBJWriter::write(
132 // "tiled.obj", displaced,
133 // collision_mesh_.edges(), collision_mesh_.faces());
134 }
135
136 void PeriodicContactForm::update_quantities(const double t, const Eigen::VectorXd &x)
137 {
139 }
140
141 double PeriodicContactForm::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
142 {
144 }
145
146 void PeriodicContactForm::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
147 {
149 }
150
151 void PeriodicContactForm::solution_changed(const Eigen::VectorXd &new_x)
152 {
154 }
155
156 void PeriodicContactForm::post_step(const polysolve::nonlinear::PostStepData &data)
157 {
159 polysolve::nonlinear::PostStepData(
160 data.iter_num, data.solver_info, single_to_tiled(data.x), single_to_tiled(data.grad)));
161 }
162
163 bool PeriodicContactForm::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const
164 {
166 }
167
168 void PeriodicContactForm::update_barrier_stiffness(const Eigen::VectorXd &x, const Eigen::MatrixXd &grad_energy)
169 {
171 }
172} // namespace polyfem::solver
std::vector< Eigen::Triplet< double > > entries
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.
virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
virtual double value_unweighted(const Eigen::VectorXd &x) const override
Compute the contact barrier potential value.
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 post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.
virtual void init(const Eigen::VectorXd &x) override
Initialize the form.
virtual 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.
virtual void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
const ipc::CollisionMesh & collision_mesh_
Collision mesh.
virtual void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
Initialize variables used during the line search.
virtual void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time-dependent fields.
virtual bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
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.
PeriodicContactForm(const ipc::CollisionMesh &periodic_collision_mesh, const Eigen::VectorXi &tiled_to_single, const double dhat, const double avg_mass, const bool use_area_weighting, const bool use_improved_max_operator, const bool use_physical_barrier, const bool use_adaptive_barrier_stiffness, const bool is_time_dependent, const bool enable_shape_derivatives, const ipc::BroadPhaseMethod broad_phase_method, const double ccd_tolerance, const int ccd_max_iterations)
Construct a new Contact Form object.
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