42template<
class Triangulation>
43Foam::DelaunayMesh<Triangulation>::DelaunayMesh(
const Time&
runTime)
52template<
class Triangulation>
53Foam::DelaunayMesh<Triangulation>::DelaunayMesh
64 Info<<
"Reading " << meshName <<
" from " <<
runTime.timeName() << endl;
72 meshName/polyMesh::meshSubDir,
74 IOobject::READ_IF_PRESENT,
79 if (
pts.typeHeaderOk<pointIOField>(
true))
121 List<Vb> pointsToInsert(
pts.size());
123 forAll(pointsToInsert, pI)
130 static_cast<indexedVertexEnum::vertexType
>(types[pI]),
137 pointsToInsert.begin(),
138 pointsToInsert.end(),
143 vertexCount_ = Triangulation::number_of_vertices();
150template<
class Triangulation>
163template<
class Triangulation>
166 Info<<
"Clearing triangulation" <<
endl;
172 Finite_vertices_iterator vit = Triangulation::finite_vertices_begin();
173 vit != Triangulation::finite_vertices_end();
199 insertPoints(vertices,
false);
201 Info<<
"Inserted " << vertexCount() <<
" fixed points" <<
endl;
205template<
class Triangulation>
208 const List<Vb>& vertices,
212 return rangeInsertWithInfo
222template<
class Triangulation>
230 return typename Gt::Less_x_3()(*(
p.first), *(q.first));
233template<
class Triangulation>
241 return typename Gt::Less_y_3()(*(
p.first), *(q.first));
244template<
class Triangulation>
252 return typename Gt::Less_z_3()(*(
p.first), *(q.first));
255template<
class Triangulation>
257Foam::DelaunayMesh<Triangulation>::Traits_for_spatial_sort::less_x_3_object()
263template<
class Triangulation>
264typename Foam::DelaunayMesh<Triangulation>::Traits_for_spatial_sort::Less_y_3
265Foam::DelaunayMesh<Triangulation>::Traits_for_spatial_sort::less_y_3_object()
271template<
class Triangulation>
272typename Foam::DelaunayMesh<Triangulation>::Traits_for_spatial_sort::Less_z_3
273Foam::DelaunayMesh<Triangulation>::Traits_for_spatial_sort::less_z_3_object()
280template<
class Triangulation>
281template<
class Po
intIterator>
294 const typename Triangulation::Point*,
297 > vectorPairPointIndex;
299 vectorPairPointIndex
points;
302 for (PointIterator it = begin; it !=
end; ++it)
306 std::make_pair(&(it->point()), count++)
310 std::shuffle(
points.begin(),
points.end(), std::default_random_engine());
316 Traits_for_spatial_sort()
321 Map<label> oldToNewIndex(
points.size());
325 typename vectorPairPointIndex::const_iterator
p =
points.begin();
330 const size_t checkInsertion = Triangulation::number_of_vertices();
332 hint = this->
insert(*(
p->first), hint);
334 const Vb& vert = *(begin +
p->second);
336 if (checkInsertion != Triangulation::number_of_vertices() - 1)
341 Triangulation::nearest_vertex(*(
p->first));
343 Pout<<
"Failed insertion : " << vert.
info()
344 <<
" nearest : " << nearV->info();
349 const label oldIndex = vert.
index();
354 oldToNewIndex.insert(oldIndex, hint->index());
357 hint->type() = vert.
type();
364 return oldToNewIndex;
CGAL::indexedVertex< K > Vb
Foam::scalar & targetCellSize()
Foam::InfoProxy< indexedVertex< Gt, Vb > > info() const noexcept
Return info proxy, used to print information to a stream.
Foam::tensor & alignment()
The vertex and cell classes must have an index defined.
Triangulation::Vertex_handle Vertex_handle
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo).
Map< label > rangeInsertWithInfo(PointIterator begin, PointIterator end, bool printErrors=false, bool reIndex=true)
Function inserting points into a triangulation and setting the.
~DelaunayMesh()
Destructor.
label getNewVertexIndex() const
Create a new unique vertex index and return.
void reset()
Clear the entire triangulation.
void append(const T &val)
Append an element at the end of the list.
A HashTable to objects of type <T> with a label key.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
iterator end() noexcept
Return an iterator to end traversing the UList.
T & last()
Access last element of the list, position [size()-1].
A HashTable to objects of type <T> with a labelPair key. The hashing is based on labelPair (FixedList...
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
PointFrompoint toPoint(const Foam::point &p)
messageStream Info
Information stream (stdout output on master, null elsewhere).
vectorIOField pointIOField
pointIOField is a vectorIOField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
pointField vertices(const blockVertexList &bvl)
IOField< label > labelIOField
IO for a Field of label.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
nonInt insert("surfaceSum(((S|magSf)*S)")
#define forAll(list, i)
Loop across all elements in list.