6 #include "GmshConfig.h"
8 #if defined(HAVE_REVOROPT)
12 #include "Revoropt/Mesh/builder_def.hpp"
13 #include "Revoropt/Mesh/sampling_def.hpp"
14 #include "Revoropt/Mesh/normals_def.hpp"
16 #include "Revoropt/RVD/rvd.hpp"
17 #include "Revoropt/RVD/rdt.hpp"
19 #include "Revoropt/CVT/minimizer.hpp"
21 #include "Revoropt/Solver/alglib_lbfgs.hpp"
36 {
GMSH_FULLRC,
"Normal computation radius",
nullptr, 0.005}};
47 return "Plugin(CVTRemesh) triangulates an input geometry using"
48 "centroïdal Voronoï tesselations. The STL mesh of the geometry"
49 "is generated and randomly sampled. An objective function derived"
50 "from centroïdal Voronoï tesselations is then defined on the"
51 "generated sampling, and optimized through LBFGS to obtain a"
52 "regular sampling of the surface. The triangulation is extracted"
53 "as the restricted Delaunay triangulation of the samples and the"
55 "If `View' < 0, the plugin is run on the current view.\n\n"
56 "Plugin(CVTRemesh) creates one new view.";
61 return sizeof(CVTRemeshOptions_Number) /
sizeof(
StringXNumber);
66 return &CVTRemeshOptions_Number[iopt];
70 class SolverCallback {
72 template <
typename Data>
void operator()(Data *data)
75 Msg::Info(
"[CVTRemesh] : iteration %d, objective function value is %f",
86 std::vector<double> vertices;
87 std::vector<unsigned int>
faces;
89 unsigned int offset = 0;
91 (*it)->buildSTLTriangulation();
92 for(std::size_t i = 0; i < (*it)->stl_vertices_uv.size(); ++i) {
93 GPoint p = (*it)->point((*it)->stl_vertices_uv[i]);
94 vertices.push_back(p.
x());
95 vertices.push_back(p.
y());
96 vertices.push_back(p.
z());
98 for(std::size_t i = 0; i < (*it)->stl_triangles.size(); ++i) {
99 faces.push_back((*it)->stl_triangles[i] + offset);
101 offset += (*it)->stl_vertices_uv.size();
104 Revoropt::MeshBuilder<3> mesh;
105 mesh.swap_vertices(vertices);
106 mesh.swap_faces(
faces);
108 double mesh_center[3];
110 Revoropt::normalize_mesh(&mesh, mesh_center, &mesh_scale);
112 double nradius = (double)CVTRemeshOptions_Number[5].def;
115 std::vector<double> normals(3 * mesh.vertices_size());
116 Revoropt::full_robust_vertex_normals(&mesh, nradius, normals.data());
119 std::vector<double> lifted_vertices(6 * mesh.vertices_size(), 0);
120 for(std::size_t vertex = 0; vertex < mesh.vertices_size(); ++vertex) {
121 std::copy(mesh.vertex(vertex), mesh.vertex(vertex) + 3,
122 lifted_vertices.data() + 6 * vertex);
123 std::copy(normals.data() + 3 * vertex, normals.data() + 3 * vertex + 3,
124 lifted_vertices.data() + 6 * vertex + 3);
128 Revoropt::ROMeshWrapper<3, 6> lifted_mesh(lifted_vertices.data(),
129 lifted_vertices.size() / 6, &mesh);
132 double twfactor = (double)CVTRemeshOptions_Number[3].def;
135 std::vector<double> triangle_weights(lifted_mesh.faces_size());
137 for(std::size_t
f = 0;
f < lifted_mesh.faces_size(); ++
f) {
139 const unsigned int *fverts = mesh.face(
f);
143 for(
int i = 0; i < 3; ++i) { x[i] = lifted_mesh.vertex(fverts[i]); }
149 typedef Eigen::Matrix<double, 3, 1> Vector3;
151 Eigen::Map<const Vector3> v0(x[0]);
152 Eigen::Map<const Vector3> v1(x[1]);
153 Eigen::Map<const Vector3> v2(x[2]);
156 Vector3 U = (v1 - v0);
157 const double U_len = U.norm();
160 Vector3 H = (v2 - v0);
161 H = H - H.dot(U) * U;
162 const double H_len = H.norm();
168 Eigen::Matrix<double, 3, 2> bar_grads;
170 bar_grads(2, 1) = 1 / H_len;
173 for(
int i = 0; i < 2; ++i) {
175 Eigen::Map<const Vector3> vi0(x[(i + 1) % 3]);
176 Eigen::Map<const Vector3> vi1(x[(i + 2) % 3]);
177 Eigen::Map<const Vector3> vi2(x[i]);
179 Vector3 Ui = (vi1 - vi0);
181 Vector3 Hi = (vi2 - vi0);
182 Hi = Hi - Hi.dot(Ui) * Ui;
183 const double Hi_invlen = 1 / Hi.norm();
185 bar_grads(i, 0) = Hi.dot(U) * Hi_invlen;
186 bar_grads(i, 1) = Hi.dot(H) * Hi_invlen;
190 Eigen::Map<const Vector3> n0(x[0] + 3);
191 Eigen::Map<const Vector3> n1(x[1] + 3);
192 Eigen::Map<const Vector3> n2(x[2] + 3);
194 Eigen::Matrix<double, 3, 2> n_grads =
195 Eigen::Matrix<double, 3, 2>::Zero();
197 n_grads = n0 * bar_grads.row(0);
198 n_grads += n1 * bar_grads.row(1);
199 n_grads += n2 * bar_grads.row(2);
202 double g_max = n_grads.row(0).dot(n_grads.row(0));
203 double g_other = n_grads.row(1).dot(n_grads.row(1));
204 g_max = g_max > g_other ? g_max : g_other;
205 g_other = n_grads.row(2).dot(n_grads.row(2));
206 g_max = g_max > g_other ? g_max : g_other;
213 triangle_weights[
f] = pow(ratio, twfactor);
218 double nfactor = (double)CVTRemeshOptions_Number[2].def;
222 for(std::size_t i = 0; i < lifted_mesh.vertices_size(); ++i) {
223 double *v = lifted_vertices.data() + 6 * i;
230 unsigned int nsites = (
unsigned int)CVTRemeshOptions_Number[0].def;
233 std::vector<double> lifted_sites(6 * nsites);
235 Revoropt::generate_random_sites<Revoropt::ROMesh<3, 6> >(
236 &lifted_mesh, nsites, lifted_sites.data(), triangle_weights.data());
239 Revoropt::generate_random_sites<Revoropt::ROMesh<3, 6> >(
240 &lifted_mesh, nsites, lifted_sites.data());
244 Revoropt::CVT::DirectMinimizer<Revoropt::ROMesh<3, 6> > cvt;
245 cvt.set_sites(lifted_sites.data(), nsites);
246 cvt.set_mesh(&lifted_mesh);
247 if(twfactor > 0) { cvt.set_triangle_weights(triangle_weights.data()); }
250 SolverCallback callback;
253 unsigned int niter = (
unsigned int)CVTRemeshOptions_Number[1].def;
255 unsigned int aniso_niter = std::min<unsigned int>(10, niter);
262 aniso_niter = std::max(aniso_niter, niter * 10 / 100);
263 cvt.set_anisotropy(1);
265 cvt.minimize<Revoropt::Solver::AlgLBFGS>(niter - aniso_niter, &callback);
271 double tanisotropy = (double)CVTRemeshOptions_Number[4].def;
275 cvt.set_anisotropy(tanisotropy);
276 status = cvt.minimize<Revoropt::Solver::AlgLBFGS>(aniso_niter, &callback);
280 std::vector<unsigned int> rdt_triangles;
281 Revoropt::RDTBuilder<Revoropt::ROMesh<3, 6> > build_rdt(rdt_triangles);
282 Revoropt::RVD<Revoropt::ROMesh<3, 6> > rvd;
283 rvd.set_sites(lifted_sites.data(), nsites);
284 rvd.set_mesh(&lifted_mesh);
285 rvd.compute(build_rdt);
291 std::vector<MVertex *> m_verts(nsites);
292 for(std::size_t i = 0; i < nsites; ++i) {
294 new MVertex(lifted_sites[6 * i] * mesh_scale + mesh_center[0],
295 lifted_sites[6 * i + 1] * mesh_scale + mesh_center[1],
296 lifted_sites[6 * i + 2] * mesh_scale + mesh_center[2]);
299 for(std::size_t i = 0; i < rdt_triangles.size() / 3; ++i) {
301 m_verts[rdt_triangles[3 * i + 1]],
302 m_verts[rdt_triangles[3 * i + 2]]));