iso surface correction

This commit is contained in:
mattijs
2009-02-24 13:45:57 +00:00
parent 5a30dd1b01
commit c8944ce200
3 changed files with 391 additions and 218 deletions

View File

@ -85,9 +85,74 @@ bool Foam::isoSurface::isEdgeOfFaceCut
} }
// Get neighbour value and position.
void Foam::isoSurface::getNeighbour
(
const labelList& boundaryRegion,
const volScalarField& cVals,
const label cellI,
const label faceI,
scalar& nbrValue,
point& nbrPoint
) const
{
const labelList& own = mesh_.faceOwner();
const labelList& nei = mesh_.faceNeighbour();
const surfaceScalarField& weights = mesh_.weights();
if (mesh_.isInternalFace(faceI))
{
label nbr = (own[faceI] == cellI ? nei[faceI] : own[faceI]);
nbrValue = cVals[nbr];
nbrPoint = mesh_.cellCentres()[nbr];
}
else
{
label bFaceI = faceI-mesh_.nInternalFaces();
label patchI = boundaryRegion[bFaceI];
const polyPatch& pp = mesh_.boundaryMesh()[patchI];
label patchFaceI = faceI-pp.start();
if (isA<emptyPolyPatch>(pp))
{
// Assume zero gradient
nbrValue = cVals[own[faceI]];
nbrPoint = mesh_.faceCentres()[faceI];
}
else if (pp.coupled())
{
if (!refCast<const coupledPolyPatch>(pp).separated())
{
// Behave as internal face:
// other side value
nbrValue = cVals.boundaryField()[patchI][patchFaceI];
// other side cell centre
nbrPoint = mesh_.C().boundaryField()[patchI][patchFaceI];
}
else
{
// Do some interpolation for now
const scalarField& w = weights.boundaryField()[patchI];
nbrPoint = mesh_.faceCentres()[faceI];
nbrValue =
(1-w[patchFaceI])*cVals[own[faceI]]
+ w[patchFaceI]*cVals.boundaryField()[patchI][patchFaceI];
}
}
else
{
nbrValue = cVals.boundaryField()[patchI][patchFaceI];
nbrPoint = mesh_.faceCentres()[faceI];
}
}
}
// Determine for every face/cell whether it (possibly) generates triangles. // Determine for every face/cell whether it (possibly) generates triangles.
void Foam::isoSurface::calcCutTypes void Foam::isoSurface::calcCutTypes
( (
const labelList& boundaryRegion,
const volScalarField& cVals, const volScalarField& cVals,
const scalarField& pVals const scalarField& pVals
) )
@ -103,7 +168,20 @@ void Foam::isoSurface::calcCutTypes
{ {
// CC edge. // CC edge.
bool ownLower = (cVals[own[faceI]] < iso_); bool ownLower = (cVals[own[faceI]] < iso_);
bool neiLower = (cVals[nei[faceI]] < iso_);
scalar nbrValue;
point nbrPoint;
getNeighbour
(
boundaryRegion,
cVals,
own[faceI],
faceI,
nbrValue,
nbrPoint
);
bool neiLower = (nbrValue < iso_);
if (ownLower != neiLower) if (ownLower != neiLower)
{ {
@ -129,15 +207,29 @@ void Foam::isoSurface::calcCutTypes
if (isA<emptyPolyPatch>(pp)) if (isA<emptyPolyPatch>(pp))
{ {
// Assume zero gradient so owner and neighbour/boundary value equal // Still needs special treatment?
forAll(pp, i) forAll(pp, i)
{ {
bool ownLower = (cVals[own[faceI]] < iso_); bool ownLower = (cVals[own[faceI]] < iso_);
scalar nbrValue;
point nbrPoint;
getNeighbour
(
boundaryRegion,
cVals,
own[faceI],
faceI,
nbrValue,
nbrPoint
);
bool neiLower = (nbrValue < iso_);
const face f = mesh_.faces()[faceI]; const face f = mesh_.faces()[faceI];
if (isEdgeOfFaceCut(pVals, f, ownLower, ownLower)) if (isEdgeOfFaceCut(pVals, f, ownLower, neiLower))
{ {
faceCutType_[faceI] = CUT; faceCutType_[faceI] = CUT;
} }
@ -150,7 +242,20 @@ void Foam::isoSurface::calcCutTypes
forAll(pp, i) forAll(pp, i)
{ {
bool ownLower = (cVals[own[faceI]] < iso_); bool ownLower = (cVals[own[faceI]] < iso_);
bool neiLower = (cVals.boundaryField()[patchI][i] < iso_);
scalar nbrValue;
point nbrPoint;
getNeighbour
(
boundaryRegion,
cVals,
own[faceI],
faceI,
nbrValue,
nbrPoint
);
bool neiLower = (nbrValue < iso_);
if (ownLower != neiLower) if (ownLower != neiLower)
{ {
@ -166,6 +271,7 @@ void Foam::isoSurface::calcCutTypes
faceCutType_[faceI] = CUT; faceCutType_[faceI] = CUT;
} }
} }
faceI++; faceI++;
} }
} }
@ -331,70 +437,6 @@ Foam::pointIndexHit Foam::isoSurface::collapseSurface
} }
// Get neighbour value and position.
void Foam::isoSurface::getNeighbour
(
const labelList& boundaryRegion,
const volScalarField& cVals,
const label cellI,
const label faceI,
scalar& nbrValue,
point& nbrPoint
) const
{
const labelList& own = mesh_.faceOwner();
const labelList& nei = mesh_.faceNeighbour();
const surfaceScalarField& weights = mesh_.weights();
if (mesh_.isInternalFace(faceI))
{
label nbr = (own[faceI] == cellI ? nei[faceI] : own[faceI]);
nbrValue = cVals[nbr];
nbrPoint = mesh_.cellCentres()[nbr];
}
else
{
label bFaceI = faceI-mesh_.nInternalFaces();
label patchI = boundaryRegion[bFaceI];
const polyPatch& pp = mesh_.boundaryMesh()[patchI];
label patchFaceI = faceI-pp.start();
if (isA<emptyPolyPatch>(pp))
{
// Assume zero gradient
nbrValue = cVals[own[faceI]];
nbrPoint = mesh_.faceCentres()[faceI];
}
else if (pp.coupled())
{
if (!refCast<const coupledPolyPatch>(pp).separated())
{
// Behave as internal face:
// other side value
nbrValue = cVals.boundaryField()[patchI][patchFaceI];
// other side cell centre
nbrPoint = mesh_.C().boundaryField()[patchI][patchFaceI];
}
else
{
// Do some interpolation for now
const scalarField& w = weights.boundaryField()[patchI];
nbrPoint = mesh_.faceCentres()[faceI];
nbrValue =
(1-w[patchFaceI])*cVals[own[faceI]]
+ w[patchFaceI]*cVals.boundaryField()[patchI][patchFaceI];
}
}
else
{
nbrValue = cVals.boundaryField()[patchI][patchFaceI];
nbrPoint = mesh_.faceCentres()[faceI];
}
}
}
// Determine per cell centre whether all the intersections get collapsed // Determine per cell centre whether all the intersections get collapsed
// to a single point // to a single point
void Foam::isoSurface::calcSnappedCc void Foam::isoSurface::calcSnappedCc
@ -618,13 +660,14 @@ void Foam::isoSurface::calcSnappedPoint
forAll(pFaces, pFaceI) forAll(pFaces, pFaceI)
{ {
// Create points for all intersections close to point
// (i.e. from pyramid edges)
label faceI = pFaces[pFaceI]; label faceI = pFaces[pFaceI];
const face& f = mesh_.faces()[faceI]; const face& f = mesh_.faces()[faceI];
label own = mesh_.faceOwner()[faceI]; label own = mesh_.faceOwner()[faceI];
// Create points for all intersections close to point // Get neighbour value
// (i.e. from pyramid edges)
scalar nbrValue; scalar nbrValue;
point nbrPoint; point nbrPoint;
getNeighbour getNeighbour
@ -843,6 +886,7 @@ Foam::triSurface Foam::isoSurface::stitchTriPoints
} }
// Determine 'flat hole' situation (see RMT paper). // Determine 'flat hole' situation (see RMT paper).
// Two unconnected triangles get connected because (some of) the edges // Two unconnected triangles get connected because (some of) the edges
// separating them get collapsed. Below only checks for duplicate triangles, // separating them get collapsed. Below only checks for duplicate triangles,
@ -870,6 +914,9 @@ Foam::triSurface Foam::isoSurface::stitchTriPoints
if (nbrTriI > triI && (tris[nbrTriI] == tri)) if (nbrTriI > triI && (tris[nbrTriI] == tri))
{ {
//Pout<< "Duplicate : " << triI << " verts:" << tri
// << " to " << nbrTriI << " verts:" << tris[nbrTriI]
// << endl;
dupTriI = nbrTriI; dupTriI = nbrTriI;
break; break;
} }
@ -1526,8 +1573,14 @@ Foam::isoSurface::isoSurface
{ {
if (debug) if (debug)
{ {
Pout<< "isoSurface :" << nl Pout<< "isoSurface:" << nl
<< " isoField : " << cVals.name() << nl << " isoField : " << cVals.name() << nl
<< " cell min/max : "
<< min(cVals.internalField()) << " / "
<< max(cVals.internalField()) << nl
<< " point min/max : "
<< min(pVals) << " / "
<< max(pVals) << nl
<< " isoValue : " << iso << nl << " isoValue : " << iso << nl
<< " regularise : " << regularise << nl << " regularise : " << regularise << nl
<< " mergeTol : " << mergeTol << nl << " mergeTol : " << mergeTol << nl
@ -1556,15 +1609,7 @@ Foam::isoSurface::isoSurface
} }
} }
// Pre-calculate patch-per-face to avoid whichPatch call.
// Determine if any cut through face/cell
calcCutTypes(cVals, pVals);
// Determine if point is on boundary. Points on boundaries are never
// snapped. Coupled boundaries are handled explicitly so not marked here.
PackedBoolList isBoundaryPoint(mesh_.nPoints());
labelList boundaryRegion(mesh_.nFaces()-mesh_.nInternalFaces()); labelList boundaryRegion(mesh_.nFaces()-mesh_.nInternalFaces());
forAll(patches, patchI) forAll(patches, patchI)
@ -1578,29 +1623,11 @@ Foam::isoSurface::isoSurface
boundaryRegion[faceI-mesh_.nInternalFaces()] = patchI; boundaryRegion[faceI-mesh_.nInternalFaces()] = patchI;
faceI++; faceI++;
} }
}
// Mark all boundary points that are not physically coupled (so anything
// but collocated coupled patches)
if
(
!pp.coupled()
|| refCast<const coupledPolyPatch>(pp).separated()
)
{
label faceI = pp.start();
forAll(pp, i) // Determine if any cut through face/cell
{ calcCutTypes(boundaryRegion, cVals, pVals);
const face& f = mesh_.faces()[faceI];
forAll(f, fp)
{
isBoundaryPoint.set(f[fp], 1);
}
faceI++;
}
}
}
DynamicList<point> snappedPoints(nCutCells_); DynamicList<point> snappedPoints(nCutCells_);
@ -1635,6 +1662,39 @@ Foam::isoSurface::isoSurface
label nCellSnaps = snappedPoints.size(); label nCellSnaps = snappedPoints.size();
// Determine if point is on boundary. Points on boundaries are never
// snapped. Coupled boundaries are handled explicitly so not marked here.
PackedBoolList isBoundaryPoint(mesh_.nPoints());
forAll(patches, patchI)
{
// Mark all boundary points that are not physically coupled (so anything
// but collocated coupled patches)
const polyPatch& pp = patches[patchI];
if
(
!pp.coupled()
|| refCast<const coupledPolyPatch>(pp).separated()
)
{
label faceI = pp.start();
forAll(pp, i)
{
const face& f = mesh_.faces()[faceI];
forAll(f, fp)
{
isBoundaryPoint.set(f[fp], 1);
}
faceI++;
}
}
}
// Per point -1 or a point inside snappedPoints. // Per point -1 or a point inside snappedPoints.
labelList snappedPoint; labelList snappedPoint;
if (regularise) if (regularise)
@ -1727,7 +1787,8 @@ Foam::isoSurface::isoSurface
if (debug) if (debug)
{ {
Pout<< "isoSurface : generated " << triMeshCells.size() Pout<< "isoSurface : generated " << triMeshCells.size()
<< " unmerged triangles." << endl; << " unmerged triangles from " << triPoints.size()
<< " unmerged points." << endl;
} }

