Graphite Version 3
An experimental 3D geometry processing program
Loading...
Searching...
No Matches
delaunay.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_DELAUNAY_DELAUNAY
41#define GEOGRAM_DELAUNAY_DELAUNAY
42
48#include <stdexcept>
49
55namespace GEO {
56
57 class Mesh;
58
59 /************************************************************************/
60
71 class GEOGRAM_API Delaunay : public Counted {
72 public:
79 struct InvalidDimension : std::logic_error {
87 coord_index_t dimension,
88 const char* name,
89 const char* expected
90 );
91
95 const char* what() const GEO_NOEXCEPT override;
96 };
97
98
104 struct InvalidInput : std::logic_error {
105
110 InvalidInput(int error_code_in);
111
117
118 ~InvalidInput() GEO_NOEXCEPT override;
119
123 const char* what() const GEO_NOEXCEPT override;
124
128 int error_code;
129
134 vector<index_t> invalid_facets;
135 };
136
155 static Delaunay* create(
156 coord_index_t dim, const std::string& name = "default"
157 );
158
159
165 static void initialize();
166
171 coord_index_t dimension() const {
172 return dimension_;
173 }
174
181 return cell_size_;
182 }
183
190 virtual void set_vertices(
191 index_t nb_vertices, const double* vertices
192 );
193
194
204 void set_reorder(bool x) {
205 do_reorder_ = x;
206 }
207
208
218 virtual void set_BRIO_levels(const vector<index_t>& levels);
219
224 const double* vertices_ptr() const {
225 return vertices_;
226 }
227
233 const double* vertex_ptr(index_t i) const {
234 geo_debug_assert(i < nb_vertices());
235 return vertices_ + vertex_stride_ * i;
236 }
237
243 return nb_vertices_;
244 }
245
252 virtual bool supports_constraints() const;
253
263 virtual void set_constraints(const Mesh* mesh) {
264 geo_assert(supports_constraints());
265 constraints_ = mesh;
266 }
267
277 void set_refine(bool x) {
278 refine_ = x;
279 }
280
287 bool get_refine() const {
288 return refine_;
289 }
290
302 void set_quality(double qual) {
303 quality_ = qual;
304 }
305
311 const Mesh* constraints() const {
312 return constraints_;
313 }
314
320 return nb_cells_;
321 }
322
332 geo_debug_assert(keeps_infinite());
333 return nb_finite_cells_;
334 }
335
340 const signed_index_t* cell_to_v() const {
341 return cell_to_v_;
342 }
343
349 return cell_to_cell_;
350 }
351
357 virtual index_t nearest_vertex(const double* p) const;
358
366 geo_debug_assert(c < nb_cells());
367 geo_debug_assert(lv < cell_size());
368 return cell_to_v_[c * cell_v_stride_ + lv];
369 }
370
380 geo_debug_assert(c < nb_cells());
381 geo_debug_assert(lf < cell_size());
382 return cell_to_cell_[c * cell_neigh_stride_ + lf];
383 }
384
392
399 bool cell_is_finite(index_t c) const {
400 return !cell_is_infinite(c);
401 }
402
412 geo_debug_assert(c < nb_cells());
413 geo_debug_assert(v < (signed_index_t) nb_vertices());
414 for(index_t iv = 0; iv < cell_size(); iv++) {
415 if(cell_vertex(c, iv) == v) {
416 return iv;
417 }
418 }
420 }
421
432 geo_debug_assert(c1 < nb_cells());
433 geo_debug_assert(c2 < nb_cells());
434 for(index_t f = 0; f < cell_size(); f++) {
435 if(cell_adjacent(c1, f) == signed_index_t(c2)) {
436 return f;
437 }
438 }
440 }
441
450 geo_debug_assert(v < nb_vertices());
451 geo_debug_assert(v < v_to_cell_.size());
452 return v_to_cell_[v];
453 }
454
455
466 geo_debug_assert(c < nb_cells());
467 geo_debug_assert(lv < cell_size());
468 return cicl_[cell_size() * c + lv];
469 }
470
482 index_t v, vector<index_t>& neighbors
483 ) const {
484 geo_debug_assert(v < nb_vertices());
485 if(store_neighbors_) {
486 neighbors_.get_array(v, neighbors);
487 } else {
488 get_neighbors_internal(v, neighbors);
489 }
490 }
491
497 void save_histogram(std::ostream& out) const;
498
507 bool stores_neighbors() const {
508 return store_neighbors_;
509 }
510
518 void set_stores_neighbors(bool x) {
519 store_neighbors_ = x;
520 if(store_neighbors_) {
521 set_stores_cicl(true);
522 }
523 }
524
531 bool stores_cicl() const {
532 return store_cicl_;
533 }
534
541 void set_stores_cicl(bool x) {
542 store_cicl_ = x;
543 }
544
545
551 bool keeps_infinite() const {
552 return keep_infinite_;
553 }
554
563 void set_keeps_infinite(bool x) {
564 keep_infinite_ = x;
565 }
566
571 bool thread_safe() const {
572 return neighbors_.thread_safe();
573 }
574
580 void set_thread_safe(bool x) {
581 neighbors_.set_thread_safe(x);
582 }
583
592 default_nb_neighbors_ = x;
593 }
594
603 return default_nb_neighbors_;
604 }
605
610 neighbors_.clear();
611 }
612
619 void set_keep_regions(bool x) {
620 keep_regions_ = x;
621 }
622
630 virtual index_t region(index_t t) const;
631
632
633 protected:
646
650 ~Delaunay() override;
651
659 index_t v, vector<index_t>& neighbors
660 ) const;
661
669 virtual void set_arrays(
670 index_t nb_cells,
671 const signed_index_t* cell_to_v, const signed_index_t* cell_to_cell
672 );
673
677 virtual void update_v_to_cell();
678
683 virtual void update_cicl();
684
688 virtual void update_neighbors();
689
697 index_t c1, index_t lv, index_t c2
698 ) {
699 geo_debug_assert(c1 < nb_cells());
700 geo_debug_assert(c2 < nb_cells());
701 geo_debug_assert(lv < cell_size());
702 cicl_[cell_size() * c1 + lv] = signed_index_t(c2);
703 }
704
705 public:
714
715 protected:
729 dimension_ = dim;
730 vertex_stride_ = dim;
731 cell_size_ = index_t(dim) + 1;
732 cell_v_stride_ = cell_size_;
733 cell_neigh_stride_ = cell_size_;
734 }
735
736 coord_index_t dimension_;
737 index_t vertex_stride_;
738 index_t cell_size_;
739 index_t cell_v_stride_;
740 index_t cell_neigh_stride_;
741
742 const double* vertices_;
743 index_t nb_vertices_;
744 index_t nb_cells_;
745 const signed_index_t* cell_to_v_;
746 const signed_index_t* cell_to_cell_;
747 vector<signed_index_t> v_to_cell_;
749 bool is_locked_;
750 PackedArrays neighbors_;
751 bool store_neighbors_;
752 index_t default_nb_neighbors_;
753
759
760 const Mesh* constraints_;
761
762 bool refine_;
763 double quality_;
764
770
776
784
785 bool keep_regions_;
786 };
787
793
805
811#define geo_register_Delaunay_creator(type, name) \
812 geo_register_creator(GEO::DelaunayFactory, type, name)
813}
814
815#endif
#define geo_assert_not_reached
Sets a non reachable point in the program.
Definition assert.h:177
#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
Base class for reference-counted objects.
Definition counted.h:71
Abstract interface for Delaunay triangulation in Nd.
Definition delaunay.h:71
index_t adjacent_index(index_t c1, index_t c2) const
Retrieves a local facet index from two adacent cell global indices.
Definition delaunay.h:431
void save_histogram(std::ostream &out) const
Saves the histogram of vertex degree (can be visualized with gnuplot).
void set_refine(bool x)
Specifies whether the mesh should be refined.
Definition delaunay.h:277
signed_index_t cell_vertex(index_t c, index_t lv) const
Gets a vertex index by cell index and local vertex index.
Definition delaunay.h:365
virtual index_t nearest_vertex(const double *p) const
Computes the nearest vertex from a query point.
void set_stores_neighbors(bool x)
Specifies whether neighbors should be stored.
Definition delaunay.h:518
virtual void get_neighbors_internal(index_t v, vector< index_t > &neighbors) const
Internal implementation for get_neighbors (with vector).
void set_stores_cicl(bool x)
Specifies whether incident tetrahedra lists should be stored.
Definition delaunay.h:541
index_t nb_cells() const
Gets the number of cells.
Definition delaunay.h:319
signed_index_t cell_adjacent(index_t c, index_t lf) const
Gets an adjacent cell index by cell index and local facet index.
Definition delaunay.h:379
void get_neighbors(index_t v, vector< index_t > &neighbors) const
Gets the one-ring neighbors of vertex v.
Definition delaunay.h:481
bool do_reorder_
If true, uses BRIO reordering (in some implementations)
Definition delaunay.h:758
void set_thread_safe(bool x)
Specifies whether thread-safe mode should be used.
Definition delaunay.h:580
const Mesh * constraints() const
Gets the constraints.
Definition delaunay.h:311
index_t nb_finite_cells_
If keep_infinite_ is true, then finite cells are 0..nb_finite_cells_-1 and infinite cells are nb_fini...
Definition delaunay.h:783
virtual bool supports_constraints() const
Tests whether constraints are supported by this Delaunay.
void set_default_nb_neighbors(index_t x)
Sets the default number of stored neighbors.
Definition delaunay.h:591
bool store_cicl_
It true, circular incident tet lists are stored.
Definition delaunay.h:769
void set_quality(double qual)
Specifies the desired quality for mesh elements when refinement is enabled (.
Definition delaunay.h:302
virtual void update_v_to_cell()
Stores for each vertex v a cell incident to v.
void set_next_around_vertex(index_t c1, index_t lv, index_t c2)
Sets the circular incident edge list.
Definition delaunay.h:696
index_t cell_size() const
Gets the number of vertices in each cell.
Definition delaunay.h:180
index_t nb_vertices() const
Gets the number of vertices.
Definition delaunay.h:242
bool stores_cicl() const
Tests whether incident tetrahedra lists are stored.
Definition delaunay.h:531
void set_reorder(bool x)
Specifies whether vertices should be reordered.
Definition delaunay.h:204
signed_index_t vertex_cell(index_t v) const
Gets an incident cell index by a vertex index.
Definition delaunay.h:449
bool get_refine() const
Tests whether mesh refinement is selected.
Definition delaunay.h:287
bool keep_infinite_
If true, infinite vertex and infinite simplices are kept.
Definition delaunay.h:775
bool thread_safe() const
Tests whether thread-safe mode is active.
Definition delaunay.h:571
const signed_index_t * cell_to_v() const
Gets a pointer to the cell-to-vertex incidence array.
Definition delaunay.h:340
virtual void set_vertices(index_t nb_vertices, const double *vertices)
Sets the vertices of this Delaunay, and recomputes the cells.
bool cell_is_infinite(index_t c) const
Tests whether a cell is infinite.
const double * vertices_ptr() const
Gets a pointer to the array of vertices.
Definition delaunay.h:224
~Delaunay() override
Delaunay destructor.
index_t nb_finite_cells() const
Gets the number of finite cells.
Definition delaunay.h:331
bool cell_is_finite(index_t c) const
Tests whether a cell is finite.
Definition delaunay.h:399
Delaunay(coord_index_t dimension)
Creates a new Delaunay triangulation.
void set_keeps_infinite(bool x)
Sets whether infinite elements should be kept.
Definition delaunay.h:563
virtual void set_arrays(index_t nb_cells, const signed_index_t *cell_to_v, const signed_index_t *cell_to_cell)
Sets the arrays that represent the combinatorics of this Delaunay.
SmartPointer< Delaunay > Delaunay_var
Smart pointer that refers to a Delaunay object.
Definition delaunay.h:792
void set_keep_regions(bool x)
Specifies whether all internal regions should be kept.
Definition delaunay.h:619
index_t index(index_t c, signed_index_t v) const
Retrieves a local vertex index from cell index and global vertex index.
Definition delaunay.h:411
virtual void update_cicl()
Updates the circular incident cell lists.
virtual index_t region(index_t t) const
Gets the region id associated with a tetrahedron.
bool keeps_infinite() const
Tests whether infinite elements are kept.
Definition delaunay.h:551
signed_index_t next_around_vertex(index_t c, index_t lv) const
Traverses the list of cells incident to a vertex.
Definition delaunay.h:465
void clear_neighbors()
Frees all memory used for neighbors storage.
Definition delaunay.h:609
Factory1< Delaunay, coord_index_t > DelaunayFactory
Delaunay Factory.
Definition delaunay.h:804
virtual void set_BRIO_levels(const vector< index_t > &levels)
Specifies the bounds of each level to be used when hierarchic ordering is specified from outside.
bool stores_neighbors() const
Tests whether neighbors are stored.
Definition delaunay.h:507
const double * vertex_ptr(index_t i) const
Gets a pointer to a vertex by its global index.
Definition delaunay.h:233
virtual void store_neighbors_CB(index_t i)
Stores the neighbors of a vertex.
void set_dimension(coord_index_t dim)
Sets the dimension of this Delaunay.
Definition delaunay.h:728
virtual void update_neighbors()
Computes the stored neighbor lists.
virtual void set_constraints(const Mesh *mesh)
Defines the constraints.
Definition delaunay.h:263
index_t default_nb_neighbors() const
Gets the default number of stored neighbors.
Definition delaunay.h:602
const signed_index_t * cell_to_cell() const
Gets a pointer to the cell-to-cell adjacency array.
Definition delaunay.h:348
Factory for types with one constructor argument.
Definition factory.h:345
Represents a mesh.
Definition mesh.h:2701
Efficient storage for array of arrays.
A smart pointer with reference-counted copy semantics.
Vector with aligned memory allocation.
Definition memory.h:660
Base class of reference-counted objects, to be used with smart pointers.
Common include file, providing basic definitions. Should be included before anything else by all head...
Generic factory mechanism.
Global Vorpaline namespace.
geo_signed_index_t signed_index_t
The type for storing and manipulating indices differences.
Definition numeric.h:343
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
Efficient storage for array of arrays.
Pointers with automatic reference counting.
Invalid dimension exception.
Definition delaunay.h:79
const char * what() const GEO_NOEXCEPT override
Gets the string identifying the exception.
InvalidDimension(coord_index_t dimension, const char *name, const char *expected)
Creates a invalid dimension exception.
Invalid input exception.
Definition delaunay.h:104
InvalidInput(int error_code_in)
InvalidInput constructor.
InvalidInput(const InvalidInput &rhs)
InvalidInput copy constructor.