indexedVertex.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 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 Class
28  CGAL::indexedVertex
29 
30 Description
31  An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep
32  track of the Delaunay vertices in the tessellation.
33 
34 SourceFiles
35  indexedVertexI.H
36  indexedVertex.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef Foam_CGAL_indexedVertex_H
41 #define Foam_CGAL_indexedVertex_H
42 
43 // Silence boost bind deprecation warnings (before CGAL-5.2.1)
44 #include "CGAL/version.h"
45 #if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1050211000)
46 #define BOOST_BIND_GLOBAL_PLACEHOLDERS
47 #endif
48 #pragma clang diagnostic ignored "-Wbitwise-instead-of-logical"
49 
50 // ------------------------------------------------------------------------- //
51 
52 #include "CGAL/Triangulation_3.h"
54 #include "tensor.H"
55 #include "triad.H"
56 #include "InfoProxy.H"
57 #include "point.H"
58 #include "indexedVertexEnum.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace CGAL
63 {
64 template<class Gt, class Vb> class indexedVertex;
65 }
66 
67 namespace Foam
68 {
69 
70 class Ostream;
71 class Istream;
72 
73 template<class Gt, class Vb> Ostream& operator<<
74 (
75  Ostream&,
77 );
78 
79 template<class Gt, class Vb> Ostream& operator<<
80 (
81  Ostream&,
83 );
84 
85 template<class Gt, class Vb> Istream& operator>>
86 (
87  Istream&,
89 );
90 
91 inline Istream& operator>>
92 (
93  Istream& is,
94  CGAL::Point_3<baseK>& p
95 );
96 
97 inline Ostream& operator<<
98 (
99  Ostream& os,
100  const CGAL::Point_3<baseK>& p
101 );
102 
103 } // End namespace Foam
104 
105 
106 namespace CGAL
107 {
108 
109 /*---------------------------------------------------------------------------*\
110  Class indexedVertex Declaration
111 \*---------------------------------------------------------------------------*/
112 
113 template<class Gt, class Vb = CGAL::Triangulation_vertex_base_3<Gt>>
114 class indexedVertex
115 :
117  public Vb
118 {
119  // Private Data
120 
121  //- Type of pair-point
122  vertexType type_;
123 
124  //- The index for this Delaunay vertex. For referred vertices, the
125  // index is negative for vertices that are the outer (slave) of point
126  // pairs
127  Foam::label index_;
128 
129  //- Number of the processor that owns this vertex
130  int processor_;
131 
132  //- Required alignment of the dual cell of this vertex
133  Foam::tensor alignment_;
134 
135  //- Target size of the dual cell of this vertex
136  Foam::scalar targetCellSize_;
137 
138  //- Specify whether the vertex is fixed or movable.
139  bool vertexFixed_;
140 
141 
142 public:
143 
144  typedef typename Vb::Triangulation_data_structure Tds;
145  typedef typename Vb::Point Point;
146  typedef typename Tds::Vertex_handle Vertex_handle;
147  typedef typename Tds::Cell_handle Cell_handle;
148 
149  template<class TDS2>
150  struct Rebind_TDS
151  {
152  typedef typename Vb::template Rebind_TDS<TDS2>::Other Vb2;
154  };
156 
157  // Generated Methods
159  //- Copy construct
160  indexedVertex(const indexedVertex&) = default;
162 
163  // Constructors
164 
165  inline indexedVertex();
166 
167  inline indexedVertex(const Point& p);
168 
169  inline indexedVertex(const Point& p, vertexType type);
170 
171  inline indexedVertex(const Foam::point& p, vertexType type);
172 
173  inline indexedVertex
174  (
175  const Point& p,
176  Foam::label index,
178  int processor
179  );
180 
181  inline indexedVertex
182  (
183  const Foam::point& p,
184  Foam::label index,
186  int processor
187  );
188 
189  inline indexedVertex(const Point& p, Cell_handle f);
190 
191  inline indexedVertex(Cell_handle f);
192 
193 
194  // Member Functions
195 
196  inline Foam::label& index();
197 
198  inline Foam::label index() const;
199 
200  inline vertexType& type();
201 
202  inline vertexType type() const;
203 
204  inline Foam::tensor& alignment();
205 
206  inline const Foam::tensor& alignment() const;
207 
208  inline Foam::scalar& targetCellSize();
209 
210  inline Foam::scalar targetCellSize() const;
211 
212  //- Is point a far-point
213  inline bool farPoint() const;
214 
215  //- Is point internal, i.e. not on boundary
216  inline bool internalPoint() const;
217 
218  //- Is this a referred vertex
219  inline bool referred() const;
220 
221  //- Is this a "real" point on this processor, i.e. is internal or part
222  // of the boundary description, and not a "far" or "referred" point
223  inline bool real() const;
224 
225  // For referred vertices, what is the original processor index
226  inline int procIndex() const;
227 
228  // For referred vertices, set the original processor index
229  inline int& procIndex();
230 
231  //- Set the point to be internal
232  inline void setInternal();
233 
234  //- Is point internal and near the boundary
235  inline bool nearBoundary() const;
236 
237  //- Set the point to be near the boundary
238  inline void setNearBoundary();
239 
240  //- Is point internal and near a proc boundary
241  inline bool nearProcBoundary() const;
242 
243  //- Set the point to be near a proc boundary
244  inline void setNearProcBoundary();
245 
246  //- Either master or slave of pointPair.
247  inline bool boundaryPoint() const;
248 
249  //- Either original internal point or master of pointPair.
250  inline bool internalOrBoundaryPoint() const;
251 
252  //- Is point near the boundary or part of the boundary definition
253  inline bool nearOrOnBoundary() const;
254 
255  //- Part of a feature point
256  inline bool featurePoint() const;
257 
258  //- Part of a feature edge
259  inline bool featureEdgePoint() const;
260 
261  //- Part of a surface point pair
262  inline bool surfacePoint() const;
263 
264  inline bool internalBoundaryPoint() const;
265  inline bool internalBaffleSurfacePoint() const;
266  inline bool internalBaffleEdgePoint() const;
267 
268  inline bool externalBoundaryPoint() const;
269  inline bool externalBaffleSurfacePoint() const;
270  inline bool externalBaffleEdgePoint() const;
271 
272  inline bool constrained() const;
273 
274  //- Is the vertex fixed or movable
275  inline bool fixed() const;
276 
277  //- Fix the vertex so that it can't be moved
278  inline bool& fixed();
279 
280 
281  // Member Operators
282 
283  //- Copy assignment
284  inline void operator=(const indexedVertex& rhs)
285  {
286  Vb::operator=(rhs);
287 
288  this->type_ = rhs.type();
289  this->index_ = rhs.index();
290  this->processor_ = rhs.procIndex();
291  this->alignment_ = rhs.alignment();
292  this->targetCellSize_ = rhs.targetCellSize();
293  this->vertexFixed_ = rhs.fixed();
294  }
295 
296  inline bool operator==(const indexedVertex& rhs) const
297  {
298  return
299  (
300  //this->point() == rhs.point()
301  this->type_ == rhs.type()
302  && this->index_ == rhs.index()
303  && this->processor_ == rhs.procIndex()
304  && this->vertexFixed_ == rhs.fixed()
305  );
306  }
307 
308  inline bool operator!=(const indexedVertex& rhs) const
309  {
310  return !(*this == rhs);
311  }
312 
313 
314  // IOstream Operators
315 
316  //- Return info proxy,
317  //- used to print information to a stream
319  {
320  return *this;
321  }
322 
323  friend Foam::Ostream& Foam::operator<< <Gt, Vb>
324  (
325  Foam::Ostream&,
327  );
328 
329  friend Foam::Ostream& Foam::operator<< <Gt, Vb>
330  (
332  const indexedVertex<Gt, Vb>&
333  );
334 
335  friend Foam::Istream& Foam::operator>> <Gt, Vb>
336  (
337  Foam::Istream&,
339  );
340 };
341 
342 
343 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344 
345 } // End namespace CGAL
346 
347 
348 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
349 
350 #ifdef CGAL_INEXACT
351 namespace Foam
352 {
353  // For inexact representations where the storage type is a double, the data
354  // is contiguous. This may not be true for exact number types.
356  template<>
357  struct is_contiguous
358  <
359  CGAL::indexedVertex<K, CGAL::Triangulation_vertex_base_3<K>>
360  > : std::true_type {};
361 
362  template<>
363  struct is_contiguous
364  <
365  CGAL::Triangulation_vertex_base_3<K>::Point
366  > : std::true_type {};
368 } // End namespace Foam
369 #endif
370 
371 
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373 
374 #include "indexedVertexI.H"
375 
376 #ifdef NoRepository
377  #include "indexedVertex.C"
378 #endif
379 
380 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
381 
382 #endif
383 
384 // ************************************************************************* //
void setNearBoundary()
Set the point to be near the boundary.
bool boundaryPoint() const
Either master or slave of pointPair.
bool farPoint() const
Is point a far-point.
Vb::Triangulation_data_structure Tds
bool internalBaffleSurfacePoint() const
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
vertexType & type()
bool internalPoint() const
Is point internal, i.e. not on boundary.
Foam::label & index()
bool internalBaffleEdgePoint() const
void setInternal()
Set the point to be internal.
Tds::Cell_handle Cell_handle
Foam::scalar & targetCellSize()
void operator=(const indexedVertex &rhs)
Copy assignment.
bool nearOrOnBoundary() const
Is point near the boundary or part of the boundary definition.
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
Definition: indexedVertex.H:59
Specializations for CGAL-related routines.
Definition: indexedCell.H:60
bool externalBoundaryPoint() const
bool constrained() const
bool surfacePoint() const
Part of a surface point pair.
Foam::InfoProxy< indexedVertex< Gt, Vb > > info() const noexcept
Return info proxy, used to print information to a stream.
bool internalBoundaryPoint() const
bool real() const
Is this a "real" point on this processor, i.e. is internal or part.
friend Foam::Ostream & Foam::operator(Foam::Ostream &, const Foam::InfoProxy< indexedVertex< Gt, Vb >> &)
bool featureEdgePoint() const
Part of a feature edge.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
bool operator!=(const indexedVertex &rhs) const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
bool operator==(const indexedVertex &rhs) const
OBJstream os(runTime.globalPath()/outputName)
labelList f(nPoints)
Foam::tensor & alignment()
CGAL::Point_3< K > Point
A helper class for outputting values to Ostream.
Definition: ensightCells.H:43
bool featurePoint() const
Part of a feature point.
bool externalBaffleEdgePoint() const
bool referred() const
Is this a referred vertex.
void setNearProcBoundary()
Set the point to be near a proc boundary.
bool fixed() const
Is the vertex fixed or movable.
Tds::Vertex_handle Vertex_handle
indexedVertex< Gt, Vb2 > Other
volScalarField & p
bool internalOrBoundaryPoint() const
Either original internal point or master of pointPair.
Tensor of scalars, i.e. Tensor<scalar>.
bool externalBaffleSurfacePoint() const
Vb::template Rebind_TDS< TDS2 >::Other Vb2
bool nearProcBoundary() const
Is point internal and near a proc boundary.
Namespace for OpenFOAM.
bool nearBoundary() const
Is point internal and near the boundary.