View File

@ -135,9 +135,20 @@ class isoSurface
const bool neiLower const bool neiLower
) const; ) const;
void getNeighbour
(
const labelList& boundaryRegion,
const volScalarField& cVals,
const label cellI,
const label faceI,
scalar& nbrValue,
point& nbrPoint
) const;
//- Set faceCutType,cellCutType. //- Set faceCutType,cellCutType.
void calcCutTypes void calcCutTypes
( (
const labelList& boundaryRegion,
const volScalarField& cVals, const volScalarField& cVals,
const scalarField& pVals const scalarField& pVals
); );
@ -156,16 +167,6 @@ class isoSurface
DynamicList<labelledTri, 64>& localTris DynamicList<labelledTri, 64>& localTris
); );
void getNeighbour
(
const labelList& boundaryRegion,
const volScalarField& cVals,
const label cellI,
const label faceI,
scalar& nbrValue,
point& nbrPoint
) const;
//- Determine per cc whether all near cuts can be snapped to single //- Determine per cc whether all near cuts can be snapped to single
// point. // point.
void calcSnappedCc void calcSnappedCc
@ -193,37 +194,39 @@ class isoSurface
template<class Type> template<class Type>
Type generatePoint Type generatePoint
( (
const DynamicList<Type>& snappedPoints,
const scalar s0, const scalar s0,
const Type& p0, const Type& p0,
const label p0Index, const bool hasSnap0,
const Type& snapP0,
const scalar s1, const scalar s1,
const Type& p1, const Type& p1,
const label p1Index const bool hasSnap1,
const Type& snapP1
) const; ) const;
template<class Type> template<class Type>
void generateTriPoints void generateTriPoints
( (
const DynamicList<Type>& snapped,
const scalar s0, const scalar s0,
const Type& p0, const Type& p0,
const label p0Index, const bool hasSnap0,
const Type& snapP0,
const scalar s1, const scalar s1,
const Type& p1, const Type& p1,
const label p1Index, const bool hasSnap1,
const Type& snapP1,
const scalar s2, const scalar s2,
const Type& p2, const Type& p2,
const label p2Index, const bool hasSnap2,
const Type& snapP2,
const scalar s3, const scalar s3,
const Type& p3, const Type& p3,
const label p3Index, const bool hasSnap3,
const Type& snapP3,
DynamicList<Type>& points DynamicList<Type>& points
) const; ) const;
@ -244,7 +247,8 @@ class isoSurface
const scalar neiVal, const scalar neiVal,
const Type& neiPt, const Type& neiPt,
const label neiSnap, const bool hasNeiSnap,
const Type& neiSnapPt,
DynamicList<Type>& triPoints, DynamicList<Type>& triPoints,
DynamicList<label>& triMeshCells DynamicList<label>& triMeshCells

