First version of cuttingPatchFringe algorithm

Similar to donorBasedOverlapFringe but better in a sense that it relies on
determining inside/outside region by creating an STL of the patch from donor
region. Note: still has bugs
This commit is contained in:
Vuko Vukcevic 2019-05-20 16:16:13 +02:00
parent db6dd0cf97
commit db859b968c
3 changed files with 936 additions and 0 deletions

View file

@ -5,6 +5,7 @@ oversetFringe/oversetFringe/newOversetFringe.C
oversetFringe/manualFringe/manualFringe.C oversetFringe/manualFringe/manualFringe.C
oversetFringe/faceCellsFringe/faceCellsFringe.C oversetFringe/faceCellsFringe/faceCellsFringe.C
oversetFringe/donorBasedLayeredOverlapFringe/donorBasedLayeredOverlapFringe.C oversetFringe/donorBasedLayeredOverlapFringe/donorBasedLayeredOverlapFringe.C
oversetFringe/cuttingPatchFringe/cuttingPatchFringe.C
oversetFringe/overlapFringe/overlapFringe/overlapFringe.C oversetFringe/overlapFringe/overlapFringe/overlapFringe.C
oversetFringe/overlapFringe/layeredOverlapFringe/layeredOverlapFringe.C oversetFringe/overlapFringe/layeredOverlapFringe/layeredOverlapFringe.C
oversetFringe/overlapFringe/adaptiveOverlapFringe/adaptiveOverlapFringe.C oversetFringe/overlapFringe/adaptiveOverlapFringe/adaptiveOverlapFringe.C

View file

