PolyFEM
Loading...
Searching...
No Matches
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 AdjointForm(variable_to_simulations), state_(state), scale_invariant_(scale_invariant), power_(power)
14 {
15 const auto &mesh = *(state_.mesh);
16 const int dim = mesh.dimension();
17 const int n_verts = mesh.n_vertices();
18 assert(mesh.is_simplicial());
19
20 surface_ids_ = std::set(surface_selections.begin(), surface_selections.end());
21
22 // collect active nodes
23 std::vector<bool> active_mask;
24 active_mask.assign(n_verts, false);
25 std::vector<Eigen::Triplet<bool>> T_adj;
26
27 for (int b = 0; b < mesh.n_boundary_elements(); b++)
28 {
29 const int boundary_id = mesh.get_boundary_id(b);
30 if (!surface_ids_.empty() && surface_ids_.find(boundary_id) == surface_ids_.end())
31 continue;
32
33 for (int lv = 0; lv < dim; lv++)
34 {
35 active_mask[mesh.boundary_element_vertex(b, lv)] = true;
36 }
37
38 for (int lv1 = 0; lv1 < dim; lv1++)
39 for (int lv2 = 0; lv2 < lv1; lv2++)
40 {
41 const int v1 = mesh.boundary_element_vertex(b, lv1);
42 const int v2 = mesh.boundary_element_vertex(b, lv2);
43 T_adj.emplace_back(v2, v1, true);
44 T_adj.emplace_back(v1, v2, true);
45 }
46 }
47
48 adj.setZero();
49 adj.resize(n_verts, n_verts);
50 adj.setFromTriplets(T_adj.begin(), T_adj.end());
51
52 std::vector<int> degrees(n_verts, 0);
53 for (int k = 0; k < adj.outerSize(); ++k)
54 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, k); it; ++it)
55 degrees[k]++;
56
57 L.setZero();
58 L.resize(n_verts, n_verts);
60 {
61 std::vector<Eigen::Triplet<double>> T_L;
62 for (int k = 0; k < adj.outerSize(); ++k)
63 {
64 if (!active_mask[k])
65 continue;
66 T_L.emplace_back(k, k, 1);
67 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, k); it; ++it)
68 {
69 assert(it.row() == k);
70 T_L.emplace_back(it.row(), it.col(), -1. / degrees[k]);
71 }
72 }
73 L.setFromTriplets(T_L.begin(), T_L.end());
74 L.prune([](int i, int j, double val) { return abs(val) > 1e-12; });
75 }
76 }
77
78 double BoundarySmoothingForm::value_unweighted(const Eigen::VectorXd &x) const
79 {
80 const auto &mesh = *(state_.mesh);
81 const int dim = mesh.dimension();
82 const int n_verts = mesh.n_vertices();
83
84 double val = 0;
86 {
87 for (int b = 0; b < adj.rows(); b++)
88 {
90 s.setZero(dim);
91 double sum_norm = 0;
92 int valence = 0;
93 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
94 {
95 assert(it.col() != b);
96 auto x = mesh.point(b) - mesh.point(it.col());
97 s += x;
98 sum_norm += x.norm();
99 valence += 1;
100 }
101 if (valence)
102 {
103 s = s / sum_norm;
104 val += pow(s.norm(), power_);
105 }
106 }
107 }
108 else
109 {
110 Eigen::MatrixXd V;
112
113 val = (L * V).eval().squaredNorm();
114 }
115
116 return val;
117 }
118
119 void BoundarySmoothingForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
120 {
121 const auto &mesh = *(state_.mesh);
122 const int dim = mesh.dimension();
123 const int n_verts = mesh.n_vertices();
124
125 Eigen::VectorXd grad;
127 {
128 grad.setZero(n_verts * dim);
129 for (int b = 0; b < adj.rows(); b++)
130 {
132 s.setZero(dim);
133 double sum_norm = 0;
134 auto sum_normalized = s;
135 int valence = 0;
136 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
137 {
138 assert(it.col() != b);
139 auto x = mesh.point(b) - mesh.point(it.col());
140 s += x;
141 sum_norm += x.norm();
142 sum_normalized += x.normalized();
143 valence += 1;
144 }
145 if (valence)
146 {
147 s = s / sum_norm;
148 const double coeff = power_ * pow(s.norm(), power_ - 2.) / sum_norm;
149
150 grad.segment(b * dim, dim) += (s * valence - s.squaredNorm() * sum_normalized) * coeff;
151 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
152 grad.segment(it.col() * dim, dim) -= (s + s.squaredNorm() * (mesh.point(it.col()) - mesh.point(b)).normalized()) * coeff;
153 }
154 }
155 }
156 else
157 {
158 Eigen::MatrixXd V;
160
161 grad = utils::flatten(2 * (L.transpose() * (L * V)));
162 }
163
165 return grad;
166 });
167 }
168} // 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:69
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
const VariableToSimulationGroup variable_to_simulations_
Eigen::SparseMatrix< bool, Eigen::RowMajor > adj
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
BoundarySmoothingForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const bool scale_invariant, const int power, const std::vector< int > &surface_selections)
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