View File

@ -33,15 +33,15 @@ License
template<class Type> template<class Type>
Type Foam::isoSurface::generatePoint Type Foam::isoSurface::generatePoint
( (
const DynamicList<Type>& snappedPoints,
const scalar s0, const scalar s0,
const Type& p0, const Type& p0,
const label p0Index, const bool hasSnap0,
const Type& snapP0,
const scalar s1, const scalar s1,
const Type& p1, const Type& p1,
const label p1Index const bool hasSnap1,
const Type& snapP1
) const ) const
{ {
scalar d = s1-s0; scalar d = s1-s0;
@ -50,13 +50,13 @@ Type Foam::isoSurface::generatePoint
{ {
scalar s = (iso_-s0)/d; scalar s = (iso_-s0)/d;
if (s >= 0.5 && s <= 1 && p1Index != -1) if (hasSnap1 && s >= 0.5 && s <= 1)
{ {
return snappedPoints[p1Index]; return snapP1;
} }
else if (s >= 0.0 && s <= 0.5 && p0Index != -1) else if (hasSnap0 && s >= 0.0 && s <= 0.5)
{ {
return snappedPoints[p0Index]; return snapP0;
} }
else else
{ {
@ -75,23 +75,25 @@ Type Foam::isoSurface::generatePoint
template<class Type> template<class Type>
void Foam::isoSurface::generateTriPoints void Foam::isoSurface::generateTriPoints
( (
const DynamicList<Type>& snapped,
const scalar s0, const scalar s0,
const Type& p0, const Type& p0,
const label p0Index, const bool hasSnap0,
const Type& snapP0,
const scalar s1, const scalar s1,
const Type& p1, const Type& p1,
const label p1Index, const bool hasSnap1,
const Type& snapP1,
const scalar s2, const scalar s2,
const Type& p2, const Type& p2,
const label p2Index, const bool hasSnap2,
const Type& snapP2,
const scalar s3, const scalar s3,
const Type& p3, const Type& p3,
const label p3Index, const bool hasSnap3,
const Type& snapP3,
DynamicList<Type>& points DynamicList<Type>& points
) const ) const
@ -123,29 +125,55 @@ void Foam::isoSurface::generateTriPoints
case 0x0E: case 0x0E:
case 0x01: case 0x01:
points.append(generatePoint(snapped,s0,p0,p0Index,s1,p1,p1Index)); points.append
points.append(generatePoint(snapped,s0,p0,p0Index,s2,p2,p2Index)); (
points.append(generatePoint(snapped,s0,p0,p0Index,s3,p3,p3Index)); generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1)
);
points.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
);
points.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
);
break; break;
case 0x0D: case 0x0D:
case 0x02: case 0x02:
points.append(generatePoint(snapped,s1,p1,p1Index,s0,p0,p0Index)); points.append
points.append(generatePoint(snapped,s1,p1,p1Index,s3,p3,p3Index)); (
points.append(generatePoint(snapped,s1,p1,p1Index,s2,p2,p2Index)); generatePoint(s1,p1,hasSnap1,snapP1,s0,p0,hasSnap0,snapP0)
);
points.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
);
points.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
);
break; break;
case 0x0C: case 0x0C:
case 0x03: case 0x03:
{ {
Type tp1 = generatePoint(snapped,s0,p0,p0Index,s2,p2,p2Index); Type tp1 =
Type tp2 = generatePoint(snapped,s1,p1,p1Index,s3,p3,p3Index); generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2);
Type tp2 =
generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3);
points.append(generatePoint(snapped,s0,p0,p0Index,s3,p3,p3Index)); points.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
);
points.append(tp1); points.append(tp1);
points.append(tp2); points.append(tp2);
points.append(tp2); points.append(tp2);
points.append(generatePoint(snapped,s1,p1,p1Index,s2,p2,p2Index)); points.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
);
points.append(tp1); points.append(tp1);
} }
break; break;
@ -153,23 +181,40 @@ void Foam::isoSurface::generateTriPoints
case 0x0B: case 0x0B:
case 0x04: case 0x04:
{ {
points.append(generatePoint(snapped,s2,p2,p2Index,s0,p0,p0Index)); points.append
points.append(generatePoint(snapped,s2,p2,p2Index,s1,p1,p1Index)); (
points.append(generatePoint(snapped,s2,p2,p2Index,s3,p3,p3Index)); generatePoint(s2,p2,hasSnap2,snapP2,s0,p0,hasSnap0,snapP0)
);
points.append
(
generatePoint(s2,p2,hasSnap2,snapP2,s1,p1,hasSnap1,snapP1)
);
points.append
(
generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3)
);
} }
break; break;
case 0x0A: case 0x0A:
case 0x05: case 0x05:
{ {
Type tp0 = generatePoint(snapped,s0,p0,p0Index,s1,p1,p1Index); Type tp0 =
Type tp1 = generatePoint(snapped,s2,p2,p2Index,s3,p3,p3Index); generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
Type tp1 =
generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
points.append(tp0); points.append(tp0);
points.append(tp1); points.append(tp1);
points.append(generatePoint(snapped,s0,p0,p0Index,s3,p3,p3Index)); points.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
);
points.append(tp0); points.append(tp0);
points.append(generatePoint(snapped,s1,p1,p1Index,s2,p2,p2Index)); points.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
);
points.append(tp1); points.append(tp1);
} }
break; break;
@ -177,23 +222,40 @@ void Foam::isoSurface::generateTriPoints
case 0x09: case 0x09:
case 0x06: case 0x06:
{ {
Type tp0 = generatePoint(snapped,s0,p0,p0Index,s1,p1,p1Index); Type tp0 =
Type tp1 = generatePoint(snapped,s2,p2,p2Index,s3,p3,p3Index); generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
Type tp1 =
generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
points.append(tp0); points.append(tp0);
points.append(generatePoint(snapped,s1,p1,p1Index,s3,p3,p3Index)); points.append
(
generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
);
points.append(tp1); points.append(tp1);
points.append(tp0); points.append(tp0);
points.append(generatePoint(snapped,s0,p0,p0Index,s2,p2,p2Index)); points.append
(
generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
);
points.append(tp1); points.append(tp1);
} }
break; break;
case 0x07: case 0x07:
case 0x08: case 0x08:
points.append(generatePoint(snapped,s3,p3,p3Index,s0,p0,p0Index)); points.append
points.append(generatePoint(snapped,s3,p3,p3Index,s2,p2,p2Index)); (
points.append(generatePoint(snapped,s3,p3,p3Index,s1,p1,p1Index)); generatePoint(s3,p3,hasSnap3,snapP3,s0,p0,hasSnap0,snapP0)
);
points.append
(
generatePoint(s3,p3,hasSnap3,snapP3,s2,p2,hasSnap2,snapP2)
);
points.append
(
generatePoint(s3,p3,hasSnap3,snapP3,s1,p1,hasSnap1,snapP1)
);
break; break;
} }
} }
@ -215,7 +277,8 @@ Foam::label Foam::isoSurface::generateFaceTriPoints
const scalar neiVal, const scalar neiVal,
const Type& neiPt, const Type& neiPt,
const label neiSnap, const bool hasNeiSnap,
const Type& neiSnapPt,
DynamicList<Type>& triPoints, DynamicList<Type>& triPoints,
DynamicList<label>& triMeshCells DynamicList<label>& triMeshCells
@ -234,23 +297,37 @@ Foam::label Foam::isoSurface::generateFaceTriPoints
generateTriPoints generateTriPoints
( (
snappedPoints,
pVals[pointI], pVals[pointI],
pCoords[pointI], pCoords[pointI],
snappedPoint[pointI], snappedPoint[pointI] != -1,
(
snappedPoint[pointI] != -1
? snappedPoints[snappedPoint[pointI]]
: pTraits<Type>::zero
),
pVals[nextPointI], pVals[nextPointI],
pCoords[nextPointI], pCoords[nextPointI],
snappedPoint[nextPointI], snappedPoint[nextPointI] != -1,
(
snappedPoint[nextPointI] != -1
? snappedPoints[snappedPoint[nextPointI]]
: pTraits<Type>::zero
),
cVals[own], cVals[own],
cCoords[own], cCoords[own],
snappedCc[own], snappedCc[own] != -1,
(
snappedCc[own] != -1
? snappedPoints[snappedCc[own]]
: pTraits<Type>::zero
),
neiVal, neiVal,
neiPt, neiPt,
neiSnap, hasNeiSnap,
neiSnapPt,
triPoints triPoints
); );
@ -311,25 +388,6 @@ void Foam::isoSurface::generateTriPoints
<< abort(FatalError); << abort(FatalError);
} }
// Determine neighbouring snap status
labelList neiSnappedCc(mesh_.nFaces()-mesh_.nInternalFaces(), -1);
forAll(patches, patchI)
{
const polyPatch& pp = patches[patchI];
if (pp.coupled())
{
label faceI = pp.start();
forAll(pp, i)
{
neiSnappedCc[faceI-mesh_.nInternalFaces()] =
snappedCc[own[faceI]];
faceI++;
}
}
}
syncTools::swapBoundaryFaceList(mesh_, neiSnappedCc, false);
// Generate triangle points // Generate triangle points
@ -356,7 +414,12 @@ void Foam::isoSurface::generateTriPoints
cVals[nei[faceI]], cVals[nei[faceI]],
cCoords[nei[faceI]], cCoords[nei[faceI]],
snappedCc[nei[faceI]], snappedCc[nei[faceI]] != -1,
(
snappedCc[nei[faceI]] != -1
? snappedPoints[snappedCc[nei[faceI]]]
: pTraits<Type>::zero
),
triPoints, triPoints,
triMeshCells triMeshCells
@ -365,6 +428,34 @@ void Foam::isoSurface::generateTriPoints
} }
// Determine neighbouring snap status
boolList neiSnapped(mesh_.nFaces()-mesh_.nInternalFaces(), false);
List<Type> neiSnappedPoint(neiSnapped.size(), pTraits<Type>::zero);
forAll(patches, patchI)
{
const polyPatch& pp = patches[patchI];
if (pp.coupled())
{
label faceI = pp.start();
forAll(pp, i)
{
label bFaceI = faceI-mesh_.nInternalFaces();
label snappedIndex = snappedCc[own[faceI]];
if (snappedIndex != -1)
{
neiSnapped[bFaceI] = true;
neiSnappedPoint[bFaceI] = snappedPoints[snappedIndex];
}
faceI++;
}
}
}
syncTools::swapBoundaryFaceList(mesh_, neiSnapped, false);
syncTools::swapBoundaryFaceList(mesh_, neiSnappedPoint, false);
forAll(patches, patchI) forAll(patches, patchI)
{ {
const polyPatch& pp = patches[patchI]; const polyPatch& pp = patches[patchI];
@ -372,9 +463,10 @@ void Foam::isoSurface::generateTriPoints
if if
( (
isA<processorPolyPatch>(pp) isA<processorPolyPatch>(pp)
&& refCast<const processorPolyPatch>(pp).owner()
&& !refCast<const processorPolyPatch>(pp).separated() && !refCast<const processorPolyPatch>(pp).separated()
) )
{
//if (refCast<const processorPolyPatch>(pp).owner())
{ {
label faceI = pp.start(); label faceI = pp.start();
@ -382,6 +474,18 @@ void Foam::isoSurface::generateTriPoints
{ {
if (faceCutType_[faceI] != NOTCUT) if (faceCutType_[faceI] != NOTCUT)
{ {
label bFaceI = faceI-mesh_.nInternalFaces();
if
(
neiSnapped[bFaceI]
&& (neiSnappedPoint[bFaceI]==pTraits<Type>::zero)
)
{
FatalErrorIn("isoSurface::generateTriPoints(..)")
<< "problem:" << abort(FatalError);
}
generateFaceTriPoints generateFaceTriPoints
( (
cVals, cVals,
@ -397,7 +501,8 @@ void Foam::isoSurface::generateTriPoints
cVals.boundaryField()[patchI][i], cVals.boundaryField()[patchI][i],
cCoords.boundaryField()[patchI][i], cCoords.boundaryField()[patchI][i],
neiSnappedCc[faceI-mesh_.nInternalFaces()], neiSnapped[faceI-mesh_.nInternalFaces()],
neiSnappedPoint[faceI-mesh_.nInternalFaces()],
triPoints, triPoints,
triMeshCells triMeshCells
@ -406,6 +511,7 @@ void Foam::isoSurface::generateTriPoints
faceI++; faceI++;
} }
} }
}
else if (isA<emptyPolyPatch>(pp)) else if (isA<emptyPolyPatch>(pp))
{ {
// Assume zero-gradient. But what about coordinates? // Assume zero-gradient. But what about coordinates?
@ -430,7 +536,8 @@ void Foam::isoSurface::generateTriPoints
cVals[own[faceI]], cVals[own[faceI]],
cCoords.boundaryField()[patchI][i], cCoords.boundaryField()[patchI][i],
-1, // fc not snapped false, // fc not snapped
pTraits<Type>::zero,
triPoints, triPoints,
triMeshCells triMeshCells
@ -462,7 +569,8 @@ void Foam::isoSurface::generateTriPoints
cVals.boundaryField()[patchI][i], cVals.boundaryField()[patchI][i],
cCoords.boundaryField()[patchI][i], cCoords.boundaryField()[patchI][i],
-1, // fc not snapped false, // fc not snapped
pTraits<Type>::zero,
triPoints, triPoints,
triMeshCells triMeshCells