@ -0,0 +1,772 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.1
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend 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.
foam-extend 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 foam-extend. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "cuttingPatchFringe.H"
#include "oversetMesh.H"
#include "oversetRegion.H"
#include "faceCellsFringe.H"
#include "oversetRegion.H"
#include "polyPatchID.H"
#include "addToRunTimeSelectionTable.H"
#include "syncTools.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(cuttingPatchFringe, 0);
addToRunTimeSelectionTable
(
oversetFringe,
cuttingPatchFringe,
dictionary
);
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::cuttingPatchFringe::init() const
{
// Set size of the list containing IDs
connectedRegionIDs_.setSize(connectedRegionNames_.size());
// Get list of all overset regions
const PtrList<oversetRegion>& allRegions =
this->region().overset().regions();
// Create list of all region names for easy lookup
wordList allRegionNames(allRegions.size());
forAll (allRegionNames, arI)
{
allRegionNames[arI] = allRegions[arI].name();
}
// Loop through all regions and check whether the overlap has been found
forAll (connectedRegionNames_, crI)
{
// Get name of this connected region
const word& crName = connectedRegionNames_[crI];
// Find this region in the list of all regions
const label regionID = findIndex(allRegionNames, crName);
if (regionID == -1)
{
FatalErrorIn("void cuttingPatchFringe::init() const")
<< "Region " << crName << " not found in list of regions."
<< "List of overset regions: " << allRegionNames
<< abort(FatalError);
}
// Collect the region index in the list
connectedRegionIDs_[crI] = regionID;
// Sanity check: if the specified connected donor region has more than 1
// donor regions, this fringe algorithm is attempted to be used for
// something that's not intended for. Issue an error
if (allRegions[regionID].donorRegions().size() != 1)
{
FatalErrorIn("void cuttingPatchFringe::init() const")
<< "Region " << crName << " specified as connected region, but"
<< " that region has "
<< allRegions[regionID].donorRegions().size()
<< " donor regions."
<< abort(FatalError);
}
// Sanity check whether the donor region of connected region is actually
// this region
if (allRegions[regionID].donorRegions()[0] != this->region().index())
{
FatalErrorIn("void cuttingPatchFringe::init() const")
<< "The donor region of region " << crName
<< " should be only region " << this->region().name()
<< abort(FatalError);
}
}
}
void Foam::cuttingPatchFringe::calcAddressing() const
{
// Make sure that either acceptorsPtr is unnalocated or if it is allocated,
// that it is empty
if (acceptorsPtr_ && !acceptorsPtr_->empty())
{
FatalErrorIn
(
"void Foam::cuttingPatchFringe::calcAddressing() const"
) << "Addressing already calculated"
<< abort(FatalError);
}
if (!isInitialized_)
{
// This is the first call, initialize the data and set flag to true
init();
isInitialized_ = true;
}
// Get list of all overset regions
const PtrList<oversetRegion>& allRegions =
this->region().overset().regions();
// Sets containing all acceptors and all holes for all connected regions
const polyMesh& mesh = this->mesh();
labelHashSet allAcceptors(0.02*mesh.nCells());
labelHashSet allFringeHoles(0.02*mesh.nCells());
if (debug)
{
Info<< "All dependent fringes are ready."
<< " Starting face cells cut patch fringe assembly..." << endl;
}
// Loop through connected regions
forAll (connectedRegionIDs_, crI)
{
// Get ID of this region
const label& regionID = connectedRegionIDs_[crI];
// Get fringe of the connected region
const oversetFringe& fringe = allRegions[regionID].fringe();
// If this is not faceCells fringe, issue an Error. This fringe
// selection algorithm is intended to work only with faceCells fringe on
// the other side. VV, 9/Apr/2019
if (!isA<faceCellsFringe>(fringe))
{
FatalErrorIn
(
"void Foam::cuttingPatchFringe::"
"updateIteration(donorAcceptorList&) const"
) << "cuttingPatch fringe is designed to work"
<< " with faceCells fringe as a connected region fringe."
<< nl
<< "Connected overset region " << allRegions[regionID].name()
<< " has " << fringe.type() << " fringe type. "
<< nl
<< "Proceed with care!"
<< abort(FatalError);
}
const faceCellsFringe& fcFringe =
refCast<const faceCellsFringe>(fringe);
// Get patch names from faceCells fringe
const wordList& fcPatchNames = fcFringe.patchNames();
// Find the patches
labelHashSet patchIDs;
forAll (fcPatchNames, nameI)
{
// Get polyPatchID and check if valid
const polyPatchID fringePatch
(
fcPatchNames[nameI],
mesh.boundaryMesh()
);
if (!fringePatch.active())
{
FatalErrorIn
(
"void cuttingPatchFringe::calcAddressing const"
) << "Fringe patch " << fcPatchNames[nameI]
<< " for region " << allRegions[regionID].name()
<< " cannot be found."
<< abort(FatalError);
}
// Store patch ID in the set
patchIDs.insert(fringePatch.index());
}
// Note: same code as in oversetRegion::calcHoleTriMesh. Consider
// refactoring. VV, 20/May/2019
// Make and invert local triSurface
triFaceList triFaces;
pointField triPoints;
// Memory management
{
triSurface ts = triSurfaceTools::triangulate
(
mesh.boundaryMesh(),
patchIDs
);
// Clean mutiple points and zero-sized triangles
ts.cleanup(false);
triFaces.setSize(ts.size());
triPoints = ts.points();
forAll (ts, tsI)
{
triFaces[tsI] = ts[tsI].reverseFace();
}
}
if (Pstream::parRun())
{
// Combine all faces and points into a single list
List<triFaceList> allTriFaces(Pstream::nProcs());
List<pointField> allTriPoints(Pstream::nProcs());
allTriFaces[Pstream::myProcNo()] = triFaces;
allTriPoints[Pstream::myProcNo()] = triPoints;
Pstream::gatherList(allTriFaces);
Pstream::scatterList(allTriFaces);
Pstream::gatherList(allTriPoints);
Pstream::scatterList(allTriPoints);
// Re-pack points and faces
label nTris = 0;
label nPoints = 0;
forAll (allTriFaces, procI)
{
nTris += allTriFaces[procI].size();
nPoints += allTriPoints[procI].size();
}
// Pack points
triPoints.setSize(nPoints);
// Prepare point renumbering array
labelListList renumberPoints(Pstream::nProcs());
nPoints = 0;
forAll (allTriPoints, procI)
{
const pointField& ptp = allTriPoints[procI];
renumberPoints[procI].setSize(ptp.size());
labelList& procRenumberPoints = renumberPoints[procI];
forAll (ptp, ptpI)
{
triPoints[nPoints] = ptp[ptpI];
procRenumberPoints[ptpI] = nPoints;
nPoints++;
}
}
// Pack triangles and renumber into complete points on the fly
triFaces.setSize(nTris);
nTris = 0;
forAll (allTriFaces, procI)
{
const triFaceList& ptf = allTriFaces[procI];
const labelList& procRenumberPoints = renumberPoints[procI];
forAll (ptf, ptfI)
{
const triFace& procFace = ptf[ptfI];
triFace& renumberFace = triFaces[nTris];
forAll (renumberFace, rfI)
{
renumberFace[rfI] =
procRenumberPoints[procFace[rfI]];
}
nTris++;
}
}
}
// Make a complete triSurface from local data
triSurface patchTriMesh(triFaces, triPoints);
// Clean up duplicate points and zero sized triangles
patchTriMesh.cleanup(false);
// Get this region
const oversetRegion& myRegion = this->region();
// Debug: write down the tri mesh
if (Pstream::master())
{
patchTriMesh.write
(
word
(
"patchTriSurface_region" + myRegion.name() +
"_connectedRegion" + allRegions[regionID].name() + ".vtk"
)
);
}
// Create the tri surface search object
triSurfaceSearch patchTriSearch(patchTriMesh);
// Get cells in this region
const labelList& myRegionCells = myRegion.regionCells();
// Get cell centres for inside-outside search using search object
vectorField myCC(myRegionCells.size());
// Cell centres from polyMesh
const vectorField& cc = mesh.cellCentres();
forAll (myCC, i)
{
myCC[i] = cc[myRegionCells[i]];
}
// Inside mask: all cells within search object will be marked
boolList insideMask(mesh.nCells(), false);
// Get inside cells for cells in my region only
boolList myRegionInsideMask = patchTriSearch.calcInside(myCC);
// Note: insideMask has the size of all mesh cells and
// myRegionInsideMask has the size of cells in this region
forAll (myRegionInsideMask, i)
{
insideMask[myRegionCells[i]] = myRegionInsideMask[i];
}
// Get necessary mesh data (from polyMesh/primitiveMesh)
const cellList& meshCells = mesh.cells();
const unallocLabelList& owner = mesh.faceOwner();
const unallocLabelList& neighbour = mesh.faceNeighbour();
// Collect all inside cells either as hole cells or acceptor
// cells. For the first iteration, acceptor cells are the cells that
// have at least one neighbour cell that is not marked
labelHashSet acceptors(myRegionCells.size()/10);
// Bool list for collecting faces with at least one unmarked
// cell (to determine the acceptors for the first iteration)
boolList hasUnmarkedCell(mesh.nFaces(), false);
// Loop through all cells
forAll (insideMask, cellI)
{
if (!insideMask[cellI])
{
// This cell is not inside (it is unmarked). Loop through
// its faces and set the flag
const cell& cFaces = meshCells[cellI];
forAll (cFaces, i)
{
// Set the mark for this global face and break out
hasUnmarkedCell[cFaces[i]] = true;
break;
}
}
}
// Sync the face list across processor boundaries
syncTools::syncFaceList
(
mesh,
hasUnmarkedCell,
orOp<bool>(),
false
);
// Mark-up for all hole faces
boolList holeFaceMask(mesh.nFaces(), false);
// Loop again through all cells and collect marked ones into
// acceptors or holes, depending on whether they have unmarked cell
// as a neighbour (indicating an acceptor)
forAll (insideMask, cellI)
{
if (insideMask[cellI])
{
// This cell is inside the covered region
const cell& cFaces = meshCells[cellI];
// Helper variable to determine whether this marked cell
// should be collected as a hole
bool isHole = true;
forAll (cFaces, i)
{
// Get global face index
const label& faceI = cFaces[i];
// Check whether this neighbour is unmarked
if (hasUnmarkedCell[faceI])
{
// This cell has unmarked neighbour, collect it into
// the acceptor list, set isHole to false and break
acceptors.insert(cellI);
isHole = false;
break;
}
}
// If this is a hole, collect it and mark its faces
if (isHole)
{
holeCellMask[cellI] = true;
// Loop through cell faces and mark them
const cell& cFaces = meshCells[cellI];
forAll (cFaces, i)
{
holeFaceMask[cFaces[i]] = true;
}
}
} // End if cell is inside
} // End for all cells
// Hash set containing new acceptors (for successive iterations)
labelHashSet newAcceptors(acceptors.size());
// Now that we have the initial set of acceptors (and holes), loop
// nLayers away from initial donors
for (label i = 0; i < nLayers_; ++i)
{
// Face markup for propagation
boolList propagateFace(mesh.nFaces(), false);
// Loop through all donors and mark faces that are around hole
// cells. This way, we make sure that we go towards the correct,
// inside direction
forAllConstIter (labelHashSet, acceptors, iter)
{
// Get the cell index and the cell
const label& cellI = iter.key();
const cell& cFaces = meshCells[cellI];
// Loop through all faces of the cell
forAll (cFaces, i)
{
// Get face index (global)
const label& faceI = cFaces[i];
if (holeFaceMask[faceI])
{
// This is a hole face, we are moving in the right
// direction. Mark the face for propagation
propagateFace[faceI] = true;
}
} // End for all faces of the cell
} // End for all donor cells
// Sync the face list across processor boundaries
syncTools::syncFaceList
(
mesh,
propagateFace,
orOp<bool>(),
false
);
// Loop through all faces and append acceptors
for (label faceI = 0; faceI < mesh.nInternalFaces(); ++faceI)
{
if (propagateFace[faceI])
{
// Face is marked, select owner or neighbour
const label& own = owner[faceI];
const label& nei = neighbour[faceI];
// Either owner or neighbour may be hole, not both
if (holeCellMask[own])
{
// Owner cell is a hole, insert it
newAcceptors.insert(own);
// Update hole mask
holeCellMask[own] = false;
}
else if (holeCellMask[nei])
{
// Neighbour cell is a hole, insert it
newAcceptors.insert(nei);
// Update hole mask
holeCellMask[nei] = false;
}
// Also update hole face mask for next iteration
holeFaceMask[faceI] = false;
}
}
// Loop through boundary faces
for
(
label faceI = mesh.nInternalFaces();
faceI < mesh.nFaces();
++faceI
)
{
if (propagateFace[faceI])
{
// Face is marked, select owner if this is the right
// side. Neighbour handled on the other side
const label& own = owner[faceI];
if (holeCellMask[own])
{
// Face cell is a hole, insert it
newAcceptors.insert(own);
// Update hole mask
holeCellMask[own] = false;
}
// Also update hole face mask for next iteration
holeFaceMask[faceI] = false;
}
}
// Transfer newAcceptors into acceptors for next iteration or
// for final assembly. Resize newAcceptors accordingly
acceptors.transfer(newAcceptors);
newAcceptors.resize(acceptors.size());
} // End for specified number of layers
// At this point, we have the final set of acceptors and we marked
// all cells that should be holes. Collect holes into hash set (could be
// optimized by using dynamic lists)
labelHashSet fringeHoles(myRegionCells.size()/10);
forAll (holeCellMask, cellI)
{
if (holeCellMask[cellI])
{
fringeHoles.insert(cellI);
}
}
// Finally, we have fringe holes and acceptors and we need to add them
// to the list containing all acceptors and holes (for all connected
// regions)
allAcceptors += acceptors;
allFringeHoles += fringeHoles;
}
// Set acceptors and holes from the data for all regions
acceptorsPtr_ = new labelList(allAcceptors.sortedToc());
fringeHolesPtr_ = new labelList(allFringeHoles.sortedToc());
if (debug)
{
Info<< "In cuttingPatchFringe::calcAddressing() const" << nl
<< "Found " << acceptorsPtr_->size() << " acceptors." << nl
<< "Found " << fringeHolesPtr_->size() << " fringe holes." << endl;
}
}
void Foam::cuttingPatchFringe::clearAddressing() const
{
deleteDemandDrivenData(fringeHolesPtr_);
deleteDemandDrivenData(acceptorsPtr_);
deleteDemandDrivenData(finalDonorAcceptorsPtr_);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
// Construct from dictionary
Foam::cuttingPatchFringe::cuttingPatchFringe
(
const fvMesh& mesh,
const oversetRegion& region,
const dictionary& dict
)
:
oversetFringe(mesh, region, dict),
connectedRegionNames_(dict.lookup("connectedRegions")),
connectedRegionIDs_(),
nLayers_(readLabel(dict.lookup("nLayers"))),
fringeHolesPtr_(nullptr),
acceptorsPtr_(nullptr),
finalDonorAcceptorsPtr_(nullptr),
isInitialized_(false)
{
// Sanity check number of layers: must be greater than 0
if (nLayers_ < 1)
{
FatalIOErrorIn
(
"cuttingPatchFringe::"
"cuttingPatchFringe\n"
"(\n"
" const fvMesh& mesh,\n"
" const oversetRegion& region,\n"
" const dictionary& dict\n"
")",
dict
) << "Invalid number of layers specified, nLayers = " << nLayers_
<< nl
<< "The number should be greater than 0."
<< abort(FatalError);
}
// Preferably, the number of layers should be at least 2
if (nLayers_ == 1)
{
WarningIn
(
"cuttingPatchFringe::"
"cuttingPatchFringe\n"
"(\n"
" const fvMesh& mesh,\n"
" const oversetRegion& region,\n"
" const dictionary& dict\n"
")"
) << "You have specified nLayers = " << nLayers_
<< nl
<< "We strongly advise to use at least 2 layers to avoid"
<< " possibility of having acceptors that cannot find decent"
<< " donors on the other side."
<< abort(FatalError);
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::cuttingPatchFringe::~cuttingPatchFringe()
{
clearAddressing();
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::cuttingPatchFringe::updateIteration
(
donorAcceptorList& donorAcceptorRegionData
) const
{
// If the donorAcceptor list has been allocated, something went wrong with
// the iteration procedure (not-updated flag): this function has been called
// more than once, which should not happen for cuttingPatchFringe
if (finalDonorAcceptorsPtr_)
{
FatalErrorIn
(
"cuttingPatchFringe::"
"updateIteration(donorAcceptorList&) const"
) << "finalDonorAcceptorPtr_ already allocated. Something went "
<< "wrong with the iteration procedure (flag was not updated)."
<< nl
<< "This should not happen for cuttingPatchFringe."
<< abort(FatalError);
}
// Allocate the list by reusing the argument list
finalDonorAcceptorsPtr_ = new donorAcceptorList
(
donorAcceptorRegionData,
true
);
// Set the flag to true and return
updateSuitableOverlapFlag(true);
return foundSuitableOverlap();
}
const Foam::labelList& Foam::cuttingPatchFringe::fringeHoles() const
{
if (!fringeHolesPtr_)
{
calcAddressing();
}
return *fringeHolesPtr_;
}
const Foam::labelList& Foam::cuttingPatchFringe::candidateAcceptors() const
{
if (!acceptorsPtr_)
{
calcAddressing();
}
return *acceptorsPtr_;
}
Foam::donorAcceptorList&
Foam::cuttingPatchFringe::finalDonorAcceptors() const
{
if (!finalDonorAcceptorsPtr_)
{
FatalErrorIn("cuttingPatchFringe::finalDonorAcceptors()")
<< "finalDonorAcceptorPtr_ not allocated. Make sure you have"
<< " called cuttingPatchFringe::updateIteration() before"
<< " asking for final set of donor/acceptor pairs."
<< abort(FatalError);
}
if (!foundSuitableOverlap())
{
FatalErrorIn("cuttingPatchFringe::finalDonorAcceptors()")
<< "Attemted to access finalDonorAcceptors but suitable overlap "
<< "has not been found. This is not allowed. "
<< abort(FatalError);
}
return *finalDonorAcceptorsPtr_;
}
void Foam::cuttingPatchFringe::update() const
{
Info<< "cuttingPatchFringe::update() const" << endl;
// Clear out
clearAddressing();
// Set flag to false and clear final donor/acceptors only
deleteDemandDrivenData(finalDonorAcceptorsPtr_);
updateSuitableOverlapFlag(false);
}
// ************************************************************************* //

View file

@ -0,0 +1,163 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.1
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend 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.
foam-extend 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 foam-extend. If not, see <http://www.gnu.org/licenses/>.
Class
cuttingPatchFringe
Description
Fringe algorithm based on collecting all the cells cut by its donor region
patch (assuming faceCells algorithm is used). Inside cells are marked and
moved nLayers towards the interior to define the acceptors as robustly as
possible.
Author
Vuko Vukcevic, Wikki Ltd. All rights reserved.
SourceFiles
cuttingPatchFringe.C
\*---------------------------------------------------------------------------*/
#ifndef cuttingPatchFringe_H
#define cuttingPatchFringe_H
#include "oversetFringe.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class cuttingPatchFringe Declaration
\*---------------------------------------------------------------------------*/
class cuttingPatchFringe
:
public oversetFringe
{
// Private data
//- Names of connected regions. Looked up on construction
wordList connectedRegionNames_;
//- Regions IDs from which the donors will be collected as a starting
// point. Note: initialized in init private member function because we
// cannot initialize it in constructor. This is because certain overset
// regions (and their fringes) may not be initialized at this point.
mutable labelList connectedRegionIDs_;
//- How many layers to move away from connected region donors to define
// acceptor (and holes)
label nLayers_;
//- Fringe hole cells
mutable labelList* fringeHolesPtr_;
//- Acceptor cells
mutable labelList* acceptorsPtr_;
//- Final donor/acceptor pairs for this region (fringe)
mutable donorAcceptorList* finalDonorAcceptorsPtr_;
//- Initialization helper
mutable bool isInitialized_;
// Private Member Functions
//- Initialization
void init() const;
//- Calculate hole and acceptor addressing
void calcAddressing() const;
//- Clear addressing
void clearAddressing() const;
public:
//- Runtime type information
TypeName("cuttingPatch");
// Constructors
//- Construct from dictionary
cuttingPatchFringe
(
const fvMesh& mesh,
const oversetRegion& region,
const dictionary& dict
);
//- Disallow default bitwise copy construct
cuttingPatchFringe
(
const cuttingPatchFringe&
) = delete;
//- Disallow default bitwise assignment
void operator=(const cuttingPatchFringe&) = delete;
// Destructor
virtual ~cuttingPatchFringe();
// Member Functions
//- Update iteration. Note: invalidates parameter
virtual bool updateIteration
(
donorAcceptorList& donorAcceptorRegionData
) const;
//- Return list of deactivated (hole) cells
// Fringe hole cells are collected in addition to geometric hole
// cells, which fall outside of all donor regions
virtual const labelList& fringeHoles() const;
//- Return list of acceptor cells
virtual const labelList& candidateAcceptors() const;
//- Return list of final donor acceptor pairs. Note: caller may
// invalidate finalDonorAcceptorsPtr_ for optimisation purposes
virtual donorAcceptorList& finalDonorAcceptors() const;
//- Update the fringe
virtual void update() const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //