Geogram Version 1.9.6-rc
A programming library of geometric algorithms
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(index_t nb_vertices, const double* vertices);
191
201 void set_reorder(bool x) {
202 do_reorder_ = x;
203 }
204
214 virtual void set_BRIO_levels(const vector<index_t>& levels);
215
220 const double* vertices_ptr() const {
221 return vertices_;
222 }
223
229 const double* vertex_ptr(index_t i) const {
230 geo_debug_assert(i < nb_vertices());
231 return vertices_ + vertex_stride_ * i;
232 }
233
239 return nb_vertices_;
240 }
241
248 virtual bool supports_constraints() const;
249
259 virtual void set_constraints(const Mesh* mesh) {
260 geo_assert(supports_constraints());
261 constraints_ = mesh;
262 }
263
273 void set_refine(bool x) {
274 refine_ = x;
275 }
276
283 bool get_refine() const {
284 return refine_;
285 }
286
298 void set_quality(double qual) {
299 quality_ = qual;
300 }
301
307 const Mesh* constraints() const {
308 return constraints_;
309 }
310
316 return nb_cells_;
317 }
318
328 geo_debug_assert(keeps_infinite());
329 return nb_finite_cells_;
330 }
331
336 const index_t* cell_to_v() const {
337 return cell_to_v_;
338 }
339
344 const index_t* cell_to_cell() const {
345 return cell_to_cell_;
346 }
347
353 virtual index_t nearest_vertex(const double* p) const;
354
362 geo_debug_assert(c < nb_cells());
363 geo_debug_assert(lv < cell_size());
364 return cell_to_v_[c * cell_v_stride_ + lv];
365 }
366
376 geo_debug_assert(c < nb_cells());
377 geo_debug_assert(lf < cell_size());
378 return cell_to_cell_[c * cell_neigh_stride_ + lf];
379 }
380
388
395 bool cell_is_finite(index_t c) const {
396 return !cell_is_infinite(c);
397 }
398
408 geo_debug_assert(c < nb_cells());
409 geo_debug_assert(v == NO_INDEX || v < nb_vertices());
410 for(index_t iv = 0; iv < cell_size(); iv++) {
411 if(cell_vertex(c, iv) == v) {
412 return iv;
413 }
414 }
416 }
417
428 geo_debug_assert(c1 < nb_cells());
429 geo_debug_assert(c2 < nb_cells());
430 for(index_t f = 0; f < cell_size(); f++) {
431 if(cell_adjacent(c1, f) == c2) {
432 return f;
433 }
434 }
436 }
437
446 geo_debug_assert(v < nb_vertices());
447 geo_debug_assert(v < v_to_cell_.size());
448 return v_to_cell_[v];
449 }
450
451
462 geo_debug_assert(c < nb_cells());
463 geo_debug_assert(lv < cell_size());
464 return cicl_[cell_size() * c + lv];
465 }
466
477 void get_neighbors(index_t v, vector<index_t>& neighbors) const {
478 geo_debug_assert(v < nb_vertices());
479 if(store_neighbors_) {
480 neighbors_.get_array(v, neighbors);
481 } else {
482 get_neighbors_internal(v, neighbors);
483 }
484 }
485
491 void save_histogram(std::ostream& out) const;
492
501 bool stores_neighbors() const {
502 return store_neighbors_;
503 }
504
512 void set_stores_neighbors(bool x) {
513 store_neighbors_ = x;
514 if(store_neighbors_) {
515 set_stores_cicl(true);
516 }
517 }
518
525 bool stores_cicl() const {
526 return store_cicl_;
527 }
528
535 void set_stores_cicl(bool x) {
536 store_cicl_ = x;
537 }
538
539
545 bool keeps_infinite() const {
546 return keep_infinite_;
547 }
548
557 void set_keeps_infinite(bool x) {
558 keep_infinite_ = x;
559 }
560
565 bool thread_safe() const {
566 return neighbors_.thread_safe();
567 }
568
574 void set_thread_safe(bool x) {
575 neighbors_.set_thread_safe(x);
576 }
577
586 default_nb_neighbors_ = x;
587 }
588
597 return default_nb_neighbors_;
598 }
599
604 neighbors_.clear();
605 }
606
613 void set_keep_regions(bool x) {
614 keep_regions_ = x;
615 }
616
624 virtual index_t region(index_t t) const;
625
626
627 protected:
640
644 ~Delaunay() override;
645
653 index_t v, vector<index_t>& neighbors
654 ) const;
655
663 virtual void set_arrays(
664 index_t nb_cells,
665 const index_t* cell_to_v, const index_t* cell_to_cell
666 );
667
671 virtual void update_v_to_cell();
672
677 virtual void update_cicl();
678
682 virtual void update_neighbors();
683
691 index_t c1, index_t lv, index_t c2
692 ) {
693 geo_debug_assert(c1 < nb_cells());
694 geo_debug_assert(c2 < nb_cells());
695 geo_debug_assert(lv < cell_size());
696 cicl_[cell_size() * c1 + lv] = c2;
697 }
698
699 public:
708
709 protected:
723 dimension_ = dim;
724 vertex_stride_ = dim;
725 cell_size_ = index_t(dim) + 1;
726 cell_v_stride_ = cell_size_;
727 cell_neigh_stride_ = cell_size_;
728 }
729
730 coord_index_t dimension_;
731 index_t vertex_stride_;
732 index_t cell_size_;
733 index_t cell_v_stride_;
734 index_t cell_neigh_stride_;
735
736 const double* vertices_;
737 index_t nb_vertices_;
738 index_t nb_cells_;
739 const index_t* cell_to_v_;
740 const index_t* cell_to_cell_;
741 vector<index_t> v_to_cell_;
742 vector<index_t> cicl_;
743 bool is_locked_;
744 PackedArrays neighbors_;
745 bool store_neighbors_;
746 index_t default_nb_neighbors_;
747
753
754 const Mesh* constraints_;
755
756 bool refine_;
757 double quality_;
758
764
770
778
779 bool keep_regions_;
780 };
781
787
799
805#define geo_register_Delaunay_creator(type, name) \
806 geo_register_creator(GEO::DelaunayFactory, type, name)
807}
808
809#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:427
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:273
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:512
virtual void get_neighbors_internal(index_t v, vector< index_t > &neighbors) const
Internal implementation for get_neighbors (with vector).
const index_t * cell_to_v() const
Gets a pointer to the cell-to-vertex incidence array.
Definition delaunay.h:336
void set_stores_cicl(bool x)
Specifies whether incident tetrahedra lists should be stored.
Definition delaunay.h:535
const index_t * cell_to_cell() const
Gets a pointer to the cell-to-cell adjacency array.
Definition delaunay.h:344
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:375
index_t nb_cells() const
Gets the number of cells.
Definition delaunay.h:315
void get_neighbors(index_t v, vector< index_t > &neighbors) const
Gets the one-ring neighbors of vertex v.
Definition delaunay.h:477
bool do_reorder_
If true, uses BRIO reordering (in some implementations)
Definition delaunay.h:752
void set_thread_safe(bool x)
Specifies whether thread-safe mode should be used.
Definition delaunay.h:574
const Mesh * constraints() const
Gets the constraints.
Definition delaunay.h:307
index_t vertex_cell(index_t v) const
Gets an incident cell index by a vertex index.
Definition delaunay.h:445
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:777
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:585
bool store_cicl_
It true, circular incident tet lists are stored.
Definition delaunay.h:763
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:361
void set_quality(double qual)
Specifies the desired quality for mesh elements when refinement is enabled (.
Definition delaunay.h:298
virtual void set_arrays(index_t nb_cells, const index_t *cell_to_v, const index_t *cell_to_cell)
Sets the arrays that represent the combinatorics of this Delaunay.
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:690
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:238
bool stores_cicl() const
Tests whether incident tetrahedra lists are stored.
Definition delaunay.h:525
void set_reorder(bool x)
Specifies whether vertices should be reordered.
Definition delaunay.h:201
bool get_refine() const
Tests whether mesh refinement is selected.
Definition delaunay.h:283
bool keep_infinite_
If true, infinite vertex and infinite simplices are kept.
Definition delaunay.h:769
bool thread_safe() const
Tests whether thread-safe mode is active.
Definition delaunay.h:565
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:220
~Delaunay() override
Delaunay destructor.
index_t nb_finite_cells() const
Gets the number of finite cells.
Definition delaunay.h:327
bool cell_is_finite(index_t c) const
Tests whether a cell is finite.
Definition delaunay.h:395
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:557
SmartPointer< Delaunay > Delaunay_var
Smart pointer that refers to a Delaunay object.
Definition delaunay.h:786
index_t next_around_vertex(index_t c, index_t lv) const
Traverses the list of cells incident to a vertex.
Definition delaunay.h:461
void set_keep_regions(bool x)
Specifies whether all internal regions should be kept.
Definition delaunay.h:613
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:545
void clear_neighbors()
Frees all memory used for neighbors storage.
Definition delaunay.h:603
Factory1< Delaunay, coord_index_t > DelaunayFactory
Delaunay Factory.
Definition delaunay.h:798
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:501
const double * vertex_ptr(index_t i) const
Gets a pointer to a vertex by its global index.
Definition delaunay.h:229
virtual void store_neighbors_CB(index_t i)
Stores the neighbors of a vertex.
index_t index(index_t c, index_t v) const
Retrieves a local vertex index from cell index and global vertex index.
Definition delaunay.h:407
void set_dimension(coord_index_t dim)
Sets the dimension of this Delaunay.
Definition delaunay.h:722
virtual void update_neighbors()
Computes the stored neighbor lists.
virtual void set_constraints(const Mesh *mesh)
Defines the constraints.
Definition delaunay.h:259
index_t default_nb_neighbors() const
Gets the default number of stored neighbors.
Definition delaunay.h:596
Factory for types with one constructor argument.
Definition factory.h:345
Represents a mesh.
Definition mesh.h:3050
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.
Generic factory mechanism.
Common include file, providing basic definitions. Should be included before anything else by all head...
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
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.