PolyFEM
Loading...
Searching...
No Matches
NormalAdhesionForm.cpp
Go to the documentation of this file.
2
7
8#include <ipc/broad_phase/create_broad_phase.hpp>
9#include <ipc/potentials/potential.hpp>
10
11#include <algorithm>
12#include <cassert>
13#include <cmath>
14
15namespace polyfem::solver
16{
17 NormalAdhesionForm::NormalAdhesionForm(const ipc::CollisionMesh &collision_mesh,
18 const double dhat_p,
19 const double dhat_a,
20 const double Y,
21 const bool is_time_dependent,
22 const bool enable_shape_derivatives,
23 const ipc::BroadPhaseMethod broad_phase_method,
24 const double ccd_tolerance,
25 const int ccd_max_iterations)
26 : collision_mesh_(collision_mesh),
27 dhat_p_(dhat_p),
28 dhat_a_(dhat_a),
29 Y_(Y),
30 is_time_dependent_(is_time_dependent),
31 enable_shape_derivatives_(enable_shape_derivatives),
32 broad_phase_method_(broad_phase_method),
33 broad_phase_(ipc::create_broad_phase(broad_phase_method)),
34 tight_inclusion_ccd_(ccd_tolerance, ccd_max_iterations),
35 normal_adhesion_potential_(dhat_p, dhat_a, Y, 1)
36 {
37 assert(dhat_p > 0);
38 assert(dhat_a > dhat_p);
39 assert(ccd_tolerance > 0);
40
41 prev_distance_ = -1;
42 collision_set_.set_enable_shape_derivatives(enable_shape_derivatives);
43 }
44
45 void NormalAdhesionForm::init(const Eigen::VectorXd &x)
46 {
48 }
49
50 void NormalAdhesionForm::update_quantities(const double t, const Eigen::VectorXd &x)
51 {
53 }
54
55 Eigen::MatrixXd NormalAdhesionForm::compute_displaced_surface(const Eigen::VectorXd &x) const
56 {
57 return collision_mesh_.displace_vertices(utils::unflatten(x, collision_mesh_.dim()));
58 }
59
60 void NormalAdhesionForm::update_collision_set(const Eigen::MatrixXd &displaced_surface)
61 {
62 // Store the previous value used to compute the constraint set to avoid duplicate computation.
63 static Eigen::MatrixXd cached_displaced_surface;
64 if (cached_displaced_surface.size() == displaced_surface.size() && cached_displaced_surface == displaced_surface)
65 return;
66
68 collision_set_.build(
69 candidates_, collision_mesh_, displaced_surface, dhat_a_);
70 else
71 collision_set_.build(
72 collision_mesh_, displaced_surface, dhat_a_, dmin_, broad_phase_);
73 cached_displaced_surface = displaced_surface;
74 }
75
80
81 Eigen::VectorXd NormalAdhesionForm::value_per_element_unweighted(const Eigen::VectorXd &x) const
82 {
83 const Eigen::MatrixXd V = compute_displaced_surface(x);
84 assert(V.rows() == collision_mesh_.num_vertices());
85
86 const size_t num_vertices = collision_mesh_.num_vertices();
87
88 if (collision_set_.empty())
89 {
90 return Eigen::VectorXd::Zero(collision_mesh_.full_num_vertices());
91 }
92
93 const Eigen::MatrixXi &E = collision_mesh_.edges();
94 const Eigen::MatrixXi &F = collision_mesh_.faces();
95
96 auto storage = utils::create_thread_storage<Eigen::VectorXd>(Eigen::VectorXd::Zero(num_vertices));
97
98 utils::maybe_parallel_for(collision_set_.size(), [&](int start, int end, int thread_id) {
99 Eigen::VectorXd &local_storage = utils::get_local_thread_storage(storage, thread_id);
100
101 for (size_t i = start; i < end; i++)
102 {
103 // Quadrature weight is premultiplied by compute_potential
104 const double potential = normal_adhesion_potential_(collision_set_[i], collision_set_[i].dof(V, E, F));
105
106 const int n_v = collision_set_[i].num_vertices();
107 const auto vis = collision_set_[i].vertex_ids(E, F);
108 for (int j = 0; j < n_v; j++)
109 {
110 assert(0 <= vis[j] && vis[j] < num_vertices);
111 local_storage[vis[j]] += potential / n_v;
112 }
113 }
114 });
115
116 Eigen::VectorXd out = Eigen::VectorXd::Zero(num_vertices);
117 for (const auto &local_potential : storage)
118 {
119 out += local_potential;
120 }
121
122 Eigen::VectorXd out_full = Eigen::VectorXd::Zero(collision_mesh_.full_num_vertices());
123 for (int i = 0; i < out.size(); i++)
124 out_full[collision_mesh_.to_full_vertex_id(i)] = out[i];
125
126 assert(std::abs(value_unweighted(x) - out_full.sum()) < std::max(1e-10 * out_full.sum(), 1e-10));
127
128 return out_full;
129 }
130
131 void NormalAdhesionForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
132 {
134 gradv = collision_mesh_.to_full_dof(gradv);
135 }
136
137 void NormalAdhesionForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
138 {
139 POLYFEM_SCOPED_TIMER("normal adhesion hessian");
140
141 ipc::PSDProjectionMethod psd_projection_method;
142
143 if (project_to_psd_)
144 {
145 psd_projection_method = ipc::PSDProjectionMethod::CLAMP;
146 }
147 else
148 {
149 psd_projection_method = ipc::PSDProjectionMethod::NONE;
150 }
151
152 hessian = normal_adhesion_potential_.hessian(collision_set_, collision_mesh_, compute_displaced_surface(x), psd_projection_method);
153 hessian = collision_mesh_.to_full_dof(hessian);
154 }
155
156 void NormalAdhesionForm::solution_changed(const Eigen::VectorXd &new_x)
157 {
159 }
160
161 void NormalAdhesionForm::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
162 {
163 candidates_.build(
167 /*inflation_radius=*/dhat_a_ / 2,
169
171 }
172
174 {
175 candidates_.clear();
177 }
178
179 void NormalAdhesionForm::post_step(const polysolve::nonlinear::PostStepData &data)
180 {
181 if (data.iter_num == 0)
182 return;
183
184 const Eigen::MatrixXd displaced_surface = compute_displaced_surface(data.x);
185
186 const double curr_distance = collision_set_.compute_minimum_distance(collision_mesh_, displaced_surface);
187
188 prev_distance_ = curr_distance;
189 }
190
191} // namespace polyfem::solver
int V
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
bool project_to_psd_
If true, the form's second derivative is projected to be positive semidefinite.
Definition Form.hpp:147
void init(const Eigen::VectorXd &x) override
Initialize the form.
bool use_cached_candidates_
If true, use the cached candidate set for the current solution.
const ipc::CollisionMesh & collision_mesh_
Collision mesh.
const std::shared_ptr< ipc::BroadPhase > broad_phase_
ipc::NormalCollisions collision_set_
Cached constraint set for the current solution.
virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
const ipc::NormalAdhesionPotential normal_adhesion_potential_
double prev_distance_
Previous minimum distance between all elements.
ipc::Candidates candidates_
Cached candidate set for the current solution.
const double dhat_a_
Adhesion activation distance.
virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
const double dmin_
Minimum distance between elements.
NormalAdhesionForm(const ipc::CollisionMesh &collision_mesh, const double dhat_p, const double dhat_a, const double Y, 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 NormalAdhesion Form object.
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time-dependent fields.
virtual double value_unweighted(const Eigen::VectorXd &x) const override
Compute the normal adhesion potential value.
void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
Initialize variables used during the line search.
Eigen::VectorXd value_per_element_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied per element.
void line_search_end() override
Clear variables used during the line search.
Eigen::MatrixXd compute_displaced_surface(const Eigen::VectorXd &x) const
Compute the displaced positions of the surface nodes.
void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
void update_collision_set(const Eigen::MatrixXd &displaced_surface)
Update the cached candidate set for the current solution.
void post_step(const polysolve::nonlinear::PostStepData &data) override
Update fields after a step in the optimization.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
void maybe_parallel_for(int size, const std::function< void(int, int, int)> &partial_for)
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24