DelaunayMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2012-2015 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).
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;
120 
121  struct Less_x_3
122  {
123  bool operator()(const Point_3& p, const Point_3& q) const;
124  };
125 
126  struct Less_y_3
127  {
128  bool operator()(const Point_3& p, const Point_3& q) const;
129  };
130 
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,
158  const List<DynamicList<label> >& patchOwners
159  ) const;
160 
161  //- Disallow default bitwise copy construct
163 
164  //- Disallow default bitwise assignment
166 
167 
168 public:
169 
170  // Constructors
171 
172  //- Construct from components
173  explicit DelaunayMesh(const Time& runTime);
174 
176  (
177  const Time& runTime,
178  const word& meshName
179  );
180 
181 
182  //- Destructor
183  ~DelaunayMesh();
184 
185 
186  // Member Functions
187 
188  // Access
189 
190  //- Return a reference to the Time object
191  inline const Time& time() const;
192 
193 
194  // Check
195 
196  //- Write the cpuTime to screen
197  inline void timeCheck
198  (
199  const string& description,
200  const bool check = true
201  ) const;
202 
203 
204  // Indexing functions
205 
206  //- Create a new unique cell index and return
207  inline label getNewCellIndex() const;
208 
209  //- Create a new unique vertex index and return
210  inline label getNewVertexIndex() const;
211 
212  //- Return the cell count (the next unique cell index)
213  inline label cellCount() const;
214 
215  //- Return the vertex count (the next unique vertex index)
216  inline label vertexCount() const;
217 
218  //- Set the cell count to zero
219  inline void resetCellCount();
220 
221  //- Set the vertex count to zero
222  inline void resetVertexCount();
223 
224 
225  // Triangulation manipulation functions
226 
227  //- Clear the entire triangulation
228  void reset();
229 
230  //- Insert the list of vertices (calls rangeInsertWithInfo)
232  (
233  const List<Vb>& vertices,
234  const bool reIndex
235  );
236 
237  //- Function inserting points into a triangulation and setting the
238  // index and type data of the point in the correct order. This is
239  // faster than inserting points individually.
240  //
241  // Adapted from a post on the CGAL lists: 2010-01/msg00004.html by
242  // Sebastien Loriot (Geometry Factory).
243  template<class PointIterator>
245  (
246  PointIterator begin,
247  PointIterator end,
248  bool printErrors = false,
249  bool reIndex = true
250  );
251 
252 
253  // Write
254 
255  //- Write mesh statistics to stream
256  void printInfo(Ostream& os) const;
257 
258  //- Write vertex statistics in the form of a table to stream
259  void printVertexInfo(Ostream& os) const;
260 
261  //- Create an fvMesh from the triangulation.
262  // The mesh is not parallel consistent - only used for viewing
264  (
265  const fileName& name,
266  labelTolabelPairHashTable& vertexMap,
267  labelList& cellMap,
268  const bool writeDelaunayData = true
269  ) const;
270 };
271 
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
275 } // End namespace Foam
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 #include "DelaunayMeshI.H"
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 #ifdef NoRepository
284  #include "DelaunayMesh.C"
285 #endif
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 #endif
290 
291 // ************************************************************************* //
Foam::DelaunayMesh::DelaunayMesh
DelaunayMesh(const DelaunayMesh< Triangulation > &)
Disallow default bitwise copy construct.
Foam::DelaunayMesh::timeCheck
void timeCheck(const string &description, const bool check=true) const
Write the cpuTime to screen.
Definition: DelaunayMeshI.H:37
Foam::DelaunayMesh::getNewVertexIndex
label getNewVertexIndex() const
Create a new unique vertex index and return.
Definition: DelaunayMeshI.H:78
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::DelaunayMesh::Traits_for_spatial_sort
Spatial sort traits to use with a pair of point pointers and an int.
Definition: DelaunayMesh.H:111
Foam::DelaunayMesh::Traits_for_spatial_sort::less_y_3_object
Less_y_3 less_y_3_object() const
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
DelaunayMeshI.H
Foam::DelaunayMesh::labelPairHashSet
HashSet< Pair< label >, FixedList< label, 2 >::Hash<> > labelPairHashSet
Definition: DelaunayMesh.H:83
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_z_3::operator()
bool operator()(const Point_3 &p, const Point_3 &q) const
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::DelaunayMesh::Traits_for_spatial_sort::less_x_3_object
Less_x_3 less_x_3_object() const
Foam::DelaunayMesh::createMesh
autoPtr< polyMesh > createMesh(const fileName &name, labelTolabelPairHashTable &vertexMap, labelList &cellMap, const bool writeDelaunayData=true) const
Create an fvMesh from the triangulation.
Foam::DelaunayMesh::resetVertexCount
void resetVertexCount()
Set the vertex count to zero.
Definition: DelaunayMeshI.H:114
Foam::DelaunayMesh::time
const Time & time() const
Return a reference to the Time object.
Definition: DelaunayMeshI.H:29
Foam::DelaunayMesh::resetCellCount
void resetCellCount()
Set the cell count to zero.
Definition: DelaunayMeshI.H:107
Foam::Map< label >
Pair.H
Foam::DelaunayMesh::addPatches
void addPatches(const label nInternalFaces, faceList &faces, labelList &owner, PtrList< dictionary > &patchDicts, const List< DynamicList< face > > &patchFaces, const List< DynamicList< label > > &patchOwners) const
Foam::DelaunayMesh::rangeInsertWithInfo
Map< label > rangeInsertWithInfo(PointIterator begin, PointIterator end, bool printErrors=false, bool reIndex=true)
Function inserting points into a triangulation and setting the.
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
Foam::DelaunayMesh::sortFaces
void sortFaces(faceList &faces, labelList &owner, labelList &neighbour) const
Foam::DelaunayMesh::operator=
void operator=(const DelaunayMesh< Triangulation > &)
Disallow default bitwise assignment.
Foam::DelaunayMesh::Traits_for_spatial_sort::less_z_3_object
Less_z_3 less_z_3_object() const
patchDicts
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
Foam::DelaunayMesh::Finite_vertices_iterator
Triangulation::Finite_vertices_iterator Finite_vertices_iterator
Definition: DelaunayMesh.H:73
Foam::DelaunayMesh::vertexCount_
label vertexCount_
Keep track of the number of vertices that have been added.
Definition: DelaunayMesh.H:99
Foam::DelaunayMesh::vertexCount
label vertexCount() const
Return the vertex count (the next unique vertex index)
Definition: DelaunayMeshI.H:100
Foam::DelaunayMesh::getNewCellIndex
label getNewCellIndex() const
Create a new unique cell index and return.
Definition: DelaunayMeshI.H:63
Foam::label
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
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_x_3::operator()
bool operator()(const Point_3 &p, const Point_3 &q) const
Foam::DelaunayMesh::cellCount
label cellCount() const
Return the cell count (the next unique cell index)
Definition: DelaunayMeshI.H:93
Foam::DelaunayMesh::insertPoints
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo)
Foam::DelaunayMesh::cellCount_
label cellCount_
Keep track of the number of cells that have been added.
Definition: DelaunayMesh.H:103
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_x_3
Definition: DelaunayMesh.H:120
Foam::DelaunayMesh::Facet
Triangulation::Facet Facet
Definition: DelaunayMesh.H:70
DelaunayMesh.C
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_y_3
Definition: DelaunayMesh.H:125
HashSet.H
Foam::DelaunayMesh::Traits_for_spatial_sort::Point_3
std::pair< const typename Triangulation::Point *, label > Point_3
Definition: DelaunayMesh.H:118
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::DelaunayMesh::~DelaunayMesh
~DelaunayMesh()
Destructor.
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
boundBox.H
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_z_3
Definition: DelaunayMesh.H:130
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::DelaunayMesh::labelTolabelPairHashTable
HashTable< label, labelPair, FixedList< label, 2 >::Hash<> > labelTolabelPairHashTable
Definition: DelaunayMesh.H:90
Foam::DelaunayMesh::Cell_handle
Triangulation::Cell_handle Cell_handle
Definition: DelaunayMesh.H:66
Foam::DelaunayMesh::Finite_cells_iterator
Triangulation::Finite_cells_iterator Finite_cells_iterator
Definition: DelaunayMesh.H:75
Foam::Pair< label >
Foam::DelaunayMesh::Edge
Triangulation::Edge Edge
Definition: DelaunayMesh.H:68
Foam::DelaunayMesh::Traits_for_spatial_sort::Gt
Triangulation::Geom_traits Gt
Definition: DelaunayMesh.H:115
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_y_3::operator()
bool operator()(const Point_3 &p, const Point_3 &q) const
Foam::DelaunayMesh::printInfo
void printInfo(Ostream &os) const
Write mesh statistics to stream.
Foam::DelaunayMesh::runTime_
const Time & runTime_
Reference to Time.
Definition: DelaunayMesh.H:106
meshName
static char meshName[]
Definition: globalFoam.H:7
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:51
Foam::DelaunayMesh
The vertex and cell classes must have an index defined.
Definition: DelaunayMesh.H:60
Foam::DelaunayMesh::printVertexInfo
void printVertexInfo(Ostream &os) const
Write vertex statistics in the form of a table to stream.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
FixedList.H
CGALTriangulation3Ddefs.H
CGAL data structures used for 3D Delaunay meshing.
Foam::labelPair
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
Foam::DelaunayMesh::Finite_facets_iterator
Triangulation::Finite_facets_iterator Finite_facets_iterator
Definition: DelaunayMesh.H:77
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::DelaunayMesh::Vertex_handle
Triangulation::Vertex_handle Vertex_handle
Definition: DelaunayMesh.H:67
Foam::DelaunayMesh::reset
void reset()
Clear the entire triangulation.
Foam::DelaunayMesh::Point
Triangulation::Point Point
Definition: DelaunayMesh.H:69
autoPtr.H