74 "\n Usage: holeToFace <faceSet> ((x0 y0 z0) (x1 y1 z1))\n\n"
75 " Select faces disconnecting the individual regions"
76 " (each indicated by a locations).\n"
302void Foam::holeToFace::writeFaces
305 const bitSet& faceFld
310 Pout<<
"Writing " << faceFld.count() <<
" faces to " << str.
name() <<
endl;
312 for (
const label facei : faceFld)
319void Foam::holeToFace::calculateDistance
322 const bitSet& isBlockedCell,
328 if (isBlockedCell.size() != mesh_.nCells())
357 for (
const label celli : isBlockedCell)
376 mesh_.globalData().nTotalCells()+1
383 if (!isBlockedCell[celli])
385 if (!cellData[celli].valid(deltaCalc.data()))
397 cellDist[celli] = cellData[celli].distance();
406 if (!faceData[facei].valid(deltaCalc.data()))
418 faceDist[facei] = faceData[facei].distance();
425Foam::bitSet Foam::holeToFace::frontFaces
427 const bitSet& isSurfaceFace,
432 const labelList& faceOwner = mesh_.faceOwner();
433 const labelList& faceNeighbour = mesh_.faceNeighbour();
436 bitSet isFrontFace(mesh_.nFaces());
437 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
439 if (!isSurfaceFace[facei])
441 const label ownHole = isHoleCell[faceOwner[facei]];
442 const label neiHole = isHoleCell[faceNeighbour[facei]];
444 if (ownHole != neiHole)
446 unsigned int masks = locationFaces[facei];
450 <<
" at:" << mesh_.faceCentres()[facei]
459 isFrontFace.set(facei);
466 bitSet isHoleNeiCell(mesh_.nBoundaryFaces());
470 label bFacei =
pp.start()-mesh_.nInternalFaces();
475 isHoleNeiCell[bFacei] = isHoleCell[celli];
483 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
485 if (!isSurfaceFace[facei])
487 const label ownHole = isHoleCell[faceOwner[facei]];
488 const label neiHole = isHoleNeiCell[facei-mesh_.nInternalFaces()];
490 if (ownHole != neiHole)
492 unsigned int masks = locationFaces[facei];
496 <<
" at:" << mesh_.faceCentres()[facei]
505 isFrontFace.set(facei);
514Foam::bitSet Foam::holeToFace::findClosure
516 const bitSet& isSurfaceFace,
517 const bitSet& isCandidateHoleCell,
521 const labelList& faceOwner = mesh_.faceOwner();
522 const labelList& faceNeighbour = mesh_.faceNeighbour();
524 if (zonePoints_.size() < 2)
530 if (zonePoints_.size() > 31)
539 bitSet isHoleCell(isCandidateHoleCell);
540 for (
const labelList& zoneCells : locationCells)
542 for (
const label celli : zoneCells)
546 isHoleCell.unset(celli);
551 bitSet notHoleCell(isHoleCell);
556 Pout<<
"holeToFace::findClosure :"
558 <<
"holeToFace::findClosure :"
559 <<
" initial blocked faces:" << isSurfaceFace.count()
560 <<
" candidate closure cells:" << isHoleCell.count()
565 labelList surfaceCellDist(mesh_.nCells(), -1);
566 labelList surfaceNeiCellDist(mesh_.nBoundaryFaces(), -1);
567 labelList surfaceFaceDist(mesh_.nFaces(), -1);
590 Pout<<
"holeToFace::findClosure :"
591 <<
" calculated topological distance to initial blocked faces."
592 <<
" max distance:" <<
gMax(surfaceCellDist)
604 forAll(locationCells, zonei)
611 const labelList& zoneCells = locationCells[zonei];
612 for (
const label celli : zoneCells)
617 const cell& cFaces = mesh_.cells()[celli];
625 bitSet isSeedFace(mesh_.nFaces(), seedFaces);
632 seedFaces = isSeedFace.toc();
649 const unsigned int mask = (1u << zonei);
652 if (faceDist[facei] >= 0)
654 locationFaces[facei] |= mask;
661 writeFaces(
"isSurfaceFace.obj", isSurfaceFace);
729 bitSet isFrontFace(frontFaces(isSurfaceFace, locationFaces, isHoleCell));
733 label surfaceDist =
gMax(surfaceCellDist);
739 while (surfaceDist >= 0)
749 for (
const label facei : isFrontFace)
751 const label own = faceOwner[facei];
754 const label ownDist = surfaceCellDist[own];
755 if (ownDist >= surfaceDist)
758 isHoleCell.unset(own);
761 const cell& cFaces = mesh_.cells()[own];
763 const unsigned int mask = locationFaces[facei];
764 for (
const label fi : cFaces)
767 locationFaces[fi] |= mask;
771 else if (mesh_.isInternalFace(facei))
773 const label nei = faceNeighbour[facei];
774 const label neiDist = surfaceCellDist[nei];
776 if (isHoleCell[nei] && neiDist >= surfaceDist)
779 isHoleCell[nei] =
false;
782 const cell& cFaces = mesh_.cells()[nei];
784 const unsigned int mask = locationFaces[facei];
785 for (
const label fi : cFaces)
788 locationFaces[fi] |= mask;
799 Pout<<
"holeToFace::findClosure :"
800 <<
" surfaceDist:" << surfaceDist
801 <<
" front:" << isFrontFace.count()
802 <<
" nChangedCells:" << nChanged
827 isFrontFace = frontFaces(isSurfaceFace, locationFaces, isHoleCell);
853 bitSet isCommonFace(mesh_.nFaces());
854 forAll(locationFaces, facei)
856 unsigned int masks = locationFaces[facei];
863 isCommonFace.set(facei);
869 for (
const label facei : isCommonFace)
871 if (surfaceFaceDist[facei] == 0)
873 isCommonFace.unset(facei);
879 Pout<<
"holeToFace::findClosure :"
880 <<
" closure faces:" << isCommonFace.count() <<
endl;
887Foam::bitSet Foam::holeToFace::erodeSet
889 const bitSet& isSurfaceFace,
895 const labelList& faceOwner = mesh_.faceOwner();
896 const labelList& faceNeighbour = mesh_.faceNeighbour();
898 bitSet isSetCell(mesh_.nCells());
899 for (
const label facei : isSetFace)
901 isSetCell.set(faceOwner[facei]);
902 if (mesh_.isInternalFace(facei))
904 isSetCell.set(faceNeighbour[facei]);
910 bitSet erodedSet(isSetFace);
911 for (
const label celli : isSetCell)
913 const cell& cFaces = mesh_.cells()[celli];
915 label nBlockedFaces = 0;
916 label nSurfaceFaces = 0;
917 for (
const label facei : cFaces)
919 if (erodedSet[facei])
923 else if (isSurfaceFace[facei])
929 if ((nSurfaceFaces + nBlockedFaces) == cFaces.size())
932 for (
const label facei : cFaces)
934 erodedSet.unset(facei);
946 for (
const label celli : isSetCell)
948 const cell& cFaces = mesh_.cells()[celli];
950 label nBlockedFaces = 0;
951 for (
const label facei : cFaces)
953 if (erodedSet[facei])
958 if (nBlockedFaces >= cFaces.size()-2)
960 for (
const label facei : cFaces)
962 erodedSet.flip(facei);
975 Pout<<
"holeToFace::erodeSet :"
976 <<
" starting set:" << isSetFace.count()
977 <<
" eroded set:" << erodedSet.count() <<
endl;
988 const bitSet& isSurfaceFace,
994 forAll(zonePoints_, zonei)
996 const pointField& zoneLocations = zonePoints_[zonei];
997 labelList& zoneCells = locationCells[zonei];
1001 const label celli = mesh_.findCell(zoneLocations[i]);
1002 zoneCells[i] = celli;
1008 const cell& cFaces = mesh_.cells()[celli];
1009 bool hasUnblocked =
false;
1010 for (
const label facei : cFaces)
1012 if (!isSurfaceFace[facei])
1014 hasUnblocked =
true;
1022 <<
"problem : location:" << zoneLocations[i]
1023 <<
" in zone:" << zonei
1024 <<
" is found in cell at:" << celli
1025 << mesh_.cellCentres()[celli]
1026 <<
" which is completely surrounded by blocked faces"
1033 bitSet isClosingFace
1045 isClosingFace = erodeSet(isSurfaceFace, isClosingFace);
1050 writeFaces(
"isClosingFace.obj", isClosingFace);
1054 for (
const label facei : isClosingFace)
1056 addOrDelete(set, facei,
add);
1061Foam::List<Foam::pointField> Foam::holeToFace::expand(
const pointField&
pts)
1067 onePt.setSize(1,
pts[i]);
1085 zonePoints_(zonePoints),
1086 blockedFaceNames_(blockedFaceNames),
1087 blockedCellNames_(blockedCellNames),
1100 blockedFaceNames_(),
1101 blockedCellNames_(),
1102 erode_(
dict.getOrDefault<bool>(
"erode", false))
1106 if (!
dict.readIfPresent(
"faceSets", blockedFaceNames_))
1108 if (dict.readEntry(
"faceSet", setName))
1110 blockedFaceNames_.resize(1);
1111 blockedFaceNames_.front() = std::move(setName);
1114 if (!
dict.readIfPresent(
"cellSets", blockedCellNames_))
1116 if (dict.readEntry(
"cellSet", setName))
1118 blockedCellNames_.resize(1);
1119 blockedCellNames_.front() = std::move(setName);
1127 const polyMesh&
mesh,
1131 topoSetFaceSource(
mesh),
1133 blockedFaceNames_(one{}, word(checkIs(is))),
1134 blockedCellNames_(),
1149 for (
const word& setName : blockedFaceNames_)
1157 if (blockedFaceNames_.size())
1159 for (
const word& setName : blockedCellNames_)
1167 isCandidateCell =
true;
1174 Info<<
" Adding all faces to disconnect regions: "
1184 Info<<
" Removing all faces to disconnect regions: "
1195 const polyMesh&
mesh,
1196 const List<pointField>& zonePoints,
1198 const globalIndex& globalBlockedFaces,
1205 if (blockedFaces.size() != globalBlockedFaces.localSize())
1208 <<
" globalBlockedFaces:" << globalBlockedFaces.localSize()
1228 faceSetSource.combine
1236 closureFaces = closureFaceSet.sortedToc();
1242 closureToBlocked.clear();
1263 const label globalBlockedi = globalBlockedFaces.toGlobal(i);
1264 const label facei = blockedFaces[i];
1269 edgeMap.insert(
edge(
f[fp],
f[nextFp]), globalBlockedi);
1279 initialEdgesInfo(2*nBndEdges);
1282 const edge&
e = edges[edgei];
1285 auto iter = edgeMap.cfind(meshE);
1291 initialEdges.append(edgei);
1292 initialEdgesInfo.append
1331 closureToBlocked.resize_nocopy(
pp.
size());
1332 closureToBlocked = -1;
1337 closureToBlocked[facei] =
allFaceInfo[facei].data();
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
void append(const T &val)
Copy append an element to the end of this list.
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
@ NO_REGISTER
Do not request registration (bool: false).
const Addr & addressing() const noexcept
The addressing used for the list.
label size() const noexcept
The number of elements in the list.
A List with indirect addressing.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
void append(const T &val)
Append an element at the end of the list.
void setSize(label n)
Alias for resize().
void clear()
Clear the list, i.e. set size to zero.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
label size() const noexcept
Number of entries.
Wave propagation of information along patch. Every iteration information goes through one layer of fa...
label nEdges() const
Number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
fileName timePath() const
Return current time path = path/timeName.
A List with indirect addressing. Like IndirectList but does not store addressing.
bool empty() const noexcept
True if List is empty (ie, size() is zero).
void size(const label n)
Older name for setAddressableSize.
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
label setMany(InputIter first, InputIter last)
Set the locations listed by the iterator range, auto-vivify entries if needed.
A collection of cell labels.
A cell is defined as a list of faces with extra functionality.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
For use with PatchEdgeFaceWave. Determines topological distance to starting edges....
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Smooth ATC in cells next to a set of patches supplied by type.
A face is a list of labels corresponding to mesh vertices.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
label toGlobal(const label proci, const label i) const
From local to global on proci.
label localSize(const label proci) const
Size of proci data.
A topoSetFaceSource to select a set of faces that closes a hole i.e. disconnects zones (specified by ...
virtual void applyToSet(const topoSetSource::setAction action, topoSet &) const
Apply specified action to the topoSet.
holeToFace(const polyMesh &mesh, const List< pointField > &zonePoints, const wordList &blockedFaceNames, const wordList &blockedCellNames, const bool erode)
Construct from components.
static autoPtr< mapDistribute > calcClosure(const polyMesh &mesh, const List< pointField > &zonePoints, const labelList &blockedFaces, const globalIndex &globalBlockedFaces, const bool erode, labelList &closureFaces, labelList &closureToBlocked)
Optional direct use to generate the set of faces and the method to.
void combine(topoSet &set, const bitSet &isBlockedFace, const bitSet &isActiveCell, const bool add) const
Optional direct use to generate a faceSet.
const Time & time() const noexcept
Return time registry.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
virtual const faceList & faces() const
Return raw faces.
virtual const pointField & points() const
Return raw points.
A patch is a list of labels that address the faces in the global face list.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
For use with FaceCellWave. Determines topological distance to starting faces. Templated on passive tr...
The topoSetFaceSource is a intermediate class for handling topoSet sources for selecting faces.
topoSetFaceSource(const polyMesh &mesh)
Construct from mesh.
Class with constructor to add usage string to table.
Base class of a source for a topoSet.
void addOrDelete(topoSet &set, const label id, const bool add) const
Add or delete id from set. Add when 'add' is true.
setAction
Enumeration defining various actions.
@ SUBTRACT
Subtract elements from current set.
@ ADD
Add elements to current set.
@ NEW
Create a new set and ADD elements to it.
bool verbose_
Output verbosity (default: true).
const polyMesh & mesh() const noexcept
Reference to the mesh.
const polyMesh & mesh_
Reference to the mesh.
static Istream & checkIs(Istream &is)
Check state of stream.
General set of labels of mesh quantity (points, cells, faces).
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
List< wallPoints > allFaceInfo(mesh_.nFaces())
const bitSet isBlockedFace(intersectedFaces())
unsigned int bit_count(UIntType x)
Count arbitrary number of bits (of an integral type).
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
List< edge > edgeList
List of edge.
List< word > wordList
List of word.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
messageStream Info
Information stream (stdout output on master, null elsewhere).
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.