Geogram Version 1.9.6-rc
A programming library of geometric algorithms
Loading...
Searching...
No Matches
CVT.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2000-2022 Inria
3 * All rights reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * * Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 * * Redistributions in binary form must reproduce the above copyright notice,
11 * this list of conditions and the following disclaimer in the documentation
12 * and/or other materials provided with the distribution.
13 * * Neither the name of the ALICE Project-Team nor the names of its
14 * contributors may be used to endorse or promote products derived from this
15 * software without specific prior written permission.
16 *
17 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
21 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27 * POSSIBILITY OF SUCH DAMAGE.
28 *
29 * Contact: Bruno Levy
30 *
31 * https://www.inria.fr/fr/bruno-levy
32 *
33 * Inria,
34 * Domaine de Voluceau,
35 * 78150 Le Chesnay - Rocquencourt
36 * FRANCE
37 *
38 */
39
40#ifndef GEOGRAM_VORONOI_CVT
41#define GEOGRAM_VORONOI_CVT
42
44#include <geogram/voronoi/RVD.h>
46#include <geogram/mesh/mesh.h>
48
49
55namespace GEO {
56
57 class RestrictedVoronoiDiagram;
58 class ProgressTask;
59
69 class GEOGRAM_API CentroidalVoronoiTesselation {
70
73
74 public:
88 Mesh* mesh,
89 coord_index_t dimension = 0,
90 const std::string& delaunay = "default"
91 );
92
109 Mesh* mesh,
110 const vector<vec3>& R3_embedding, coord_index_t dimension = 0,
111 const std::string& delaunay = "default"
112 );
113
118
129 bool compute_initial_sampling(index_t nb_samples, bool verbose=false);
130
139 void set_points(index_t nb_points, const double* points);
140
146 void resize_points(index_t nb_points);
147
157 virtual void Lloyd_iterations(index_t nb_iter);
158
164 virtual void Newton_iterations(index_t nb_iter, index_t m = 7);
165
172 void compute_surface(Mesh* mesh, bool multinerve = true);
173
179 void compute_volume(Mesh* mesh);
180
185 void set_show_iterations(bool x) {
186 show_iterations_ = x;
187 }
188
196 use_RVC_centroids_ = x;
197 }
198
204 void set_constrained_cvt(bool x) {
205 constrained_cvt_ = x;
206 }
207
212 return mesh_;
213 }
214
219 return delaunay_;
220 }
221
226 return RVD_;
227 }
228
236 void set_facets_range(index_t facets_begin, index_t facets_end) {
237 RVD_->set_facets_range(facets_begin, facets_end);
238 }
239
250 geo_assert(instance_ == nullptr);
251 instance_ = this;
252 }
253
264 geo_assert(instance_ == this);
265 instance_ = nullptr;
266 }
267
268 public:
277 static void funcgrad_CB(
278 index_t n, double* x, double& f, double* g
279 );
280
290 static void newiteration_CB(
291 index_t n, const double* x, double f, const double* g, double gnorm
292 );
293
299 progress_ = progress;
300 }
301
307 return dimension_;
308 }
309
314 return index_t(points_.size() / dimension_);
315 }
316
323 const vec3& R3_embedding(index_t p) const {
324 return RVD_->R3_embedding(p);
325 }
326
333 double* embedding(index_t p) {
334 geo_debug_assert(p < nb_points());
335 return &(points_[0]) + dimension_ * p;
336 }
337
341 bool volumetric() const {
342 return RVD_->volumetric();
343 }
344
350 void set_volumetric(bool x) {
351 RVD_->set_volumetric(x);
352 }
353
361 bool point_is_locked(index_t i) const {
363 point_is_locked_.size() == 0 || i < point_is_locked_.size()
364 );
365 return point_is_locked_.size() != 0 && point_is_locked_[i];
366 }
367
376 geo_debug_assert(i < nb_points());
377 if(point_is_locked_.size() != nb_points()) {
378 point_is_locked_.resize(nb_points(), false);
379 }
380 point_is_locked_[i] = true;
381 }
382
391 geo_debug_assert(i < nb_points());
392 if(
393 point_is_locked_.size() != nb_points()
394 ) {
395 point_is_locked_.resize(nb_points(), false);
396 }
397 point_is_locked_[i] = false;
398 }
399
406 point_is_locked_.clear();
407 }
408
409 protected:
414 virtual void newiteration();
415
423 virtual void funcgrad(index_t n, double* x, double& f, double* g);
424
431 void constrain_points(double* g) const;
432
439
440 static CentroidalVoronoiTesselation* instance_;
441 bool show_iterations_;
442 coord_index_t dimension_;
443 Delaunay_var delaunay_;
445 Mesh* mesh_;
446
447 vector<double> points_;
448 vector<vec3> points_R3_;
449 vector<bool> point_is_locked_;
450
451 ProgressTask* progress_;
452 index_t cur_iter_;
453 index_t nb_iter_;
454
456 bool constrained_cvt_;
457 bool use_RVC_centroids_;
458
462 private:
465
467 thisclass& operator= (const thisclass& rhs);
468 };
469}
470
471#endif
Class and functions to compute restricted Voronoi diagrams and extract information from them.
#define geo_assert(x)
Verifies that a condition is met.
Definition assert.h:149
#define geo_debug_assert(x)
Verifies that a condition is met.
Definition assert.h:196
CentroidalVoronoiTesselation is the main component of the remeshing algorithm.
Definition CVT.h:69
virtual ~CentroidalVoronoiTesselation()
Destructor.
void done_current()
Resets the current CentroidalVoronoiTesselation to nullptr.
Definition CVT.h:263
virtual void funcgrad(index_t n, double *x, double &f, double *g)
Computes the objective function and its gradient.
void set_volumetric(bool x)
Sets volumetric mode.
Definition CVT.h:350
double * embedding(index_t p)
Returns the representation of a point in embedding space.
Definition CVT.h:333
void resize_points(index_t nb_points)
Changes the number of points.
void unlock_point(index_t i)
Unlocks a point.
Definition CVT.h:390
CentroidalVoronoiTesselation(Mesh *mesh, const vector< vec3 > &R3_embedding, coord_index_t dimension=0, const std::string &delaunay="default")
Constructs a new CentroidalVoronoiTesselation.
void constrain_points(double *g) const
Constrains the locked points.
bool compute_initial_sampling(index_t nb_samples, bool verbose=false)
Computes a random initial sampling of the surface in nD.
void set_use_RVC_centroids(bool x)
Specifies whether centroids of Voronoi cells should be used.
Definition CVT.h:195
void compute_volume(Mesh *mesh)
Computes the volumetric mesh (using the current points).
const vec3 & R3_embedding(index_t p) const
Gets the representation of a point in R3.
Definition CVT.h:323
CentroidalVoronoiTesselation(Mesh *mesh, coord_index_t dimension=0, const std::string &delaunay="default")
Constructs a new CentroidalVoronoiTesselation.
void make_current()
Makes this CentroidalVoronoiTesselation the current one.
Definition CVT.h:249
void compute_R3_embedding()
Computes the 3d representation of the Nd points.
virtual void newiteration()
Callback for the numerical solver.
static void funcgrad_CB(index_t n, double *x, double &f, double *g)
Callback for the numerical solver.
coord_index_t dimension() const
Gets the dimension of the points.
Definition CVT.h:306
virtual void Newton_iterations(index_t nb_iter, index_t m=7)
Relaxes the points with Newton-Lloyd's algorithm.
virtual void Lloyd_iterations(index_t nb_iter)
Relaxes the points with Lloyd's algorithm.
bool point_is_locked(index_t i) const
Tests whether a point is locked.
Definition CVT.h:361
IntegrationSimplex_var simplex_func_
Integration simplex used by custom codes, e.g. LpCVT.
Definition CVT.h:459
void set_points(index_t nb_points, const double *points)
Initializes the points with a user-specified vector.
index_t nb_points() const
Gets the number of points to be optimized.
Definition CVT.h:313
void compute_surface(Mesh *mesh, bool multinerve=true)
Computes the surfacic mesh (using the current points).
void lock_point(index_t i)
Locks a point.
Definition CVT.h:375
static void newiteration_CB(index_t n, const double *x, double f, const double *g, double gnorm)
Callback for the numerical solver.
void unlock_all_points()
Unlocks all the points.
Definition CVT.h:405
void set_constrained_cvt(bool x)
Specifies whether constrained mode should be used.
Definition CVT.h:204
bool volumetric() const
Tests whether volumetric mode is used.
Definition CVT.h:341
void set_show_iterations(bool x)
Specifies whether a progress bar should be used.
Definition CVT.h:185
RestrictedVoronoiDiagram * RVD()
Definition CVT.h:225
void set_progress_logger(ProgressTask *progress)
Sets a client for the progress bars.
Definition CVT.h:298
void set_facets_range(index_t facets_begin, index_t facets_end)
Restricts computation to a part of the input mesh.
Definition CVT.h:236
Abstract interface for Delaunay triangulation in Nd.
Definition delaunay.h:71
Represents a mesh.
Definition mesh.h:3050
Tracks the progress of a task.
Definition progress.h:240
Computes a Restricted Voronoi Diagram (RVD).
Definition RVD.h:89
Vector with aligned memory allocation.
Definition memory.h:660
Abstract interface for Delaunay.
Common include file, providing basic definitions. Should be included before anything else by all head...
base classes for computing integrals over the cells of a restricted Voronoi diagram
The class that represents a mesh.
Global Vorpaline namespace.
Definition basic.h:55
geo_index_t index_t
The type for storing and manipulating indices.
Definition numeric.h:329
geo_coord_index_t coord_index_t
The type for storing coordinate indices, and iterating on the coordinates of a point.
Definition numeric.h:363