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 {
9 const auto &mesh = *(state_.mesh);
10 const int dim = mesh.dimension();
11 const int n_verts = mesh.n_vertices();
12 assert(mesh.is_simplicial());
13
14 std::vector<int> tmp = {}; // args_["surface_selection"];
15 std::set<int> surface_ids = std::set(tmp.begin(), tmp.end());
16
17 // collect active nodes
18 std::vector<bool> active_mask;
19 active_mask.assign(n_verts, false);
20 std::vector<Eigen::Triplet<bool>> T_adj;
21
22 for (int b = 0; b < mesh.n_boundary_elements(); b++)
23 {
24 const int boundary_id = mesh.get_boundary_id(b);
25 if (!surface_ids.empty() && surface_ids.find(boundary_id) == surface_ids.end())
26 continue;
27
28 for (int lv = 0; lv < dim; lv++)
29 {
30 active_mask[mesh.boundary_element_vertex(b, lv)] = true;
31 }
32
33 for (int lv1 = 0; lv1 < dim; lv1++)
34 for (int lv2 = 0; lv2 < lv1; lv2++)
35 {
36 T_adj.emplace_back(mesh.boundary_element_vertex(b, lv2), mesh.boundary_element_vertex(b, lv1), true);
37 T_adj.emplace_back(mesh.boundary_element_vertex(b, lv1), mesh.boundary_element_vertex(b, lv2), true);
38 }
39 }
40
41 adj.setZero();
42 adj.resize(n_verts, n_verts);
43 adj.setFromTriplets(T_adj.begin(), T_adj.end());
44
45 std::vector<int> degrees(n_verts, 0);
46 for (int k = 0; k < adj.outerSize(); ++k)
47 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, k); it; ++it)
48 degrees[k]++;
49
50 L.setZero();
51 L.resize(n_verts, n_verts);
53 {
54 std::vector<Eigen::Triplet<double>> T_L;
55 for (int k = 0; k < adj.outerSize(); ++k)
56 {
57 if (!active_mask[k])
58 continue;
59 T_L.emplace_back(k, k, degrees[k]);
60 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, k); it; ++it)
61 {
62 assert(it.row() == k);
63 T_L.emplace_back(it.row(), it.col(), -1);
64 }
65 }
66 L.setFromTriplets(T_L.begin(), T_L.end());
67 L.prune([](int i, int j, double val) { return abs(val) > 1e-12; });
68 }
69 }
70
71 double BoundarySmoothingForm::value_unweighted(const Eigen::VectorXd &x) const
72 {
73 const auto &mesh = *(state_.mesh);
74 const int dim = mesh.dimension();
75 const int n_verts = mesh.n_vertices();
76
77 double val = 0;
79 {
80 for (int b = 0; b < adj.rows(); b++)
81 {
83 s.setZero(dim);
84 double sum_norm = 0;
85 int valence = 0;
86 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
87 {
88 assert(it.col() != b);
89 auto x = mesh.point(b) - mesh.point(it.col());
90 s += x;
91 sum_norm += x.norm();
92 valence += 1;
93 }
94 if (valence)
95 {
96 s = s / sum_norm;
97 val += pow(s.norm(), power_);
98 }
99 }
100 }
101 else
102 {
103 Eigen::MatrixXd V;
105
106 val = (L * V).eval().squaredNorm();
107 }
108
109 return val;
110 }
111
112 void BoundarySmoothingForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
113 {
114 const auto &mesh = *(state_.mesh);
115 const int dim = mesh.dimension();
116 const int n_verts = mesh.n_vertices();
117
118 Eigen::VectorXd grad;
120 {
121 grad.setZero(n_verts * dim);
122 for (int b = 0; b < adj.rows(); b++)
123 {
125 s.setZero(dim);
126 double sum_norm = 0;
127 auto sum_normalized = s;
128 int valence = 0;
129 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
130 {
131 assert(it.col() != b);
132 auto x = mesh.point(b) - mesh.point(it.col());
133 s += x;
134 sum_norm += x.norm();
135 sum_normalized += x.normalized();
136 valence += 1;
137 }
138 if (valence)
139 {
140 s = s / sum_norm;
141 const double coeff = power_ * pow(s.norm(), power_ - 2.) / sum_norm;
142
143 grad.segment(b * dim, dim) += (s * valence - s.squaredNorm() * sum_normalized) * coeff;
144 for (Eigen::SparseMatrix<bool, Eigen::RowMajor>::InnerIterator it(adj, b); it; ++it)
145 grad.segment(it.col() * dim, dim) -= (s + s.squaredNorm() * (mesh.point(it.col()) - mesh.point(b)).normalized()) * coeff;
146 }
147 }
148 }
149 else
150 {
151 Eigen::MatrixXd V;
153
154 grad = utils::flatten(2 * (L.transpose() * (L * V)));
155 }
156
158 return grad;
159 });
160 }
161} // namespace polyfem::solver
int V
double val
Definition Assembler.cpp:86
int x
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
virtual double weight() const
Get the form's multiplicative constant weight.
Definition Form.hpp:126
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:11