/*---------------------------------------------------------------------------*\
========= |
\\ / 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) 2018-2022 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 "isoSurfacePoint.H"
#include "polyMesh.H"
#include "syncTools.H"
#include "surfaceFields.H"
#include "OFstream.H"
#include "meshTools.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template
Foam::tmp>
Foam::isoSurfacePoint::adaptPatchFields
(
const VolumeField& fld
) const
{
auto tslice = tmp>::New
(
IOobject
(
fld.name(),
fld.instance(),
fld.db(),
IOobject::NO_READ,
IOobject::NO_WRITE,
IOobject::NO_REGISTER
),
fld, // internal field
true // preserveCouples
);
auto& sliceFld = tslice.ref();
const fvMesh& mesh = fld.mesh();
const polyBoundaryMesh& patches = mesh.boundaryMesh();
auto& sliceFldBf = sliceFld.boundaryFieldRef();
forAll(patches, patchi)
{
const polyPatch& pp = patches[patchi];
if
(
isA(pp)
&& pp.size() != sliceFldBf[patchi].size()
)
{
// Clear old value. Cannot resize it since is a slice.
// Set new value we can change
sliceFldBf.release(patchi);
sliceFldBf.set
(
patchi,
new calculatedFvPatchField
(
mesh.boundary()[patchi],
sliceFld
)
);
// Note: cannot use patchInternalField since uses emptyFvPatch::size
// Do our own internalField instead.
const labelUList& faceCells =
mesh.boundary()[patchi].patch().faceCells();
Field& pfld = sliceFldBf[patchi];
pfld.setSize(faceCells.size());
forAll(faceCells, i)
{
pfld[i] = sliceFld[faceCells[i]];
}
}
else if (isA(pp))
{
// Already has interpolate as value
}
else if (isA(pp))
{
fvPatchField& pfld = const_cast&>
(
sliceFldBf[patchi]
);
tmp> f = lerp
(
pfld.patchNeighbourField(),
pfld.patchInternalField(),
mesh.weights().boundaryField()[patchi]
);
bitSet isCollocated
(
collocatedFaces(refCast(pp))
);
forAll(isCollocated, i)
{
if (!isCollocated[i])
{
pfld[i] = f()[i];
}
}
}
}
return tslice;
}
template
Type Foam::isoSurfacePoint::generatePoint
(
const scalar s0,
const Type& p0,
const bool hasSnap0,
const Type& snapP0,
const scalar s1,
const Type& p1,
const bool hasSnap1,
const Type& snapP1
) const
{
const scalar d = s1-s0;
if (mag(d) > VSMALL)
{
const scalar s = (iso_-s0)/d;
if (hasSnap1 && s >= 0.5 && s <= 1)
{
return snapP1;
}
else if (hasSnap0 && s >= 0.0 && s <= 0.5)
{
return snapP0;
}
else
{
return s*p1 + (1.0-s)*p0;
}
}
else
{
constexpr scalar s = 0.4999;
return s*p1 + (1.0-s)*p0;
}
}
template
void Foam::isoSurfacePoint::generateTriPoints
(
const scalar s0,
const Type& p0,
const bool hasSnap0,
const Type& snapP0,
const scalar s1,
const Type& p1,
const bool hasSnap1,
const Type& snapP1,
const scalar s2,
const Type& p2,
const bool hasSnap2,
const Type& snapP2,
const scalar s3,
const Type& p3,
const bool hasSnap3,
const Type& snapP3,
DynamicList& pts
) const
{
// Note: cannot use simpler isoSurfaceCell::generateTriPoints since
// the need here to sometimes pass in remote 'snappedPoints'
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(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1)
);
pts.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
);
pts.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
);
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(s1,p1,hasSnap1,snapP1,s0,p0,hasSnap0,snapP0)
);
pts.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
);
pts.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
);
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(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2);
Type p1p3 =
generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3);
pts.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
);
pts.append(p1p3);
pts.append(p0p2);
pts.append(p1p3);
pts.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
);
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(s2,p2,hasSnap2,snapP2,s0,p0,hasSnap0,snapP0)
);
pts.append
(
generatePoint(s2,p2,hasSnap2,snapP2,s1,p1,hasSnap1,snapP1)
);
pts.append
(
generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3)
);
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(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
Type p2p3 =
generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
pts.append(p0p1);
pts.append(p2p3);
pts.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
);
pts.append(p0p1);
pts.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
);
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(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
Type p2p3 =
generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
pts.append(p0p1);
pts.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
);
pts.append(p2p3);
pts.append(p0p1);
pts.append(p2p3);
pts.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
);
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(s3,p3,hasSnap3,snapP3,s0,p0,hasSnap0,snapP0)
);
pts.append
(
generatePoint(s3,p3,hasSnap3,snapP3,s2,p2,hasSnap2,snapP2)
);
pts.append
(
generatePoint(s3,p3,hasSnap3,snapP3,s1,p1,hasSnap1,snapP1)
);
if (triIndex == 0x07)
{
// Flip normals
const label sz = pts.size();
std::swap(pts[sz-2], pts[sz-1]);
}
}
break;
}
}
template
Foam::label Foam::isoSurfacePoint::generateFaceTriPoints
(
const volScalarField& cVals,
const scalarField& pVals,
const VolumeField& cCoords,
const Field& pCoords,
const UList& snappedPoints,
const labelList& snappedCc,
const labelList& snappedPoint,
const label facei,
const scalar neiVal,
const Type& neiPt,
const bool hasNeiSnap,
const Type& neiSnapPt,
DynamicList& triPoints,
DynamicList