40#ifndef Foam_CGAL_indexedVertex_H
41#define Foam_CGAL_indexedVertex_H
44#include "CGAL/version.h"
45#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1050211000)
46#define BOOST_BIND_GLOBAL_PLACEHOLDERS
51#include "CGAL/Triangulation_3.h"
63template<
class Gt,
class Vb>
class indexedVertex;
72template<
class Gt,
class Vb>
Ostream&
operator<<
78template<
class Gt,
class Vb>
Ostream&
operator<<
84template<
class Gt,
class Vb>
Istream&
operator>>
93 CGAL::Point_3<baseK>&
p
99 const CGAL::Point_3<baseK>&
p
112template<
class Gt,
class Vb = CGAL::Triangulation_vertex_base_3<Gt>>
135 Foam::scalar targetCellSize_;
143 typedef typename Vb::Triangulation_data_structure
Tds;
195 inline Foam::label&
index();
197 inline Foam::label
index()
const;
222 inline bool real()
const;
274 inline bool fixed()
const;
277 inline bool&
fixed();
287 this->type_ = rhs.type();
288 this->index_ = rhs.index();
289 this->processor_ = rhs.procIndex();
290 this->alignment_ = rhs.alignment();
291 this->targetCellSize_ = rhs.targetCellSize();
292 this->vertexFixed_ = rhs.fixed();
300 this->type_ == rhs.type()
301 && this->index_ == rhs.index()
302 && this->processor_ == rhs.procIndex()
303 && this->vertexFixed_ == rhs.fixed()
309 return !(*
this == rhs);
357 CGAL::indexedVertex<K, CGAL::Triangulation_vertex_base_3<K>>
358 > : std::true_type {};
364 > : std::true_type {};
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.
bool internalBaffleEdgePoint() const
bool externalBoundaryPoint() const
bool fixed() const
Is the vertex fixed or movable.
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.
friend Foam::Ostream & Foam::operator(Foam::Ostream &, const Foam::InfoProxy< indexedVertex< Gt, Vb > > &)
bool referred() const
Is this a referred vertex.
bool farPoint() const
Is point a far-point.
Tds::Cell_handle Cell_handle
bool internalOrBoundaryPoint() const
Either original internal point or master of pointPair.
Tds::Vertex_handle Vertex_handle
bool nearOrOnBoundary() const
Is point near the boundary or part of the boundary definition.
bool externalBaffleSurfacePoint() const
Foam::InfoProxy< indexedVertex< Gt, Vb > > info() const
Info proxy, to print information to a stream.
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.
bool internalBoundaryPoint() const
bool surfacePoint() const
Part of a surface point pair.
void operator=(const indexedVertex &rhs)
Copy assignment.
bool internalBaffleSurfacePoint() const
Vb::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.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
OBJstream os(runTime.globalPath()/outputName)
Vb::template Rebind_TDS< TDS2 >::Other Vb2
indexedVertex< Gt, Vb2 > Other