All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DelaunayMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::DelaunayMesh
26 
27 Description
28  The vertex and cell classes must have an index defined
29 
30 SourceFiles
31  DelaunayMeshI.H
32  DelaunayMesh.C
33  DelaunayMeshIO.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef DelaunayMesh_H
38 #define DelaunayMesh_H
39 
40 #include "Pair.H"
41 #include "HashSet.H"
42 #include "FixedList.H"
43 #include "boundBox.H"
44 #include "indexedVertex.H"
46 #include "Time.H"
47 #include "autoPtr.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 class fvMesh;
55 
56 /*---------------------------------------------------------------------------*\
57  Class DelaunayMesh Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class Triangulation>
61 class DelaunayMesh
62 :
63  public Triangulation
64 {
65 public:
66 
67  typedef typename Triangulation::Cell_handle Cell_handle;
68  typedef typename Triangulation::Vertex_handle Vertex_handle;
69  typedef typename Triangulation::Edge Edge;
70  typedef typename Triangulation::Point Point;
71  typedef typename Triangulation::Facet Facet;
72 
73  typedef typename Triangulation::Finite_vertices_iterator
75  typedef typename Triangulation::Finite_cells_iterator
77  typedef typename Triangulation::Finite_facets_iterator
79 
80  typedef HashSet
81  <
85 
86  typedef HashTable
87  <
88  label,
89  labelPair,
92 
93 
94 private:
95 
96  // Private Data
97 
98  //- Keep track of the number of vertices that have been added.
99  // This allows a unique index to be assigned to each vertex.
100  mutable label vertexCount_;
101 
102  //- Keep track of the number of cells that have been added.
103  // This allows a unique index to be assigned to each cell.
104  mutable label cellCount_;
105 
106  //- Reference to Time
107  const Time& runTime_;
108 
109  //- Spatial sort traits to use with a pair of point pointers and an int.
110  // Taken from a post on the CGAL lists: 2010-01/msg00004.html by
111  // Sebastien Loriot (Geometry Factory).
112  struct Traits_for_spatial_sort
113  :
114  public Triangulation::Geom_traits
115  {
116  typedef typename Triangulation::Geom_traits Gt;
117 
118  typedef std::pair<const typename Triangulation::Point*, label>
119  Point_3;
121  struct Less_x_3
122  {
123  bool operator()(const Point_3& p, const Point_3& q) const;
124  };
126  struct Less_y_3
127  {
128  bool operator()(const Point_3& p, const Point_3& q) const;
129  };
131  struct Less_z_3
132  {
133  bool operator()(const Point_3& p, const Point_3& q) const;
134  };
135 
136  Less_x_3 less_x_3_object() const;
137  Less_y_3 less_y_3_object() const;
138  Less_z_3 less_z_3_object() const;
139  };
140 
141 
142  // Private Member Functions
143 
144  void sortFaces
145  (
146  faceList& faces,
147  labelList& owner,
148  labelList& neighbour
149  ) const;
150 
151  void addPatches
152  (
153  const label nInternalFaces,
154  faceList& faces,
155  labelList& owner,
157  const List<DynamicList<face>>& patchFaces,
158  const List<DynamicList<label>>& patchOwners
159  ) const;
160 
161 
162 public:
163 
164  // Constructors
165 
166  //- Construct from components
167  explicit DelaunayMesh(const Time& runTime);
168 
170  (
171  const Time& runTime,
172  const word& meshName
173  );
174 
175  //- Disallow default bitwise copy construction
176  DelaunayMesh(const DelaunayMesh<Triangulation>&) = delete;
177 
178 
179  //- Destructor
180  ~DelaunayMesh();
181 
182 
183  // Member Functions
184 
185  // Access
186 
187  //- Return a reference to the Time object
188  inline const Time& time() const;
189 
190 
191  // Check
192 
193  //- Write the cpuTime to screen
194  inline void timeCheck
195  (
196  const string& description,
197  const bool check = true
198  ) const;
199 
200 
201  // Indexing functions
202 
203  //- Create a new unique cell index and return
204  inline label getNewCellIndex() const;
205 
206  //- Create a new unique vertex index and return
207  inline label getNewVertexIndex() const;
208 
209  //- Return the cell count (the next unique cell index)
210  inline label cellCount() const;
211 
212  //- Return the vertex count (the next unique vertex index)
213  inline label vertexCount() const;
214 
215  //- Set the cell count to zero
216  inline void resetCellCount();
217 
218  //- Set the vertex count to zero
219  inline void resetVertexCount();
220 
221 
222  // Triangulation manipulation functions
223 
224  //- Clear the entire triangulation
225  void reset();
226 
227  //- Insert the list of vertices (calls rangeInsertWithInfo)
229  (
230  const List<Vb>& vertices,
231  const bool reIndex
232  );
233 
234  //- Function inserting points into a triangulation and setting the
235  // index and type data of the point in the correct order. This is
236  // faster than inserting points individually.
237  //
238  // Adapted from a post on the CGAL lists: 2010-01/msg00004.html by
239  // Sebastien Loriot (Geometry Factory).
240  template<class PointIterator>
242  (
243  PointIterator begin,
244  PointIterator end,
245  bool printErrors = false,
246  bool reIndex = true
247  );
248 
249 
250  // Write
251 
252  //- Write mesh statistics to stream
253  void printInfo(Ostream& os) const;
254 
255  //- Write vertex statistics in the form of a table to stream
256  void printVertexInfo(Ostream& os) const;
257 
258  //- Create an fvMesh from the triangulation.
259  // The mesh is not parallel consistent - only used for viewing
261  (
262  const fileName& name,
263  labelTolabelPairHashTable& vertexMap,
264  labelList& cellMap,
265  const bool writeDelaunayData = true
266  ) const;
267 
268 
269  // Member Operators
270 
271  //- Disallow default bitwise assignment
272  void operator=(const DelaunayMesh<Triangulation>&) = delete;
273 };
274 
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 #include "DelaunayMeshI.H"
283 
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 
286 #ifdef NoRepository
287  #include "DelaunayMesh.C"
288 #endif
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 #endif
293 
294 // ************************************************************************* //
Map< label > rangeInsertWithInfo(PointIterator begin, PointIterator end, bool printErrors=false, bool reIndex=true)
Function inserting points into a triangulation and setting the.
A HashTable with keys but without contents.
Definition: HashSet.H:59
label cellCount() const
Return the cell count (the next unique cell index)
Definition: DelaunayMeshI.H:93
Triangulation::Facet Facet
Definition: DelaunayMesh.H:70
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A class for handling file names.
Definition: fileName.H:79
The vertex and cell classes must have an index defined.
Definition: DelaunayMesh.H:60
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
void operator=(const DelaunayMesh< Triangulation > &)=delete
Disallow default bitwise assignment.
Triangulation::Cell_handle Cell_handle
Definition: DelaunayMesh.H:66
~DelaunayMesh()
Destructor.
HashTable< label, labelPair, FixedList< label, 2 >::Hash<> > labelTolabelPairHashTable
Definition: DelaunayMesh.H:90
engineTime & runTime
label getNewVertexIndex() const
Create a new unique vertex index and return.
Definition: DelaunayMeshI.H:78
void resetVertexCount()
Set the vertex count to zero.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
void printVertexInfo(Ostream &os) const
Write vertex statistics in the form of a table to stream.
Triangulation::Finite_cells_iterator Finite_cells_iterator
Definition: DelaunayMesh.H:75
Triangulation::Edge Edge
Definition: DelaunayMesh.H:68
void printInfo(Ostream &os) const
Write mesh statistics to stream.
pointField vertices(const blockVertexList &bvl)
DelaunayMesh(const Time &runTime)
Construct from components.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
static char meshName[]
Definition: globalFoam.H:7
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A class for handling words, derived from string.
Definition: word.H:59
void timeCheck(const string &description, const bool check=true) const
Write the cpuTime to screen.
Definition: DelaunayMeshI.H:37
Triangulation::Finite_facets_iterator Finite_facets_iterator
Definition: DelaunayMesh.H:77
void reset()
Clear the entire triangulation.
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
An STL-conforming hash table.
Definition: HashTable.H:61
const Time & time() const
Return a reference to the Time object.
Definition: DelaunayMeshI.H:29
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
HashSet< Pair< label >, FixedList< label, 2 >::Hash<> > labelPairHashSet
Definition: DelaunayMesh.H:83
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void resetCellCount()
Set the cell count to zero.
label vertexCount() const
Return the vertex count (the next unique vertex index)
Triangulation::Vertex_handle Vertex_handle
Definition: DelaunayMesh.H:67
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
bool operator()(const Point_3 &p, const Point_3 &q) const
label getNewCellIndex() const
Create a new unique cell index and return.
Definition: DelaunayMeshI.H:63
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
CGAL data structures used for 3D Delaunay meshing.
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo)
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
volScalarField & p
Triangulation::Point Point
Definition: DelaunayMesh.H:69
Triangulation::Finite_vertices_iterator Finite_vertices_iterator
Definition: DelaunayMesh.H:73
autoPtr< polyMesh > createMesh(const fileName &name, labelTolabelPairHashTable &vertexMap, labelList &cellMap, const bool writeDelaunayData=true) const
Create an fvMesh from the triangulation.
Namespace for OpenFOAM.