DelaunayMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2015 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "DelaunayMesh.H"
30 #include "polyMesh.H"
31 #include "labelPairHashes.H"
32 #include "PrintTable.H"
33 #include "pointIOField.H"
34 #include "scalarIOField.H"
35 #include "labelIOField.H"
36 #include "pointConversion.H"
37 #include <algorithm>
38 #include <random>
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class Triangulation>
44 :
45  Triangulation(),
46  vertexCount_(0),
47  cellCount_(0),
48  runTime_(runTime)
49 {}
50 
51 
52 template<class Triangulation>
54 (
55  const Time& runTime,
56  const word& meshName
57 )
58 :
59  Triangulation(),
60  vertexCount_(0),
61  cellCount_(0),
62  runTime_(runTime)
63 {
64  Info<< "Reading " << meshName << " from " << runTime.timeName() << endl;
65 
66  pointIOField pts
67  (
68  IOobject
69  (
70  "points",
71  runTime.timeName(),
72  meshName/polyMesh::meshSubDir,
73  runTime,
76  )
77  );
78 
79  if (pts.typeHeaderOk<pointIOField>(true))
80  {
81  labelIOField types
82  (
83  IOobject
84  (
85  "types",
86  runTime.timeName(),
87  meshName,
88  runTime,
91  )
92  );
93 
94 // Do not read in indices
95 // labelIOField indices
96 // (
97 // IOobject
98 // (
99 // "indices",
100 // runTime.timeName(),
101 // meshName,
102 // runTime,
103 // IOobject::MUST_READ,
104 // IOobject::NO_WRITE
105 // )
106 // );
107 
108  labelIOField processorIndices
109  (
110  IOobject
111  (
112  "processorIndices",
113  runTime.timeName(),
114  meshName,
115  runTime,
118  )
119  );
120 
121  List<Vb> pointsToInsert(pts.size());
122 
123  forAll(pointsToInsert, pI)
124  {
125  pointsToInsert[pI] =
126  Vb
127  (
128  toPoint(pts[pI]),
129  pI,
130  static_cast<indexedVertexEnum::vertexType>(types[pI]),
131  processorIndices[pI]
132  );
133  }
134 
135  rangeInsertWithInfo
136  (
137  pointsToInsert.begin(),
138  pointsToInsert.end(),
139  false,
140  false
141  );
142 
143  vertexCount_ = Triangulation::number_of_vertices();
144  }
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
149 
150 template<class Triangulation>
152 {}
153 
154 
155 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
156 
157 
158 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
159 
160 
161 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
162 
163 template<class Triangulation>
165 {
166  Info<< "Clearing triangulation" << endl;
167 
168  DynamicList<Vb> vertices;
169 
170  for
171  (
172  Finite_vertices_iterator vit = Triangulation::finite_vertices_begin();
173  vit != Triangulation::finite_vertices_end();
174  ++vit
175  )
176  {
177  if (vit->fixed())
178  {
179  vertices.append
180  (
181  Vb
182  (
183  vit->point(),
184  vit->index(),
185  vit->type(),
186  vit->procIndex()
187  )
188  );
189 
190  vertices.last().fixed() = vit->fixed();
191  }
192  }
193 
194  this->clear();
195 
196  resetVertexCount();
197  resetCellCount();
198 
199  insertPoints(vertices, false);
200 
201  Info<< "Inserted " << vertexCount() << " fixed points" << endl;
202 }
203 
204 
205 template<class Triangulation>
207 (
208  const List<Vb>& vertices,
209  const bool reIndex
210 )
211 {
212  return rangeInsertWithInfo
213  (
214  vertices.begin(),
215  vertices.end(),
216  false,
217  reIndex
218  );
219 }
220 
221 
222 template<class Triangulation>
224 operator()
225 (
226  const Point_3& p,
227  const Point_3& q
228 ) const
229 {
230  return typename Gt::Less_x_3()(*(p.first), *(q.first));
231 }
232 
233 template<class Triangulation>
235 operator()
236 (
237  const Point_3& p,
238  const Point_3& q
239 ) const
240 {
241  return typename Gt::Less_y_3()(*(p.first), *(q.first));
242 }
243 
244 template<class Triangulation>
246 operator()
247 (
248  const Point_3& p,
249  const Point_3& q
250 ) const
251 {
252  return typename Gt::Less_z_3()(*(p.first), *(q.first));
253 }
254 
255 template<class Triangulation>
258 const
259 {
260  return Less_x_3();
261 }
262 
263 template<class Triangulation>
266 const
267 {
268  return Less_y_3();
269 }
270 
271 template<class Triangulation>
274 const
275 {
276  return Less_z_3();
277 }
278 
279 
280 template<class Triangulation>
281 template<class PointIterator>
283 (
284  PointIterator begin,
285  PointIterator end,
286  bool printErrors,
287  bool reIndex
288 )
289 {
290  typedef DynamicList
291  <
292  std::pair
293  <
294  const typename Triangulation::Point*,
295  label
296  >
297  > vectorPairPointIndex;
298 
299  vectorPairPointIndex points;
300 
301  label count = 0;
302  for (PointIterator it = begin; it != end; ++it)
303  {
304  points.append
305  (
306  std::make_pair(&(it->point()), count++)
307  );
308  }
309 
310  std::shuffle(points.begin(), points.end(), std::default_random_engine());
311 
312  spatial_sort
313  (
314  points.begin(),
315  points.end(),
316  Traits_for_spatial_sort()
317  );
318 
319  Vertex_handle hint;
320 
321  Map<label> oldToNewIndex(points.size());
322 
323  for
324  (
325  typename vectorPairPointIndex::const_iterator p = points.begin();
326  p != points.end();
327  ++p
328  )
329  {
330  const size_t checkInsertion = Triangulation::number_of_vertices();
331 
332  hint = this->insert(*(p->first), hint);
333 
334  const Vb& vert = *(begin + p->second);
335 
336  if (checkInsertion != Triangulation::number_of_vertices() - 1)
337  {
338  if (printErrors)
339  {
340  Vertex_handle nearV =
341  Triangulation::nearest_vertex(*(p->first));
342 
343  Pout<< "Failed insertion : " << vert.info()
344  << " nearest : " << nearV->info();
345  }
346  }
347  else
348  {
349  const label oldIndex = vert.index();
350  hint->index() = getNewVertexIndex();
351 
352  if (reIndex)
353  {
354  oldToNewIndex.insert(oldIndex, hint->index());
355  }
356 
357  hint->type() = vert.type();
358  hint->procIndex() = vert.procIndex();
359  hint->targetCellSize() = vert.targetCellSize();
360  hint->alignment() = vert.alignment();
361  }
362  }
363 
364  return oldToNewIndex;
365 }
366 
367 
368 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369 
370 #include "DelaunayMeshIO.C"
371 
372 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
insert
srcOptions insert("case", fileName(rootDirSource/caseDirSource))
runTime
engineTime & runTime
Definition: createEngineTime.H:13
CGAL::indexedVertex::alignment
Foam::tensor & alignment()
Definition: indexedVertexI.H:189
PrintTable.H
p
volScalarField & p
Definition: createFieldRefs.H:8
pointConversion.H
stdFoam::begin
constexpr auto begin(C &c) -> decltype(c.begin())
Return iterator to the beginning of the container c.
Definition: stdFoam.H:97
Vb
CGAL::indexedVertex< K > Vb
Definition: CGALTriangulation3Ddefs.H:48
CGAL::indexedVertex::type
vertexType & type()
Definition: indexedVertexI.H:174
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:321
Foam::labelIOField
IOField< label > labelIOField
labelField with IO.
Definition: labelIOField.H:44
CGAL::indexedVertex
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
Definition: indexedVertex.H:54
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: lumpedPointController.H:69
pointIOField.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
scalarIOField.H
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.
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
Foam::indexedVertexEnum::vertexType
vertexType
Definition: indexedVertexEnum.H:52
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
labelPairHashes.H
A HashTable to objects of type <T> with a labelPair key. The hashing is based on labelPair (FixedList...
Foam::pointIOField
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:44
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DelaunayMeshIO.C
CGAL::indexedVertex::index
Foam::label & index()
Definition: indexedVertexI.H:159
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::DelaunayMesh::insertPoints
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo)
CGAL::indexedVertex::info
Foam::InfoProxy< indexedVertex< Gt, Vb > > info() const
Info proxy, to print information to a stream.
Definition: indexedVertex.H:307
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_x_3
Definition: DelaunayMesh.H:106
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_y_3
Definition: DelaunayMesh.H:111
CGAL::indexedVertex::targetCellSize
Foam::scalar & targetCellSize()
Definition: indexedVertexI.H:203
Foam::vertices
pointField vertices(const blockVertexList &bvl)
Definition: blockVertexList.H:49
Foam::toPoint
PointFrompoint toPoint(const Foam::point &p)
Definition: pointConversion.H:82
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::DelaunayMesh::~DelaunayMesh
~DelaunayMesh()
Destructor.
Foam::DelaunayMesh::Traits_for_spatial_sort::Less_z_3
Definition: DelaunayMesh.H:116
CGAL::indexedVertex::procIndex
int procIndex() const
Definition: indexedVertexI.H:251
clear
patchWriters clear()
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::shuffle
void shuffle(UList< T > &a)
Definition: UList.C:289
points
const pointField & points
Definition: gmvOutputHeader.H:1
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:53
Foam::DelaunayMesh
The vertex and cell classes must have an index defined.
Definition: DelaunayMesh.H:60
DelaunayMesh.H
Foam::DelaunayMesh::reset
void reset()
Clear the entire triangulation.
labelIOField.H
Foam::IOobject::MUST_READ
Definition: IOobject.H:185