Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SmoothingForms.cpp
Go to the documentation of this file.
1#include "SmoothingForms.hpp"
2#include <polyfem/State.hpp>
4
5namespace polyfem::solver
6{
8 const VariableToSimulationGroup &variable_to_simulations,
9 const State &state,
10 const bool scale_invariant,
11 const int power,
12 const std::vector<int> &surface_selections,
13 const std::vector<int> &active_dims) :
14 AdjointForm(variable_to_simulations), state_(state), scale_invariant_(scale_invariant), power_(power), active_dims_(active_dims)
15 {
16 const auto &mesh = *(state_.mesh);
17 const int dim = mesh.dimension();
18 const int n_verts = mesh.n_vertices();
19 assert(mesh.is_simplicial());
20
21 surface_ids_ = std::set(surface_selections.begin(), surface_selections.end());
22
23 // collect active nodes
24 std::vector<bool> active_mask;
25 active_mask.assign(n_verts, false);
26 std::vector<Eigen::Triplet<bool>> T_adj;
27
28 for (int b = 0; b < mesh.n_boundary_elements(); b++)
29 {
30 const int boundary_id = mesh.get_boundary_id(b);
31 if (!surface_ids_.empty() && surface_ids_.find(boundary_id) == surface_ids_.end())
32 continue;
33
34 for (int lv = 0; lv < dim; lv++)
35 {
36 active_mask[mesh.boundary_element_vertex(b, lv)] = true;
37 }
38
39 for (int lv1 = 0; lv1 < dim; lv1++)
40 for (int lv2 = 0; lv2 < lv1; lv2++)
41 {
42 const int v1 = mesh.boundary_element_vertex(b, lv1);
43 const int v2 = mesh.boundary_element_vertex(b, lv2);
44 T_adj.emplace_back(v2, v1, true);
45 T_adj.emplace_back(v1, v2, true);
46 }
47 }
48
49 adj.setZero();
50 adj.resize(n_verts, n_verts);
51 adj.setFromTriplets(T_adj.begin(), T_adj.end());
52
53 std::vector<int> degrees(n_verts, 0);
54 for (int k = 0; k < adj.outerSize(); ++k)
55 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, k); it; ++it)
56 degrees[k]++;
57
58 L.setZero();
59 L.resize(n_verts, n_verts);
61 {
62 std::vector<Eigen::Triplet<double>> T_L;
63 for (int k = 0; k < adj.outerSize(); ++k)
64 {
65 if (!active_mask[k])
66 continue;
67 T_L.emplace_back(k, k, 1);
68 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, k); it; ++it)
69 {
70 assert(it.row() == k);
71 T_L.emplace_back(it.row(), it.col(), -1. / degrees[k]);
72 }
73 }
74 L.setFromTriplets(T_L.begin(), T_L.end());
75 L.prune([](int i, int j, double val) { return abs(val) > 1e-12; });
76 }
77 }
78
79 double BoundarySmoothingForm::value_unweighted(const Eigen::VectorXd &x) const
80 {
81 const auto &mesh = *(state_.mesh);
82 const int dim = mesh.dimension();
83 const int n_verts = mesh.n_vertices();
84
85 double val = 0;
87 {
88 for (int b = 0; b < adj.rows(); b++)
89 {
91 s.setZero(dim);
92 double sum_norm = 0;
93 int valence = 0;
94 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
95 {
96 assert(it.col() != b);
97 auto x = mesh.point(b) - mesh.point(it.col());
98 s += x;
99 sum_norm += x.norm();
100 valence += 1;
101 }
102 if (valence)
103 {
104 s = s / sum_norm;
105 val += pow(s.norm(), power_);
106 }
107 }
108 }
109 else
110 {
111 Eigen::MatrixXd V;
113
114 val = (L * V(Eigen::all, active_dims_)).squaredNorm();
115 }
116
117 return val;
118 }
119
120 void BoundarySmoothingForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
121 {
122 const auto &mesh = *(state_.mesh);
123 const int dim = mesh.dimension();
124 const int n_verts = mesh.n_vertices();
125
126 Eigen::VectorXd grad;
128 {
129 grad.setZero(n_verts * dim);
130 for (int b = 0; b < adj.rows(); b++)
131 {
133 s.setZero(dim);
134 double sum_norm = 0;
135 auto sum_normalized = s;
136 int valence = 0;
137 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
138 {
139 assert(it.col() != b);
140 auto x = mesh.point(b) - mesh.point(it.col());
141 s += x;
142 sum_norm += x.norm();
143 sum_normalized += x.normalized();
144 valence += 1;
145 }
146 if (valence)
147 {
148 s = s / sum_norm;
149 const double coeff = power_ * pow(s.norm(), power_ - 2.) / sum_norm;
150
151 grad.segment(b * dim, dim) += (s * valence - s.squaredNorm() * sum_normalized) * coeff;
152 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
153 grad.segment(it.col() * dim, dim) -= (s + s.squaredNorm() * (mesh.point(it.col()) - mesh.point(b)).normalized()) * coeff;
154 }
155 }
156 }
157 else
158 {
159 Eigen::MatrixXd V;
161
162 Eigen::MatrixXd grad_mat = 2 * (L.transpose() * (L * V));
163 for (int d = 0; d < dim; d++)
164 if (std::find(active_dims_.begin(), active_dims_.end(), d) == active_dims_.end())
165 grad_mat.col(d).setZero();
166 grad = utils::flatten(grad_mat);
167 }
168
170 return grad;
171 });
172 }
173} // namespace polyfem::solver
int V
double val
Definition Assembler.cpp:86
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
void get_vertices(Eigen::MatrixXd &vertices) const
Definition StateDiff.cpp:72
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:471
const VariableToSimulationGroup variable_to_simulations_
Eigen::SparseMatrix< bool, Eigen::RowMajor > adj
BoundarySmoothingForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const bool scale_invariant, const int power, const std::vector< int > &surface_selections, const std::vector< int > &active_dims)
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Eigen::SparseMatrix< double, Eigen::RowMajor > L
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:128
A collection of VariableToSimulation.
virtual Eigen::VectorXd apply_parametrization_jacobian(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, const std::function< Eigen::VectorXd()> &grad) const
Maps the partial gradient wrt.
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13