87void insertDuplicateMerge
101 label otherFacei = duplicates[bFacei];
103 if (otherFacei != -1 && otherFacei > bFacei)
107 label face0 = boundaryFaces[bFacei];
108 label face1 = boundaryFaces[otherFacei];
110 label own0 = faceOwner[face0];
111 label own1 = faceOwner[face1];
116 label zoneID = faceZones.
whichZone(face0);
117 bool zoneFlip =
false;
121 const faceZone& fZone = faceZones[zoneID];
145 label zoneID = faceZones.
whichZone(face1);
146 bool zoneFlip =
false;
150 const faceZone& fZone = faceZones[zoneID];
202 faceIDs[sz++] =
pp.start()+ppi;
206 if (faceIDs.size() != sz)
228 forAll(duplicates, bFacei)
230 if (duplicates[bFacei] != -1)
232 label facei = boundaryFaces[bFacei];
233 label patchi =
patches.whichPatch(facei);
238 <<
"Duplicate face " << facei
239 <<
" is on a processorPolyPatch."
240 <<
"This is not allowed." <<
nl
242 <<
" is on patch:" <<
patches[patchi].name()
255 mesh.nBoundaryFaces()/256
258 forAll(duplicates, bFacei)
260 label otherFacei = duplicates[bFacei];
262 if (otherFacei != -1 && otherFacei > bFacei)
264 duplicateSet.insert(boundaryFaces[bFacei]);
265 duplicateSet.insert(boundaryFaces[otherFacei]);
270 <<
" duplicate faces to faceSet " << duplicateSet.objectPath()
272 duplicateSet.write();
279int main(
int argc,
char *argv[])
283 "Detect faces that share points (baffles).\n"
284 "Merge them or duplicate the points."
293 "Specify a dictionary to read actions from"
298 "Find baffles only, but do not merge or split them"
303 "Topologically split duplicate surfaces"
312 const word oldInstance =
mesh.pointsInstance();
315 const bool readDict =
args.found(
"dict");
317 const bool overwrite =
args.found(
"overwrite");
318 const bool detectOnly =
args.found(
"detectOnly");
320 if (readDict && (
split || detectOnly))
323 <<
"Use of dictionary for settings not compatible with"
324 <<
" using command line arguments for \"split\""
341 if (
dict.found(
"detect"))
343 detectPatchIDs =
patches.patchSet
348 Info<<
"Detecting baffles on " << detectPatchIDs.
size()
353 if (
dict.found(
"merge"))
355 mergePatchIDs =
patches.patchSet
360 Info<<
"Merge baffles on " << mergePatchIDs.
size()
365 if (
dict.found(
"split"))
367 splitPatchIDs =
patches.patchSet
372 Info<<
"Split baffles on " << splitPatchIDs.
size()
395 if (detectPatchIDs.
size())
397 findBaffles(
mesh, patchFaces(
mesh, detectPatchIDs));
446 if (mergePatchIDs.
size())
448 Info<<
"Merging duplicate faces" <<
nl <<
endl;
453 const labelList boundaryFaces(patchFaces(
mesh, mergePatchIDs));
459 insertDuplicateMerge(
mesh, boundaryFaces, duplicates, meshMod);
470 mesh.updateMesh(map());
473 if (map().hasMotionPoints())
475 mesh.movePoints(map().preMotionPoints());
480 mesh.setInstance(oldInstance);
487 if (splitPatchIDs.
size())
489 Info<<
"Topologically splitting duplicate surfaces"
490 <<
", i.e. duplicating points internal to duplicate surfaces"
499 sz +=
patches[splitPatchIDs[i]].nPoints();
503 bitSet isCandidate(
mesh.nPoints());
509 label pointi =
mp[mpi];
510 if (isCandidate.set(pointi))
512 candidates.
append(pointi);
529 pointDuplicator.setRefinement(
regionSide, meshMod);
540 mesh.updateMesh(map());
543 if (map().hasMotionPoints())
545 mesh.movePoints(map().preMotionPoints());
550 mesh.setInstance(oldInstance);
559 const labelList& pointMap = map().pointMap();
567 label oldPointi = pointMap[pointi];
569 nDupPerPoint[oldPointi]++;
571 if (nDupPerPoint[oldPointi] > 1)
573 dupPoints.insert(map().reversePointMap()[oldPointi]);
574 dupPoints.insert(pointi);
579 <<
" duplicated points to pointSet "
580 << dupPoints.objectPath() <<
nl <<
endl;
Field reading functions for post-processing utilities.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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.
void setCapacity(const label len)
Alter the size of the underlying storage.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
void size(const label n)
Older name for setAddressableSize.
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
static void addOption(const word &optName, const string ¶m="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
static void addNote(const string ¬e)
Add extra notes for the usage information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A subset of mesh faces organised as a primitive patch.
label whichFace(const label meshFaceID) const
The local index of the given mesh face, -1 if not in the zone.
const boolList & flipMap() const noexcept
Return face flip map.
Takes mesh with 'baffles' (= boundary faces sharing points). Determines for selected points on bounda...
static labelList findDuplicateFaces(const primitiveMesh &, const labelList &)
Helper routine to find baffles (two boundary faces using the.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Class describing modification of a face.
A patch is a list of labels that address the faces in the global face list.
Class containing data for face removal.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
Determines the 'side' for every face and connected to a singly-connected (through edges) region of fa...
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
const polyBoundaryMesh & patches
static bool split(const std::string &line, std::string &key, std::string &val)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const word dictName("faMeshDefinition")
return returnReduce(nRefine-oldNRefine, sumOp< label >())
const dimensionedScalar mp
Proton mass.
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
List< label > labelList
A List of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a).
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.