/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2016 OpenFOAM Foundation
Copyright (C) 2020-2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "isoSurfaceCell.H"
#include "isoSurfacePoint.H"
#include "polyMesh.H"
#include "tetMatcher.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template
Type Foam::isoSurfaceCell::generatePoint
(
const UList& snappedPoints,
const scalar s0,
const Type& p0,
const label p0Index,
const scalar s1,
const Type& p1,
const label p1Index
) const
{
const scalar d = s1-s0;
if (mag(d) > VSMALL)
{
const scalar s = (iso_-s0)/d;
if (s >= 0.5 && s <= 1 && p1Index != -1)
{
return snappedPoints[p1Index];
}
else if (s >= 0.0 && s <= 0.5 && p0Index != -1)
{
return snappedPoints[p0Index];
}
else
{
return s*p1 + (1.0-s)*p0;
}
}
else
{
constexpr scalar s = 0.4999;
return s*p1 + (1.0-s)*p0;
}
}
template
void Foam::isoSurfaceCell::generateTriPoints
(
const UList& snapped,
const scalar s0,
const Type& p0,
const label p0Index,
const scalar s1,
const Type& p1,
const label p1Index,
const scalar s2,
const Type& p2,
const label p2Index,
const scalar s3,
const Type& p3,
const label p3Index,
DynamicList& pts
) const
{
int triIndex = 0;
if (s0 < iso_)
{
triIndex |= 1;
}
if (s1 < iso_)
{
triIndex |= 2;
}
if (s2 < iso_)
{
triIndex |= 4;
}
if (s3 < iso_)
{
triIndex |= 8;
}
// Form the vertices of the triangles for each case
switch (triIndex)
{
case 0x00:
case 0x0F:
break;
case 0x01:
case 0x0E:
{
pts.append(generatePoint(snapped,s0,p0,p0Index,s1,p1,p1Index));
pts.append(generatePoint(snapped,s0,p0,p0Index,s2,p2,p2Index));
pts.append(generatePoint(snapped,s0,p0,p0Index,s3,p3,p3Index));
if (triIndex == 0x0E)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
case 0x02:
case 0x0D:
{
pts.append(generatePoint(snapped,s1,p1,p1Index,s0,p0,p0Index));
pts.append(generatePoint(snapped,s1,p1,p1Index,s3,p3,p3Index));
pts.append(generatePoint(snapped,s1,p1,p1Index,s2,p2,p2Index));
if (triIndex == 0x0D)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
case 0x03:
case 0x0C:
{
Type p0p2 = generatePoint(snapped,s0,p0,p0Index,s2,p2,p2Index);
Type p1p3 = generatePoint(snapped,s1,p1,p1Index,s3,p3,p3Index);
pts.append(generatePoint(snapped,s0,p0,p0Index,s3,p3,p3Index));
pts.append(p1p3);
pts.append(p0p2);
pts.append(p1p3);
pts.append
(
generatePoint(snapped,s1,p1,p1Index,s2,p2,p2Index)
);
pts.append(p0p2);
if (triIndex == 0x0C)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-5], pts[sz-4]);
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
case 0x04:
case 0x0B:
{
pts.append(generatePoint(snapped,s2,p2,p2Index,s0,p0,p0Index));
pts.append(generatePoint(snapped,s2,p2,p2Index,s1,p1,p1Index));
pts.append(generatePoint(snapped,s2,p2,p2Index,s3,p3,p3Index));
if (triIndex == 0x0B)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
case 0x05:
case 0x0A:
{
Type p0p1 = generatePoint(snapped,s0,p0,p0Index,s1,p1,p1Index);
Type p2p3 = generatePoint(snapped,s2,p2,p2Index,s3,p3,p3Index);
pts.append(p0p1);
pts.append(p2p3);
pts.append(generatePoint(snapped,s0,p0,p0Index,s3,p3,p3Index));
pts.append(p0p1);
pts.append(generatePoint(snapped,s1,p1,p1Index,s2,p2,p2Index));
pts.append(p2p3);
if (triIndex == 0x0A)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-5], pts[sz-4]);
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
case 0x06:
case 0x09:
{
Type p0p1 = generatePoint(snapped,s0,p0,p0Index,s1,p1,p1Index);
Type p2p3 = generatePoint(snapped,s2,p2,p2Index,s3,p3,p3Index);
pts.append(p0p1);
pts.append(generatePoint(snapped,s1,p1,p1Index,s3,p3,p3Index));
pts.append(p2p3);
pts.append(p0p1);
pts.append(p2p3);
pts.append(generatePoint(snapped,s0,p0,p0Index,s2,p2,p2Index));
if (triIndex == 0x09)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-5], pts[sz-4]);
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
case 0x08:
case 0x07:
{
pts.append(generatePoint(snapped,s3,p3,p3Index,s0,p0,p0Index));
pts.append(generatePoint(snapped,s3,p3,p3Index,s2,p2,p2Index));
pts.append(generatePoint(snapped,s3,p3,p3Index,s1,p1,p1Index));
if (triIndex == 0x07)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
}
}
template
void Foam::isoSurfaceCell::generateTriPoints
(
const scalarField& cVals,
const scalarField& pVals,
const Field& cCoords,
const Field& pCoords,
const UList& snappedPoints,
const labelList& snappedCc,
const labelList& snappedPoint,
DynamicList& triPoints,
DynamicList