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