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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "DelaunayMesh.H"
29 #include "polyMesh.H"
30 #include "labelPairHashes.H"
31 #include "PrintTable.H"
32 #include "pointIOField.H"
33 #include "scalarIOField.H"
34 #include "labelIOField.H"
35 #include "pointConversion.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class Triangulation>
41 :
42  Triangulation(),
43  vertexCount_(0),
44  cellCount_(0),
45  runTime_(runTime)
46 {}
47 
48 
49 template<class Triangulation>
51 (
52  const Time& runTime,
53  const word& meshName
54 )
55 :
56  Triangulation(),
57  vertexCount_(0),
58  cellCount_(0),
59  runTime_(runTime)
60 {
61  Info<< "Reading " << meshName << " from " << runTime.timeName() << endl;
62 
63  pointIOField pts
64  (
65  IOobject
66  (
67  "points",
68  runTime.timeName(),
69  meshName/polyMesh::meshSubDir,
70  runTime,
73  )
74  );
75 
76  if (pts.typeHeaderOk<pointIOField>(true))
77  {
78  labelIOField types
79  (
80  IOobject
81  (
82  "types",
83  runTime.timeName(),
84  meshName,
85  runTime,
88  )
89  );
90 
91 // Do not read in indices
92 // labelIOField indices
93 // (
94 // IOobject
95 // (
96 // "indices",
97 // runTime.timeName(),
98 // meshName,
99 // runTime,
100 // IOobject::MUST_READ,
101 // IOobject::NO_WRITE
102 // )
103 // );
104 
105  labelIOField processorIndices
106  (
107  IOobject
108  (
109  "processorIndices",
110  runTime.timeName(),
111  meshName,
112  runTime,
115  )
116  );
117 
118  List<Vb> pointsToInsert(pts.size());
119 
120  forAll(pointsToInsert, pI)
121  {
122  pointsToInsert[pI] =
123  Vb
124  (
125  toPoint(pts[pI]),
126  pI,
127  static_cast<indexedVertexEnum::vertexType>(types[pI]),
128  processorIndices[pI]
129  );
130  }
131 
132  rangeInsertWithInfo
133  (
134  pointsToInsert.begin(),
135  pointsToInsert.end(),
136  false,
137  false
138  );
139 
140  vertexCount_ = Triangulation::number_of_vertices();
141  }
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146 
147 template<class Triangulation>
149 {}
150 
151 
152 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
153 
154 
155 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
156 
157 
158 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
159 
160 template<class Triangulation>
162 {
163  Info<< "Clearing triangulation" << endl;
164 
165  DynamicList<Vb> vertices;
166 
167  for
168  (
169  Finite_vertices_iterator vit = Triangulation::finite_vertices_begin();
170  vit != Triangulation::finite_vertices_end();
171  ++vit
172  )
173  {
174  if (vit->fixed())
175  {
176  vertices.append
177  (
178  Vb
179  (
180  vit->point(),
181  vit->index(),
182  vit->type(),
183  vit->procIndex()
184  )
185  );
186 
187  vertices.last().fixed() = vit->fixed();
188  }
189  }
190 
191  this->clear();
192 
193  resetVertexCount();
194  resetCellCount();
195 
196  insertPoints(vertices, false);
197 
198  Info<< "Inserted " << vertexCount() << " fixed points" << endl;
199 }
200 
201 
202 template<class Triangulation>
204 (
205  const List<Vb>& vertices,
206  const bool reIndex
207 )
208 {
209  return rangeInsertWithInfo
210  (
211  vertices.begin(),
212  vertices.end(),
213  false,
214  reIndex
215  );
216 }
217 
218 
219 template<class Triangulation>
221 operator()
222 (
223  const Point_3& p,
224  const Point_3& q
225 ) const
226 {
227  return typename Gt::Less_x_3()(*(p.first), *(q.first));
228 }
229 
230 template<class Triangulation>
232 operator()
233 (
234  const Point_3& p,
235  const Point_3& q
236 ) const
237 {
238  return typename Gt::Less_y_3()(*(p.first), *(q.first));
239 }
240 
241 template<class Triangulation>
243 operator()
244 (
245  const Point_3& p,
246  const Point_3& q
247 ) const
248 {
249  return typename Gt::Less_z_3()(*(p.first), *(q.first));
250 }
251 
252 template<class Triangulation>
255 const
256 {
257  return Less_x_3();
258 }
259 
260 template<class Triangulation>
263 const
264 {
265  return Less_y_3();
266 }
267 
268 template<class Triangulation>
271 const
272 {
273  return Less_z_3();
274 }
275 
276 
277 template<class Triangulation>
278 template<class PointIterator>
280 (
281  PointIterator begin,
282  PointIterator end,
283  bool printErrors,
284  bool reIndex
285 )
286 {
287  typedef DynamicList
288  <
289  std::pair
290  <
291  const typename Triangulation::Point*,
292  label
293  >
294  > vectorPairPointIndex;
295 
296  vectorPairPointIndex points;
297 
298  label count = 0;
299  for (PointIterator it = begin; it != end; ++it)
300  {
301  points.append
302  (
303  std::make_pair(&(it->point()), count++)
304  );
305  }
306 
307  std::random_shuffle(points.begin(), points.end());
308 
309  spatial_sort
310  (
311  points.begin(),
312  points.end(),
313  Traits_for_spatial_sort()
314  );
315 
316  Vertex_handle hint;
317 
318  Map<label> oldToNewIndex(points.size());
319 
320  for
321  (
322  typename vectorPairPointIndex::const_iterator p = points.begin();
323  p != points.end();
324  ++p
325  )
326  {
327  const size_t checkInsertion = Triangulation::number_of_vertices();
328 
329  hint = this->insert(*(p->first), hint);
330 
331  const Vb& vert = *(begin + p->second);
332 
333  if (checkInsertion != Triangulation::number_of_vertices() - 1)
334  {
335  if (printErrors)
336  {
337  Vertex_handle nearV =
338  Triangulation::nearest_vertex(*(p->first));
339 
340  Pout<< "Failed insertion : " << vert.info()
341  << " nearest : " << nearV->info();
342  }
343  }
344  else
345  {
346  const label oldIndex = vert.index();
347  hint->index() = getNewVertexIndex();
348 
349  if (reIndex)
350  {
351  oldToNewIndex.insert(oldIndex, hint->index());
352  }
353 
354  hint->type() = vert.type();
355  hint->procIndex() = vert.procIndex();
356  hint->targetCellSize() = vert.targetCellSize();
357  hint->alignment() = vert.alignment();
358  }
359  }
360 
361  return oldToNewIndex;
362 }
363 
364 
365 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
366 
367 #include "DelaunayMeshIO.C"
368 
369 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
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:350
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 (uses stdout - output is on the master only)
DelaunayMeshIO.C
CGAL::indexedVertex::index
Foam::label & index()
Definition: indexedVertexI.H:159
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
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
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:120