Go to the documentation of this file.
79 bool zoneFlip =
false;
122 if (!isA<processorPolyPatch>(pp))
133 addedPatchNames.found(pp.
name())
135 || isA<coupledPolyPatch>(pp)
151 Info<<
"Removing zero-sized patch " << pp.
name()
152 <<
" type " << pp.type()
153 <<
" at position " << patchi <<
endl;
162 if (isA<processorPolyPatch>(pp))
179 Info<<
"Removing empty processor patch " << pp.
name()
180 <<
" at position " << patchi <<
endl;
186 if (nAllPatches != nOldPatches)
195 Info<<
"No patches removed." <<
endl;
198 delete allPatches[i];
211 if (cpp && cpp->
owner())
213 const auto& cycPatch = *cpp;
218 OFstream str(prefix+cycPatch.name()+
".obj");
219 Pout<<
"Dumping " << cycPatch.
name()
220 <<
" faces to " << str.
name() <<
endl;
230 OFstream str(prefix+nbrPatch.name()+
".obj");
231 Pout<<
"Dumping " << nbrPatch.
name()
232 <<
" faces to " << str.
name() <<
endl;
243 OFstream str(prefix+cycPatch.name()+nbrPatch.name()+
"_match.obj");
246 Pout<<
"Dumping cyclic match as lines between face centres to "
258 str<<
"l " << vertI-1 <<
' ' << vertI <<
nl;
271 if (separation.size() == 1)
277 field[i] += separation[0];
280 else if (separation.size() ==
field.size())
284 field[i] += separation[i];
290 <<
"Sizes of field and transformation not equal. field:"
291 <<
field.size() <<
" transformation:" << separation.size()
298 template<
class CombineOp>
303 const CombineOp& cop,
304 const point& nullValue
310 <<
"Number of values " <<
points.size()
311 <<
" is not equal to the number of points in the mesh "
318 bool hasTransformation =
false;
325 for (
const label patchi : procPatches)
328 const auto& procPatch = refCast<const processorPolyPatch>(pp);
330 if (pp.
nPoints() && procPatch.owner())
333 pointField patchInfo(procPatch.nPoints(), nullValue);
335 const labelList& meshPts = procPatch.meshPoints();
336 const labelList& nbrPts = procPatch.neighbPoints();
340 label nbrPointi = nbrPts[pointi];
341 if (nbrPointi >= 0 && nbrPointi < patchInfo.size())
343 patchInfo[nbrPointi] =
points[meshPts[pointi]];
350 procPatch.neighbProcNo()
359 for (
const label patchi : procPatches)
362 const auto& procPatch = refCast<const processorPolyPatch>(pp);
364 if (pp.
nPoints() && !procPatch.owner())
373 procPatch.neighbProcNo()
375 fromNbr >> nbrPatchInfo;
378 nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
380 if (!procPatch.parallel())
382 hasTransformation =
true;
385 else if (procPatch.separated())
387 hasTransformation =
true;
388 separateList(-procPatch.separation(), nbrPatchInfo);
391 const labelList& meshPts = procPatch.meshPoints();
395 label meshPointi = meshPts[pointi];
396 points[meshPointi] = nbrPatchInfo[pointi];
407 if (cpp && cpp->
owner())
411 const auto& cycPatch = *cpp;
414 const edgeList& coupledPoints = cycPatch.coupledPoints();
415 const labelList& meshPts = cycPatch.meshPoints();
416 const labelList& nbrMeshPts = nbrPatch.meshPoints();
422 const edge&
e = coupledPoints[i];
423 label point0 = meshPts[
e[0]];
424 half0Values[i] =
points[point0];
427 if (!cycPatch.parallel())
429 hasTransformation =
true;
432 else if (cycPatch.separated())
434 hasTransformation =
true;
435 separateList(cycPatch.separation(), half0Values);
440 const edge&
e = coupledPoints[i];
441 label point1 = nbrMeshPts[
e[1]];
442 points[point1] = half0Values[i];
456 if (hasTransformation)
459 <<
"There are decomposed cyclics in this mesh with"
460 <<
" transformations." <<
endl
461 <<
"This is not supported. The result will be incorrect"
492 int main(
int argc,
char *argv[])
496 "Create patches out of selected boundary faces, which are either"
497 " from existing patches or from a faceSet"
506 "Write obj files showing the cyclic matching process"
514 const word meshRegionName =
517 const bool overwrite =
args.
found(
"overwrite");
521 const bool writeObj =
args.
found(
"writeObj");
532 const bool pointSync(
dict.
get<
bool>(
"pointSync"));
542 dumpCyclicMatch(
"initial_",
mesh);
558 if (patchSources.size())
570 if (!isA<processorPolyPatch>(pp))
582 startFacei += pp.size();
592 if (destPatchi == -1)
596 destPatchi = allPatches.size();
598 Info<<
"Adding new patch " << patchName
599 <<
" as patch " << destPatchi
600 <<
" from " << patchDict <<
endl;
602 patchDict.set(
"nFaces", 0);
603 patchDict.set(
"startFace", startFacei);
619 Info<<
"Patch '" << patchName <<
"' already exists. Only "
620 <<
"moving patch faces - type will remain the same" <<
endl;
629 if (isA<processorPolyPatch>(pp))
641 startFacei += pp.size();
668 if (destPatchi == -1)
671 <<
"patch " << patchName <<
" not added. Problem."
677 if (sourceType ==
"patches")
685 for (
const label patchi : patchSources)
689 Info<<
"Moving faces from patch " << pp.
name()
690 <<
" to patch " << destPatchi <<
endl;
702 if (!isRepatchedBoundary.set(pp.
offset()+i))
705 <<
"Face " << pp.
start() + i <<
" from patch "
706 << pp.
name() <<
" is already marked to be moved"
707 <<
" to patch " << meshMod.
region()[pp.
start() + i]
713 else if (sourceType ==
"set")
720 <<
" faces from faceSet " <<
set.name() <<
endl;
725 for (
const label facei : faceLabels)
730 <<
"Face " << facei <<
" specified in set "
732 <<
" is not an external face of the mesh." <<
endl
733 <<
"This application can only repatch existing boundary"
740 <<
"Face " << facei <<
" from set "
741 <<
set.name() <<
" is already marked to be moved"
742 <<
" to patch " << meshMod.
region()[facei]
758 <<
"Invalid source type " << sourceType <<
endl
767 Info<<
"Doing topology modification to order faces." <<
nl <<
endl;
773 dumpCyclicMatch(
"coupled_",
mesh);
779 Info<<
"Not synchronising points." <<
nl <<
endl;
796 if (pp.size() && isA<coupledPolyPatch>(pp))
799 refCast<const coupledPolyPatch>(pp);
803 Info<<
"On coupled patch " << pp.
name()
804 <<
" separation[0] was "
807 if (isA<cyclicPolyPatch>(pp) && pp.size())
810 refCast<const cyclicPolyPatch>(pp);
815 Info<<
"On cyclic translation patch " << pp.
name()
816 <<
" forcing uniform separation of "
833 Info<<
"On coupled patch " << pp.
name()
834 <<
" forcing uniform separation of "
839 Info<<
"On coupled patch " << pp.
name()
840 <<
" forcing uniform rotation of "
852 Info<<
"On coupled patch " << pp.
name()
853 <<
" forcing uniform rotation of "
859 Info<<
"Synchronising points." <<
endl;
868 point(GREAT, GREAT, GREAT)
881 Info<<
"Removing patches with no faces in them." <<
nl<<
endl;
882 filterPatches(
mesh, addedPatchNames);
887 dumpCyclicMatch(
"final_",
mesh);
void addPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
vectorField pointField
pointField is a vectorField.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
virtual const pointField & points() const
points setSize(newPointi)
virtual const vectorField & separation() const
void set(List< bool > &bools, const labelRange &range)
A class for handling words, derived from Foam::string.
A class for handling file names.
const boolList & flipMap() const noexcept
virtual bool write(const bool valid=true) const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
A PtrList of objects of type <T> with automated input and output.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
static word defaultRegion
T getOrDefault(const word &optName, const T &deflt) const
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Type gAverage(const FieldField< Field, Type > &f)
Output inter-processor communications stream.
const cyclicPolyPatch & neighbPatch() const
virtual const tensorField & forwardT() const
static std::string name(const std::string &str)
const word dictName("faMeshDefinition")
static word timeName(const scalar t, const int precision=precision_)
Direct mesh changes based on v1.3 polyTopoChange syntax.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
virtual const tensorField & reverseT() const
static void addNote(const string ¬e)
Class describing modification of a face.
virtual tmp< scalarField > movePoints(const pointField &)
const polyBoundaryMesh & boundaryMesh() const
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
Ostream & endl(Ostream &os)
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
A HashTable with keys but without contents that is similar to std::unordered_set.
Mesh consisting of general polyhedral cells.
const labelList & processorPatches() const noexcept
label nPoints() const noexcept
virtual bool parallel() const
scalar diff(const triad &A, const triad &B)
const fileName & pointsInstance() const
virtual bool separated() const
bool checkParallelSync(const bool report=false) const
void transformList(const tensor &rotTensor, UList< T > &field)
static void noFunctionObjects(bool addWithOption=false)
label nGlobalPoints() const
Generic templated field type.
A subset of mesh faces organised as a primitive patch.
A patch is a list of labels that address the faces in the global face list.
virtual const fileName & name() const
virtual const labelList & faceOwner() const
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
label nBoundaryFaces() const noexcept
const faceZoneMesh & faceZones() const noexcept
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
label findPatchID(const word &patchName, const bool allowNotFound=true) const
label max(const labelHashSet &set, label maxValue=labelMin)
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
label whichZone(const label objectIndex) const
virtual bool owner() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
const labelIOList & zoneID
label whichFace(const label globalCellID) const
errorManip< error > abort(error &err)
virtual transformType transform() const
errorManipArg< error, int > exit(error &err, const int errNo=1)
Output to file stream, using an OSstream.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
const word & name() const noexcept
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
static unsigned int defaultPrecision() noexcept
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
bool isInternalFace(const label faceIndex) const noexcept
virtual const faceList & faces() const
#define FatalErrorInFunction
label setAction(const topoAction &action)
label nInternalFaces() const noexcept
static bool & parRun() noexcept
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const vectorField & faceCentres() const
static void removeFiles(const polyMesh &mesh)
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const dimensionedScalar e
A List of wordRe with additional matching capabilities.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
bool insert(const Key &key)
const polyBoundaryMesh & patches
const fileName & instance() const noexcept
virtual autoPtr< polyPatch > clone(const labelList &faceCells) const
const labelList & sharedPointLabels() const
const labelList & sharedPointAddr() const
Input inter-processor communications stream.
const word & name() const noexcept
const DynamicList< label > & region() const
const globalMeshData & globalData() const
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
vector point
Point is a vector.
static void addOption(const word &optName, const string ¶m="", const string &usage="", bool advanced=false)
Foam::argList args(argc, argv)
static void removeFiles(const polyMesh &)
#define WarningInFunction
const vector & separationVector() const
Type gMax(const FieldField< Field, Type > &f)
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
bool found(const word &optName) const