/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
\\/ 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 "triSurfaceTools.H"
#include "triSurface.H"
#include "OFstream.H"
#include "mergePoints.H"
#include "polyMesh.H"
#include "plane.H"
#include "geompack.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
const Foam::label Foam::triSurfaceTools::ANYEDGE = -1;
const Foam::label Foam::triSurfaceTools::NOEDGE = -2;
const Foam::label Foam::triSurfaceTools::COLLAPSED = -3;
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
/*
Refine by splitting all three edges of triangle ('red' refinement).
Neighbouring triangles (which are not marked for refinement get split
in half ('green') refinement. (R. Verfuerth, "A review of a posteriori
error estimation and adaptive mesh refinement techniques",
Wiley-Teubner, 1996)
*/
// FaceI gets refined ('red'). Propagate edge refinement.
void Foam::triSurfaceTools::calcRefineStatus
(
const triSurface& surf,
const label faceI,
List& refine
)
{
if (refine[faceI] == RED)
{
// Already marked for refinement. Do nothing.
}
else
{
// Not marked or marked for 'green' refinement. Refine.
refine[faceI] = RED;
const labelList& myNeighbours = surf.faceFaces()[faceI];
forAll(myNeighbours, myNeighbourI)
{
label neighbourFaceI = myNeighbours[myNeighbourI];
if (refine[neighbourFaceI] == GREEN)
{
// Change to red refinement and propagate
calcRefineStatus(surf, neighbourFaceI, refine);
}
else if (refine[neighbourFaceI] == NONE)
{
refine[neighbourFaceI] = GREEN;
}
}
}
}
// Split faceI along edgeI at position newPointI
void Foam::triSurfaceTools::greenRefine
(
const triSurface& surf,
const label faceI,
const label edgeI,
const label newPointI,
DynamicList& newFaces
)
{
const labelledTri& f = surf.localFaces()[faceI];
const edge& e = surf.edges()[edgeI];
// Find index of edge in face.
label fp0 = findIndex(f, e[0]);
label fp1 = f.fcIndex(fp0);
label fp2 = f.fcIndex(fp1);
if (f[fp1] == e[1])
{
// Edge oriented like face
newFaces.append
(
labelledTri
(
f[fp0],
newPointI,
f[fp2],
f.region()
)
);
newFaces.append
(
labelledTri
(
newPointI,
f[fp1],
f[fp2],
f.region()
)
);
}
else
{
newFaces.append
(
labelledTri
(
f[fp2],
newPointI,
f[fp1],
f.region()
)
);
newFaces.append
(
labelledTri
(
newPointI,
f[fp0],
f[fp1],
f.region()
)
);
}
}
// Refine all triangles marked for refinement.
Foam::triSurface Foam::triSurfaceTools::doRefine
(
const triSurface& surf,
const List& refineStatus
)
{
// Storage for new points. (start after old points)
DynamicList newPoints(surf.nPoints());
forAll(surf.localPoints(), pointI)
{
newPoints.append(surf.localPoints()[pointI]);
}
label newVertI = surf.nPoints();
// Storage for new faces
DynamicList newFaces(surf.size());
// Point index for midpoint on edge
labelList edgeMid(surf.nEdges(), -1);
forAll(refineStatus, faceI)
{
if (refineStatus[faceI] == RED)
{
// Create new vertices on all edges to be refined.
const labelList& fEdges = surf.faceEdges()[faceI];
forAll(fEdges, i)
{
label edgeI = fEdges[i];
if (edgeMid[edgeI] == -1)
{
const edge& e = surf.edges()[edgeI];
// Create new point on mid of edge
newPoints.append
(
0.5
* (
surf.localPoints()[e.start()]
+ surf.localPoints()[e.end()]
)
);
edgeMid[edgeI] = newVertI++;
}
}
// Now we have new mid edge vertices for all edges on face
// so create triangles for RED rerfinement.
const edgeList& edges = surf.edges();
// Corner triangles
newFaces.append
(
labelledTri
(
edgeMid[fEdges[0]],
edges[fEdges[0]].commonVertex(edges[fEdges[1]]),
edgeMid[fEdges[1]],
surf[faceI].region()
)
);
newFaces.append
(
labelledTri
(
edgeMid[fEdges[1]],
edges[fEdges[1]].commonVertex(edges[fEdges[2]]),
edgeMid[fEdges[2]],
surf[faceI].region()
)
);
newFaces.append
(
labelledTri
(
edgeMid[fEdges[2]],
edges[fEdges[2]].commonVertex(edges[fEdges[0]]),
edgeMid[fEdges[0]],
surf[faceI].region()
)
);
// Inner triangle
newFaces.append
(
labelledTri
(
edgeMid[fEdges[0]],
edgeMid[fEdges[1]],
edgeMid[fEdges[2]],
surf[faceI].region()
)
);
// Create triangles for GREEN refinement.
forAll(fEdges, i)
{
const label edgeI = fEdges[i];
label otherFaceI = otherFace(surf, faceI, edgeI);
if ((otherFaceI != -1) && (refineStatus[otherFaceI] == GREEN))
{
greenRefine
(
surf,
otherFaceI,
edgeI,
edgeMid[edgeI],
newFaces
);
}
}
}
}
// Copy unmarked triangles since keep original vertex numbering.
forAll(refineStatus, faceI)
{
if (refineStatus[faceI] == NONE)
{
newFaces.append(surf.localFaces()[faceI]);
}
}
newFaces.shrink();
newPoints.shrink();
// Transfer DynamicLists to straight ones.
pointField allPoints;
allPoints.transfer(newPoints);
newPoints.clear();
return triSurface(newFaces, surf.patches(), allPoints, true);
}
// Angle between two neighbouring triangles,
// angle between shared-edge end points and left and right face end points
Foam::scalar Foam::triSurfaceTools::faceCosAngle
(
const point& pStart,
const point& pEnd,
const point& pLeft,
const point& pRight
)
{
const vector common(pEnd - pStart);
const vector base0(pLeft - pStart);
const vector base1(pRight - pStart);
vector n0(common ^ base0);
n0 /= Foam::mag(n0);
vector n1(base1 ^ common);
n1 /= Foam::mag(n1);
return n0 & n1;
}
// Protect edges around vertex from collapsing.
// Moving vertI to a different position
// - affects obviously the cost of the faces using it
// - but also their neighbours since the edge between the faces changes
void Foam::triSurfaceTools::protectNeighbours
(
const triSurface& surf,
const label vertI,
labelList& faceStatus
)
{
// const labelList& myFaces = surf.pointFaces()[vertI];
// forAll(myFaces, i)
// {
// label faceI = myFaces[i];
//
// if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
// {
// faceStatus[faceI] = NOEDGE;
// }
// }
const labelList& myEdges = surf.pointEdges()[vertI];
forAll(myEdges, i)
{
const labelList& myFaces = surf.edgeFaces()[myEdges[i]];
forAll(myFaces, myFaceI)
{
label faceI = myFaces[myFaceI];
if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
{
faceStatus[faceI] = NOEDGE;
}
}
}
}
//
// Edge collapse helper functions
//
// Get all faces that will get collapsed if edgeI collapses.
Foam::labelHashSet Foam::triSurfaceTools::getCollapsedFaces
(
const triSurface& surf,
label edgeI
)
{
const edge& e = surf.edges()[edgeI];
label v1 = e.start();
label v2 = e.end();
// Faces using edge will certainly get collapsed.
const labelList& myFaces = surf.edgeFaces()[edgeI];
labelHashSet facesToBeCollapsed(2*myFaces.size());
forAll(myFaces, myFaceI)
{
facesToBeCollapsed.insert(myFaces[myFaceI]);
}
// From faces using v1 check if they share an edge with faces
// using v2.
// - share edge: are part of 'splay' tree and will collapse if edge
// collapses
const labelList& v1Faces = surf.pointFaces()[v1];
forAll(v1Faces, v1FaceI)
{
label face1I = v1Faces[v1FaceI];
label otherEdgeI = oppositeEdge(surf, face1I, v1);
// Step across edge to other face
label face2I = otherFace(surf, face1I, otherEdgeI);
if (face2I != -1)
{
// found face on other side of edge. Now check if uses v2.
if (oppositeVertex(surf, face2I, otherEdgeI) == v2)
{
// triangles face1I and face2I form splay tree and will
// collapse.
facesToBeCollapsed.insert(face1I);
facesToBeCollapsed.insert(face2I);
}
}
}
return facesToBeCollapsed;
}
// Return value of faceUsed for faces using vertI
Foam::label Foam::triSurfaceTools::vertexUsesFace
(
const triSurface& surf,
const labelHashSet& faceUsed,
const label vertI
)
{
const labelList& myFaces = surf.pointFaces()[vertI];
forAll(myFaces, myFaceI)
{
label face1I = myFaces[myFaceI];
if (faceUsed.found(face1I))
{
return face1I;
}
}
return -1;
}
// Calculate new edge-face addressing resulting from edge collapse
void Foam::triSurfaceTools::getMergedEdges
(
const triSurface& surf,
const label edgeI,
const labelHashSet& collapsedFaces,
HashTable