44#include "indexedVertex.H"
60template<
class Triangulation>
86 mutable label vertexCount_;
90 mutable label cellCount_;
98 struct Traits_for_spatial_sort
100 public Triangulation::Geom_traits
102 typedef typename Triangulation::Geom_traits Gt;
104 typedef std::pair<const typename Triangulation::Point*, label>
109 bool operator()(
const Point_3&
p,
const Point_3& q)
const;
114 bool operator()(
const Point_3&
p,
const Point_3& q)
const;
119 bool operator()(
const Point_3&
p,
const Point_3& q)
const;
139 const label nInternalFaces,
185 const string& description,
186 const bool check =
true
229 template<
class Po
intIterator>
234 bool printErrors =
false,
254 const bool writeDelaunayData =
true
CGAL data structures used for 3D Delaunay meshing.
The vertex and cell classes must have an index defined.
Triangulation::Vertex_handle Vertex_handle
void resetVertexCount()
Set the vertex count to zero.
Triangulation::Facet Facet
const Time & time() const
Return a reference to the Time object.
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo)
Triangulation::Finite_vertices_iterator Finite_vertices_iterator
Map< label > rangeInsertWithInfo(PointIterator begin, PointIterator end, bool printErrors=false, bool reIndex=true)
Function inserting points into a triangulation and setting the.
label vertexCount() const
Return the vertex count (the next unique vertex index)
label getNewCellIndex() const
Create a new unique cell index and return.
void printVertexInfo(Ostream &os) const
Write vertex statistics in the form of a table to stream.
Triangulation::Finite_facets_iterator Finite_facets_iterator
DelaunayMesh(const Time &runTime, const word &meshName)
Triangulation::Cell_handle Cell_handle
~DelaunayMesh()
Destructor.
autoPtr< polyMesh > createMesh(const fileName &name, labelPairLookup &vertexMap, labelList &cellMap, const bool writeDelaunayData=true) const
Create an fvMesh from the triangulation.
label getNewVertexIndex() const
Create a new unique vertex index and return.
Triangulation::Finite_cells_iterator Finite_cells_iterator
DelaunayMesh(const Time &runTime)
Construct from components.
void reset()
Clear the entire triangulation.
label cellCount() const
Return the cell count (the next unique cell index)
void resetCellCount()
Set the cell count to zero.
void printInfo(Ostream &os) const
Write mesh statistics to stream.
Triangulation::Point Point
void timeCheck(const string &description, const bool check=true) const
Write the cpuTime to screen.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
A HashTable to objects of type <T> with a label key.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
CellSizeDelaunay::Cell_handle Cell_handle
CellSizeDelaunay::Vertex_handle Vertex_handle
CellSizeDelaunay::Point Point
A class for handling file names.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
A HashTable to objects of type <T> with a labelPair key. The hashing is based on labelPair (FixedList...
static void check(const int retVal, const char *what)
pointField vertices(const blockVertexList &bvl)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
PtrList< dictionary > patchDicts
bool operator()(const Point_3 &p, const Point_3 &q) const
bool operator()(const Point_3 &p, const Point_3 &q) const
bool operator()(const Point_3 &p, const Point_3 &q) const