/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see .
\*---------------------------------------------------------------------------*/
#include "CV2D.H"
#include "OFstream.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::CV2D::writePoints(const fileName& fName, bool internalOnly) const
{
Info<< "Writing points to " << fName << nl << endl;
OFstream str(fName);
for
(
Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
vit != finite_vertices_end();
++vit
)
{
if (!internalOnly || vit->internalOrBoundaryPoint())
{
meshTools::writeOBJ(str, toPoint3D(vit->point()));
}
}
}
void Foam::CV2D::writeTriangles(const fileName& fName, bool internalOnly) const
{
Info<< "Writing triangles to " << fName << nl << endl;
OFstream str(fName);
labelList vertexMap(number_of_vertices(), -2);
label verti = 0;
for
(
Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
vit != finite_vertices_end();
++vit
)
{
if (!internalOnly || !vit->farPoint())
{
vertexMap[vit->index()] = verti++;
meshTools::writeOBJ(str, toPoint3D(vit->point()));
}
}
for
(
Triangulation::Finite_faces_iterator fit = finite_faces_begin();
fit != finite_faces_end();
++fit
)
{
if
(
!internalOnly
|| (
fit->vertex(0)->internalOrBoundaryPoint()
|| fit->vertex(1)->internalOrBoundaryPoint()
|| fit->vertex(2)->internalOrBoundaryPoint()
)
)
{
str << "f";
for (label i = 0; i < 3; ++i)
{
str << " " << vertexMap[fit->vertex(i)->index()] + 1;
}
str << nl;
}
}
}
void Foam::CV2D::writeFaces(const fileName& fName, bool internalOnly) const
{
Info<< "Writing dual faces to " << fName << nl << endl;
OFstream str(fName);
label dualVerti = 0;
for
(
Triangulation::Finite_faces_iterator fit = finite_faces_begin();
fit != finite_faces_end();
++fit
)
{
if
(
!internalOnly
|| (
fit->vertex(0)->internalOrBoundaryPoint()
|| fit->vertex(1)->internalOrBoundaryPoint()
|| fit->vertex(2)->internalOrBoundaryPoint()
)
)
{
fit->faceIndex() = dualVerti++;
meshTools::writeOBJ(str, toPoint3D(circumcenter(fit)));
}
else
{
fit->faceIndex() = -1;
}
}
for
(
Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
vit != finite_vertices_end();
++vit
)
{
if (!internalOnly || vit->internalOrBoundaryPoint())
{
Face_circulator fcStart = incident_faces(vit);
Face_circulator fc = fcStart;
str<< 'f';
do
{
if (!is_infinite(fc))
{
if (fc->faceIndex() < 0)
{
FatalErrorIn
(
"Foam::CV2D::writeFaces"
"(const fileName& fName, bool internalOnly)"
)<< "Dual face uses vertex defined by a triangle"
" defined by an external point"
<< exit(FatalError);
}
str<< ' ' << fc->faceIndex() + 1;
}
} while (++fc != fcStart);
str<< nl;
}
}
}
void Foam::CV2D::calcDual
(
point2DField& dualPoints,
faceList& dualFaces,
wordList& patchNames,
labelList& patchSizes,
EdgeMap