createPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  createPatch
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Create patches out of selected boundary faces, which are either
35  from existing patches or from a faceSet.
36 
37  More specifically it:
38  - creates new patches (from selected boundary faces).
39  Synchronise faces on coupled patches.
40  - synchronises points on coupled boundaries
41  - remove patches with 0 faces in them
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #include "cyclicPolyPatch.H"
46 #include "syncTools.H"
47 #include "argList.H"
48 #include "Time.H"
49 #include "OFstream.H"
50 #include "meshTools.H"
51 #include "faceSet.H"
52 #include "IOPtrList.H"
53 #include "polyTopoChange.H"
54 #include "polyModifyFace.H"
55 #include "polyAddFace.H"
56 #include "wordRes.H"
57 #include "processorMeshes.H"
58 #include "IOdictionary.H"
59 #include "regionProperties.H"
60 #include "faceAreaWeightAMI2D.H"
61 #include "fvMeshTools.H"
62 #include "ReadFields.H"
63 
64 using namespace Foam;
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
71 }
72 
73 
74 word patchName(const word& name, const fvMesh& mesh0, const fvMesh& mesh1)
75 {
76  word pName(name != "none" ? name : word::null);
77  pName += mesh0.name();
78  pName += "_to_";
79  pName += mesh1.name();
80  return pName;
81 }
82 
83 
84 void matchPatchFaces
85 (
86  const word& entryName,
87  const word& AMIMethod,
88  const dictionary& AMIDict,
89  const PtrList<fvMesh>& meshes,
90 
91  const label meshi,
92  const label nSourcei,
93  const label sourcei,
94  const labelList& patchesi,
95 
96  const label meshj,
97  const label nSourcej,
98  const label sourcej,
99  const labelList& patchesj,
100 
101  DynamicList<labelList>& interfaceMesh0,
102  DynamicList<label>& interfaceSource0,
103  DynamicList<labelList>& interfacePatch0,
104  DynamicList<wordList>& interfaceNames0,
105  DynamicList<List<DynamicList<label>>>& interfaceFaces0,
106 
107  DynamicList<labelList>& interfaceMesh1,
108  DynamicList<label>& interfaceSource1,
109  DynamicList<labelList>& interfacePatch1,
110  DynamicList<wordList>& interfaceNames1,
111  DynamicList<List<DynamicList<label>>>& interfaceFaces1
112 )
113 {
114  // Now we have:
115  // - meshi, sourcei, patchesi
116  // - meshj, sourcej, patchesj
117 
118 
119  // Attempt to match patches
120  forAll(patchesi, i)
121  {
122  const auto& ppi = meshes[meshi].boundaryMesh()[patchesi[i]];
123 
124  forAll(patchesj, j)
125  {
126  const auto& ppj = meshes[meshj].boundaryMesh()[patchesj[j]];
127 
128  // Use AMI to try and find matches
129  auto AMPtr(AMIInterpolation::New(AMIMethod, AMIDict));
130 
131  AMPtr->calculate(ppi, ppj, nullptr);
132  if
133  (
134  gAverage(AMPtr->tgtWeightsSum()) > SMALL
135  || gAverage(AMPtr->srcWeightsSum()) > SMALL
136  )
137  {
138  const label inti = interfaceMesh0.size();
139 
140  Info<< "Introducing interface " << inti << " between"
141  << " mesh " << meshes[meshi].name()
142  //<< " source:" << sourcei
143  << " patch " << ppi.name()
144  << " and mesh " << meshes[meshj].name()
145  //<< " source:" << sourcej
146  << " patch " << ppj.name()
147  << endl;
148 
149  // Mesh 0
150  //~~~~~~~
151 
152  auto& intMesh0 = interfaceMesh0.emplace_back(nSourcei, -1);
153  intMesh0[sourcei] = meshi;
154 
155  interfaceSource0.push_back(sourcei);
156 
157  auto& intPatch0 = interfacePatch0.emplace_back(nSourcei, -1);
158  intPatch0[sourcei] = ppi.index();
159 
160  auto& intNames0 = interfaceNames0.emplace_back(nSourcei);
161  intNames0[sourcei] =
162  patchName(entryName, meshes[meshi], meshes[meshj]);
163 
164 
165  // Mesh 1
166  //~~~~~~~
167 
168  auto& intMesh1 = interfaceMesh1.emplace_back(nSourcej, -1);
169  intMesh1[sourcej] = meshj;
170 
171  interfaceSource1.push_back(sourcej);
172 
173  auto& intPatch1 = interfacePatch1.emplace_back(nSourcej, -1);
174  intPatch1[sourcej] = ppj.index();
175 
176  auto& intNames1 = interfaceNames1.emplace_back(nSourcej);
177  intNames1[sourcej] =
178  patchName(entryName, meshes[meshj], meshes[meshi]);
179 
180  auto& intFaces0 = interfaceFaces0.emplace_back(nSourcei);
181  DynamicList<label>& faces0 = intFaces0[sourcei];
182  faces0.setCapacity(ppi.size());
183 
184  auto& intFaces1 = interfaceFaces1.emplace_back(nSourcej);
185  DynamicList<label>& faces1 = intFaces1[sourcej];
186  faces1.setCapacity(ppj.size());
187 
188 
189  // Mark any mesh faces:
190  // - that have a weight > 0.5.
191  // - contribute to these faces (even if < 0.5)
192 
193  faces0.clear();
194  faces1.clear();
195 
196  // Marked as donor of face with high weight
197  scalarField targetMask;
198  {
199  const scalarField& weights = AMPtr->srcWeightsSum();
200  scalarField mask(weights.size(), Zero);
201  forAll(weights, facei)
202  {
203  const scalar sum = weights[facei];
204  if (sum > 0.5)
205  {
206  mask[facei] = sum;
207  }
208  }
209 
210  // Push source field mask to target
211  targetMask = AMPtr->interpolateToTarget(mask);
212  }
213 
214  scalarField sourceMask;
215  {
216  const scalarField& weights = AMPtr->tgtWeightsSum();
217  scalarField mask(weights.size(), Zero);
218  forAll(weights, facei)
219  {
220  const scalar sum = weights[facei];
221  if (sum > 0.5)
222  {
223  mask[facei] = sum;
224  }
225  }
226 
227  // Push target mask back to source
228  sourceMask = AMPtr->interpolateToSource(mask);
229  }
230 
231  {
232  const scalarField& weights = AMPtr->srcWeightsSum();
233  forAll(weights, facei)
234  {
235  if (weights[facei] > 0.5 || sourceMask[facei] > SMALL)
236  {
237  faces0.push_back(ppi.start()+facei);
238  }
239  }
240  }
241  {
242  const scalarField& weights = AMPtr->tgtWeightsSum();
243  forAll(weights, facei)
244  {
245  if (weights[facei] > 0.5 || targetMask[facei] > SMALL)
246  {
247  faces1.push_back(ppj.start()+facei);
248  }
249  }
250  }
251 
252  faces0.shrink();
253  faces1.shrink();
254  }
255  }
256  }
257 }
258 
259 
260 void matchPatchFaces
261 (
262  const bool includeOwn,
263  const PtrList<fvMesh>& meshes,
264  const List<DynamicList<label>>& interRegionSources,
265  const List<wordList>& patchNames,
266  const labelListListList& matchPatchIDs,
267  const List<wordList>& AMIMethods,
268  List<PtrList<dictionary>> patchInfoDicts,
269 
270  DynamicList<labelList>& interfaceMesh0,
271  DynamicList<label>& interfaceSource0,
272  DynamicList<labelList>& interfacePatch0,
273  DynamicList<List<DynamicList<label>>>& interfaceFaces0,
274  DynamicList<wordList>& interfaceNames0,
275 
276  DynamicList<labelList>& interfaceMesh1,
277  DynamicList<label>& interfaceSource1,
278  DynamicList<labelList>& interfacePatch1,
279  DynamicList<List<DynamicList<label>>>& interfaceFaces1,
280  DynamicList<wordList>& interfaceNames1
281 )
282 {
283  // Add approximate matches
284  forAll(meshes, meshi)
285  {
286  const labelListList& allPatchesi = matchPatchIDs[meshi];
287 
288  const label meshjStart(includeOwn ? meshi : meshi+1);
289 
290  for (label meshj = meshjStart; meshj < meshes.size(); meshj++)
291  {
292  const labelListList& allPatchesj = matchPatchIDs[meshj];
293 
294  for (const label sourcei : interRegionSources[meshi])
295  {
296  const labelList& patchesi = allPatchesi[sourcei];
297 
298  for (const label sourcej : interRegionSources[meshj])
299  {
300  const labelList& patchesj = allPatchesj[sourcej];
301 
302  // Now we have:
303  // - meshi, sourcei, patchesi
304  // - meshj, sourcej, patchesj
305 
306  matchPatchFaces
307  (
308  patchNames[meshi][sourcei],
309  (
310  AMIMethods[meshi][sourcei].size()
311  ? AMIMethods[meshi][sourcei]
312  : faceAreaWeightAMI2D::typeName
313  ),
314  patchInfoDicts[meshi][sourcei],
315  meshes,
316 
317  meshi,
318  allPatchesi.size(), // nSourcei,
319  sourcei,
320  patchesi,
321 
322  meshj,
323  allPatchesj.size(), // nSourcej,
324  sourcej,
325  patchesj,
326 
327  interfaceMesh0,
328  interfaceSource0,
329  interfacePatch0,
330  interfaceNames0,
331  interfaceFaces0,
332 
333  interfaceMesh1,
334  interfaceSource1,
335  interfacePatch1,
336  interfaceNames1,
337  interfaceFaces1
338  );
339  }
340  }
341  }
342  }
343 }
344 
345 
346 void changePatchID
347 (
348  const bool modify,
349  const fvMesh& mesh,
350  const label faceID,
351  const label patchID,
352  polyTopoChange& meshMod
353 )
354 {
355  const label zoneID = mesh.faceZones().whichZone(faceID);
356 
357  bool zoneFlip = false;
358 
359  if (zoneID >= 0)
360  {
361  const faceZone& fZone = mesh.faceZones()[zoneID];
362 
363  zoneFlip = fZone.flipMap()[fZone.whichFace(faceID)];
364  }
365 
366  if (modify)
367  {
368  meshMod.setAction
369  (
371  (
372  mesh.faces()[faceID], // face
373  faceID, // face ID
374  mesh.faceOwner()[faceID], // owner
375  -1, // neighbour
376  false, // flip flux
377  patchID, // patch ID
378  false, // remove from zone
379  zoneID, // zone ID
380  zoneFlip // zone flip
381  )
382  );
383  }
384  else
385  {
386  meshMod.setAction
387  (
389  (
390  mesh.faces()[faceID], // modified face
391  mesh.faceOwner()[faceID], // owner
392  -1, // neighbour
393  -1, // master point
394  -1, // master edge
395  faceID, // master face
396  false, // face flip
397  patchID, // patch for face
398  zoneID, // zone for face
399  zoneFlip // face flip in zone
400  )
401  );
402  }
403 }
404 
405 
406 void changePatchID
407 (
408  const fvMesh& mesh,
409  const labelList& faceLabels,
410  const label patchID,
411  bitSet& isRepatchedBoundary,
412  polyTopoChange& meshMod
413 )
414 {
415  for (const label facei : faceLabels)
416  {
417  if (mesh.isInternalFace(facei))
418  {
420  << "Face " << facei
421  << " is not an external face of the mesh." << endl
422  << "This application can only repatch"
423  << " existing boundary faces." << exit(FatalError);
424  }
425 
426  const bool isFirst =
427  isRepatchedBoundary.set(facei-mesh.nInternalFaces());
428  if (!isFirst)
429  {
430  static label nWarnings = 0;
431  if (nWarnings == 0)
432  {
433  const label newPatchi = meshMod.region()[facei];
434  //FatalErrorInFunction
436  << "Face " << facei
437  << " at " << mesh.faceCentres()[facei]
438  << " marked for patch " << patchID
439  << " name " << mesh.boundaryMesh()[patchID].name()
440  << " is already marked for patch " << newPatchi
441  << " name " << mesh.boundaryMesh()[newPatchi].name()
442  << ". Creating duplicate face. Suppressing further warnings"
443  //<< exit(FatalError);
444  << endl;
445  }
446  nWarnings++;
447  }
448 
449  changePatchID(isFirst, mesh, facei, patchID, meshMod);
450  }
451 }
452 
453 
454 // Dump for all patches the current match
455 void dumpCyclicMatch(const fileName& prefix, const polyMesh& mesh)
456 {
457  for (const polyPatch& pp : mesh.boundaryMesh())
458  {
459  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(pp);
460 
461  if (cpp && cpp->owner())
462  {
463  const auto& cycPatch = *cpp;
464  const auto& nbrPatch = cycPatch.neighbPatch();
465 
466  // Dump patches
467  {
468  OFstream str(prefix+cycPatch.name()+".obj");
469  Pout<< "Dumping " << cycPatch.name()
470  << " faces to " << str.name() << endl;
472  (
473  str,
474  cycPatch,
475  cycPatch.points()
476  );
477  }
478 
479  {
480  OFstream str(prefix+nbrPatch.name()+".obj");
481  Pout<< "Dumping " << nbrPatch.name()
482  << " faces to " << str.name() << endl;
484  (
485  str,
486  nbrPatch,
487  nbrPatch.points()
488  );
489  }
490 
491 
492  // Lines between corresponding face centres
493  OFstream str(prefix+cycPatch.name()+nbrPatch.name()+"_match.obj");
494  label vertI = 0;
495 
496  Pout<< "Dumping cyclic match as lines between face centres to "
497  << str.name() << endl;
498 
499  forAll(cycPatch, facei)
500  {
501  const point& fc0 = mesh.faceCentres()[cycPatch.start()+facei];
502  meshTools::writeOBJ(str, fc0);
503  vertI++;
504  const point& fc1 = mesh.faceCentres()[nbrPatch.start()+facei];
505  meshTools::writeOBJ(str, fc1);
506  vertI++;
507 
508  str<< "l " << vertI-1 << ' ' << vertI << nl;
509  }
510  }
511  }
512 }
513 
514 
515 void separateList
516 (
517  const vectorField& separation,
519 )
520 {
521  if (separation.size() == 1)
522  {
523  // Single value for all.
524 
525  forAll(field, i)
526  {
527  field[i] += separation[0];
528  }
529  }
530  else if (separation.size() == field.size())
531  {
532  forAll(field, i)
533  {
534  field[i] += separation[i];
535  }
536  }
537  else
538  {
540  << "Sizes of field and transformation not equal. field:"
541  << field.size() << " transformation:" << separation.size()
542  << abort(FatalError);
543  }
544 }
545 
546 
547 // Synchronise points on both sides of coupled boundaries.
548 template<class CombineOp>
549 void syncPoints
550 (
551  const polyMesh& mesh,
553  const CombineOp& cop,
554  const point& nullValue
555 )
556 {
557  if (points.size() != mesh.nPoints())
558  {
560  << "Number of values " << points.size()
561  << " is not equal to the number of points in the mesh "
562  << mesh.nPoints() << abort(FatalError);
563  }
564 
566 
567  // Is there any coupled patch with transformation?
568  bool hasTransformation = false;
569 
570  if (UPstream::parRun())
571  {
572  const labelList& procPatches = mesh.globalData().processorPatches();
573 
574  // Send
575  for (const label patchi : procPatches)
576  {
577  const polyPatch& pp = patches[patchi];
578  const auto& procPatch = refCast<const processorPolyPatch>(pp);
579 
580  if (pp.nPoints() && procPatch.owner())
581  {
582  // Get data per patchPoint in neighbouring point numbers.
583  pointField patchInfo(procPatch.nPoints(), nullValue);
584 
585  const labelList& meshPts = procPatch.meshPoints();
586  const labelList& nbrPts = procPatch.neighbPoints();
587 
588  forAll(nbrPts, pointi)
589  {
590  label nbrPointi = nbrPts[pointi];
591  if (nbrPointi >= 0 && nbrPointi < patchInfo.size())
592  {
593  patchInfo[nbrPointi] = points[meshPts[pointi]];
594  }
595  }
596 
597  OPstream toNbr
598  (
600  procPatch.neighbProcNo()
601  );
602  toNbr << patchInfo;
603  }
604  }
605 
606 
607  // Receive and set.
608 
609  for (const label patchi : procPatches)
610  {
611  const polyPatch& pp = patches[patchi];
612  const auto& procPatch = refCast<const processorPolyPatch>(pp);
613 
614  if (pp.nPoints() && !procPatch.owner())
615  {
616  // We do not know the number of points on the other side
617  // so cannot use UIPstream::read
618 
619  pointField nbrPatchInfo;
620  IPstream::recv(nbrPatchInfo, procPatch.neighbProcNo());
621 
622  // Null any value which is not on neighbouring processor
623  nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
624 
625  if (!procPatch.parallel())
626  {
627  hasTransformation = true;
628  transformList(procPatch.forwardT(), nbrPatchInfo);
629  }
630  else if (procPatch.separated())
631  {
632  hasTransformation = true;
633  separateList(-procPatch.separation(), nbrPatchInfo);
634  }
635 
636  const labelList& meshPts = procPatch.meshPoints();
637 
638  forAll(meshPts, pointi)
639  {
640  label meshPointi = meshPts[pointi];
641  points[meshPointi] = nbrPatchInfo[pointi];
642  }
643  }
644  }
645  }
646 
647  // Do the cyclics.
648  for (const polyPatch& pp : patches)
649  {
650  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(pp);
651 
652  if (cpp && cpp->owner())
653  {
654  // Owner does all.
655 
656  const auto& cycPatch = *cpp;
657  const auto& nbrPatch = cycPatch.neighbPatch();
658 
659  const edgeList& coupledPoints = cycPatch.coupledPoints();
660  const labelList& meshPts = cycPatch.meshPoints();
661  const labelList& nbrMeshPts = nbrPatch.meshPoints();
662 
663  pointField half0Values(coupledPoints.size());
664 
665  forAll(coupledPoints, i)
666  {
667  const edge& e = coupledPoints[i];
668  label point0 = meshPts[e[0]];
669  half0Values[i] = points[point0];
670  }
671 
672  if (!cycPatch.parallel())
673  {
674  hasTransformation = true;
675  transformList(cycPatch.reverseT(), half0Values);
676  }
677  else if (cycPatch.separated())
678  {
679  hasTransformation = true;
680  separateList(cycPatch.separation(), half0Values);
681  }
682 
683  forAll(coupledPoints, i)
684  {
685  const edge& e = coupledPoints[i];
686  label point1 = nbrMeshPts[e[1]];
687  points[point1] = half0Values[i];
688  }
689  }
690  }
691 
692  //- Note: hasTransformation is only used for warning messages so
693  // reduction not strictly necessary.
694  //Pstream::reduceOr(hasTransformation);
695 
696  // Synchronize multiple shared points.
697  const globalMeshData& pd = mesh.globalData();
698 
699  if (pd.nGlobalPoints() > 0)
700  {
701  if (hasTransformation)
702  {
704  << "There are decomposed cyclics in this mesh with"
705  << " transformations." << endl
706  << "This is not supported. The result will be incorrect"
707  << endl;
708  }
709 
710 
711  // Values on shared points.
712  pointField sharedPts(pd.nGlobalPoints(), nullValue);
713 
714  forAll(pd.sharedPointLabels(), i)
715  {
716  label meshPointi = pd.sharedPointLabels()[i];
717  // Fill my entries in the shared points
718  sharedPts[pd.sharedPointAddr()[i]] = points[meshPointi];
719  }
720 
721  // Combine - globally consistent
722  Pstream::listCombineReduce(sharedPts, cop);
723 
724  // Now we will all have the same information. Merge it back with
725  // my local information.
726  forAll(pd.sharedPointLabels(), i)
727  {
728  label meshPointi = pd.sharedPointLabels()[i];
729  points[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
730  }
731  }
732 }
733 
734 
735 
736 int main(int argc, char *argv[])
737 {
739  (
740  "Create patches out of selected boundary faces, which are either"
741  " from existing patches or from a faceSet"
742  );
743 
744  #include "addOverwriteOption.H"
745  #include "addAllRegionOptions.H"
746 
747  argList::addOption("dict", "file", "Alternative createPatchDict");
749  (
750  "writeObj",
751  "Write obj files showing the cyclic matching process"
752  );
753 
754  argList::noFunctionObjects(); // Never use function objects
755 
756  #include "setRootCase.H"
757  #include "createTime.H"
758  #include "getAllRegionOptions.H"
759 
760  const bool overwrite = args.found("overwrite");
761 
762  #include "createNamedMeshes.H"
763 
764  const bool writeObj = args.found("writeObj");
765 
766 
767  // Read dictionaries and extract various
768  wordList oldInstances(meshes.size());
771  List<labelListList> matchPatchIDs(meshes.size());
772  List<wordList> matchMethods(meshes.size());
773  List<PtrList<dictionary>> patchInfoDicts(meshes.size());
774 
775  List<DynamicList<label>> interRegionSources(meshes.size());
776  forAll(meshes, meshi)
777  {
778  fvMesh& mesh = meshes[meshi];
780 
781  // If running parallel check same patches everywhere
783 
784  oldInstances[meshi] = mesh.pointsInstance();
785 
786  const word dictName("createPatchDict");
787  #include "setSystemMeshDictionaryIO.H"
788  Info<< "Reading " << dictIO.instance()/dictIO.name() << nl << endl;
789 
790  dicts.set(meshi, new IOdictionary(dictIO));
791 
792  const auto& dict = dicts[meshi];
793  PtrList<dictionary> patchSources(dict.lookup("patches"));
794  patchNames[meshi].setSize(patchSources.size());
795  matchPatchIDs[meshi].setSize(patchSources.size());
796  matchMethods[meshi].setSize(patchSources.size());
797  patchInfoDicts[meshi].setSize(patchSources.size());
798 
799  // Read patch construct info from dictionary
800  forAll(patchSources, sourcei)
801  {
802  const auto& pDict = patchSources[sourcei];
803  patchNames[meshi][sourcei] = pDict.getOrDefault<word>
804  (
805  "name",
806  word::null,
808  );
809 
810  patchInfoDicts[meshi].set
811  (
812  sourcei,
813  new dictionary(pDict.subDict("patchInfo"))
814  );
815  dictionary& patchDict = patchInfoDicts[meshi][sourcei];
816  if (patchDict.found("AMIMethod"))
817  {
818  matchMethods[meshi][sourcei] = patchDict.get<word>("AMIMethod");
819  // Disable full matching since we're trying to use AMIMethod to
820  // find out actual overlap
821  patchDict.add("requireMatch", false);
822  }
823 
824  wordRes matchNames;
825  if (pDict.readIfPresent("patches", matchNames, keyType::LITERAL))
826  {
827  matchPatchIDs[meshi][sourcei] =
828  patches.patchSet(matchNames).sortedToc();
829  }
830 
831  if (pDict.get<word>("constructFrom") == "autoPatch")
832  {
833  interRegionSources[meshi].push_back(sourcei);
834  }
835  }
836  }
837 
838 
839 
840  // Determine matching faces. This is an interface between two
841  // sets of faces:
842  // - originating mesh
843  // - originating patch
844  // - originating (mesh)facelabels
845  // It matches all mesh against each other. Lower numbered mesh gets
846  // postfix 0, higher numbered mesh postfix 1.
847 
848  // Per interface, per patchSource:
849  // 1. the lower numbered mesh
850  DynamicList<labelList> interfaceMesh0;
851  // 1b. the source index (i.e. the patch dictionary)
852  DynamicList<label> interfaceSource0;
853  // 2. the patch on the interfaceMesh0
854  DynamicList<labelList> interfacePatch0;
855  // 3. the facelabels on the interfaceMesh0
856  DynamicList<List<DynamicList<label>>> interfaceFaces0;
857  // 4. generated interface name
858  DynamicList<wordList> interfaceNames0;
859 
860  // Same for the higher numbered mesh
861  DynamicList<labelList> interfaceMesh1;
862  DynamicList<label> interfaceSource1;
863  DynamicList<labelList> interfacePatch1;
864  DynamicList<List<DynamicList<label>>> interfaceFaces1;
865  DynamicList<wordList> interfaceNames1;
866  {
867  // Whether to match to patches in own mesh
868  const bool includeOwn = (meshes.size() == 1);
869 
870  matchPatchFaces
871  (
872  includeOwn,
873  meshes,
874  interRegionSources,
875  patchNames,
876  matchPatchIDs,
877  matchMethods,
878  patchInfoDicts,
879 
880  interfaceMesh0,
881  interfaceSource0,
882  interfacePatch0,
883  interfaceFaces0,
884  interfaceNames0,
885 
886  interfaceMesh1,
887  interfaceSource1,
888  interfacePatch1,
889  interfaceFaces1,
890  interfaceNames1
891  );
892  }
893 
894 
895 
896  // Read fields
907 
908  {
909  forAll(meshes, meshi)
910  {
911  const fvMesh& mesh = meshes[meshi];
912 
913  bool noFields = true;
914  for (const auto& d : patchInfoDicts[meshi])
915  {
916  if (d.found("patchFields"))
917  {
918  noFields = false;
919  }
920  }
921 
922  if (!noFields)
923  {
924  // Read objects in time directory
925  IOobjectList objects(mesh, runTime.timeName());
926 
927  // Read volume fields.
928  ReadFields(mesh, objects, vsFlds[meshi]);
929  ReadFields(mesh, objects, vvFlds[meshi]);
930  ReadFields(mesh, objects, vstFlds[meshi]);
931  ReadFields(mesh, objects, vsymtFlds[meshi]);
932  ReadFields(mesh, objects, vtFlds[meshi]);
933 
934  // Read surface fields.
935  ReadFields(mesh, objects, ssFlds[meshi]);
936  ReadFields(mesh, objects, svFlds[meshi]);
937  ReadFields(mesh, objects, sstFlds[meshi]);
938  ReadFields(mesh, objects, ssymtFlds[meshi]);
939  ReadFields(mesh, objects, stFlds[meshi]);
940  }
941  }
942  }
943 
944 
945 
946  // Loop over all regions
947 
948  forAll(meshes, meshi)
949  {
950  fvMesh& mesh = meshes[meshi];
951 
952  Info<< "\n\nAdding patches to mesh " << mesh.name() << nl << endl;
953 
955  const dictionary& dict = dicts[meshi];
956 
957  if (writeObj)
958  {
959  dumpCyclicMatch("initial_", mesh);
960  }
961 
962  // Read patch construct info from dictionary
963  PtrList<dictionary> patchSources(dict.lookup("patches"));
964 
965 
966  // 1. Add all new patches
967  // ~~~~~~~~~~~~~~~~~~~~~~
968 
969  forAll(patchSources, sourcei)
970  {
971  const dictionary& dict = patchSources[sourcei];
972  const word sourceType(dict.get<word>("constructFrom"));
973  const word patchName
974  (
976  (
977  "name",
978  word::null,
980  )
981  );
982 
983  dictionary patchDict(patchInfoDicts[meshi][sourcei]);
984  patchDict.set("nFaces", 0);
985  patchDict.set("startFace", 0); // Gets overwritten
986 
987 
988  if (sourceType == "autoPatch")
989  {
990  // Special : automatically create the necessary inter-region
991  // coupling patches
992 
993  forAll(interfaceMesh0, inti)
994  {
995  const labelList& allMeshes0 = interfaceMesh0[inti];
996  const wordList& allNames0 = interfaceNames0[inti];
997  const labelList& allMeshes1 = interfaceMesh1[inti];
998  const wordList& allNames1 = interfaceNames1[inti];
999 
1000  if
1001  (
1002  interfaceSource0[inti] == sourcei
1003  && allMeshes0[sourcei] == meshi
1004  )
1005  {
1006  // Current mesh is mesh0. mesh1 is the remote mesh.
1007  const label sourcej = interfaceSource1[inti];
1008  const word& patchName = allNames0[sourcei];
1009  if (patches.findPatchID(patchName) == -1)
1010  {
1011  dictionary allDict(patchDict);
1012  const auto& mesh1 = meshes[allMeshes1[sourcej]];
1013  allDict.set("sampleRegion", mesh1.name());
1014  const auto& destPatch = allNames1[sourcej];
1015  allDict.set("samplePatch", destPatch);
1016  allDict.set("neighbourPatch", destPatch);
1017 
1018  Info<< "Adding new patch " << patchName
1019  << " from " << allDict << endl;
1020 
1021  autoPtr<polyPatch> ppPtr
1022  (
1024  (
1025  patchName,
1026  allDict,
1027  0, // overwritten
1028  patches
1029  )
1030  );
1032  (
1033  mesh,
1034  ppPtr(),
1035  patchDict.subOrEmptyDict("patchFields"),
1037  true
1038  );
1039  }
1040  }
1041  }
1042 
1043  forAll(interfaceMesh1, inti)
1044  {
1045  const labelList& allMeshes0 = interfaceMesh0[inti];
1046  const wordList& allNames0 = interfaceNames0[inti];
1047  const labelList& allMeshes1 = interfaceMesh1[inti];
1048  const wordList& allNames1 = interfaceNames1[inti];
1049 
1050  if
1051  (
1052  interfaceSource1[inti] == sourcei
1053  && allMeshes1[sourcei] == meshi
1054  )
1055  {
1056  // Current mesh is mesh1. mesh0 is the remote mesh.
1057 
1058  const label sourcej = interfaceSource0[inti];
1059  const word& patchName = allNames1[sourcei];
1060  if (patches.findPatchID(patchName) == -1)
1061  {
1062  dictionary allDict(patchDict);
1063  const auto& destPatch = allNames0[sourcej];
1064  const auto& mesh0 = meshes[allMeshes0[sourcej]];
1065  allDict.set("sampleRegion", mesh0.name());
1066  allDict.set("samplePatch", destPatch);
1067  allDict.set("neighbourPatch", destPatch);
1068 
1069  Info<< "Adding new patch " << patchName
1070  << " from " << allDict << endl;
1071 
1072  autoPtr<polyPatch> ppPtr
1073  (
1075  (
1076  patchName,
1077  allDict,
1078  0, // overwritten
1079  patches
1080  )
1081  );
1083  (
1084  mesh,
1085  ppPtr(),
1086  patchDict.subOrEmptyDict("patchFields"),
1088  true
1089  );
1090  }
1091  }
1092  }
1093  }
1094  else
1095  {
1096  if (patches.findPatchID(patchName) == -1)
1097  {
1098  Info<< "Adding new patch " << patchName
1099  << " from " << patchDict << endl;
1100 
1101  autoPtr<polyPatch> ppPtr
1102  (
1104  (
1105  patchName,
1106  patchDict,
1107  0, // overwritten
1108  patches
1109  )
1110  );
1112  (
1113  mesh,
1114  ppPtr(),
1115  patchDict.subOrEmptyDict("patchFields"),
1117  true
1118  );
1119  }
1120  }
1121  }
1122 
1123  Info<< endl;
1124  }
1125 
1126 
1127  // 2. Repatch faces
1128  // ~~~~~~~~~~~~~~~~
1129 
1130  forAll(meshes, meshi)
1131  {
1132  fvMesh& mesh = meshes[meshi];
1133 
1134  Info<< "\n\nRepatching mesh " << mesh.name() << nl << endl;
1135 
1136 
1138  const dictionary& dict = dicts[meshi];
1139 
1140  // Whether to synchronise points
1141  const bool pointSync(dict.get<bool>("pointSync"));
1142 
1143  // Read patch construct info from dictionary
1144  PtrList<dictionary> patchSources(dict.lookup("patches"));
1145 
1146 
1147  polyTopoChange meshMod(mesh);
1148 
1149  // Mark all repatched faces. This makes sure that the faces to repatch
1150  // do not overlap
1151  bitSet isRepatchedBoundary(mesh.nBoundaryFaces());
1152 
1153  forAll(patchSources, sourcei)
1154  {
1155  const dictionary& dict = patchSources[sourcei];
1156  const word patchName
1157  (
1159  (
1160  "name",
1161  word::null,
1163  )
1164  );
1165  const word sourceType(dict.get<word>("constructFrom"));
1166 
1167  if (sourceType == "autoPatch")
1168  {
1169  forAll(interfaceMesh0, inti)
1170  {
1171  const labelList& allMeshes0 = interfaceMesh0[inti];
1172  const wordList& allNames0 = interfaceNames0[inti];
1173  const auto& allFaces0 = interfaceFaces0[inti];
1174 
1175  const labelList& allMeshes1 = interfaceMesh1[inti];
1176  const wordList& allNames1 = interfaceNames1[inti];
1177  const auto& allFaces1 = interfaceFaces1[inti];
1178 
1179  if
1180  (
1181  interfaceSource0[inti] == sourcei
1182  && allMeshes0[sourcei] == meshi
1183  )
1184  {
1185  // Current mesh is mesh0. mesh1 is the remote mesh.
1186 
1187  const label destPatchi =
1188  patches.findPatchID(allNames0[sourcei], false);
1189 
1190  //const auto& mesh1 =
1191  // meshes[allMeshes1[interfaceSource1[inti]]];
1192  //Pout<< "Matched mesh:" << mesh.name()
1193  // << " to mesh:" << mesh1.name()
1194  // << " through:" << allNames0[sourcei] << endl;
1195 
1196  changePatchID
1197  (
1198  mesh,
1199  allFaces0[sourcei],
1200  destPatchi,
1201  isRepatchedBoundary,
1202  meshMod
1203  );
1204  }
1205  if
1206  (
1207  interfaceSource1[inti] == sourcei
1208  && allMeshes1[sourcei] == meshi
1209  )
1210  {
1211  // Current mesh is mesh1. mesh0 is the remote mesh.
1212 
1213  const label destPatchi =
1214  patches.findPatchID(allNames1[sourcei], false);
1215 
1216  //const auto& mesh0 =
1217  // meshes[allMeshes0[interfaceSource0[inti]]];
1218  //Pout<< "Matched mesh:" << mesh.name()
1219  // << " to mesh:" << mesh0.name()
1220  // << " through:" << allNames1[sourcei] << endl;
1221 
1222  changePatchID
1223  (
1224  mesh,
1225  allFaces1[sourcei],
1226  destPatchi,
1227  isRepatchedBoundary,
1228  meshMod
1229  );
1230  }
1231  }
1232  }
1233  else if (sourceType == "patches")
1234  {
1235  const label destPatchi = patches.findPatchID(patchName, false);
1236  labelHashSet patchSources
1237  (
1238  patches.patchSet(dict.get<wordRes>("patches"))
1239  );
1240 
1241  // Repatch faces of the patches.
1242  for (const label patchi : patchSources)
1243  {
1244  const polyPatch& pp = patches[patchi];
1245 
1246  Info<< "Moving faces from patch " << pp.name()
1247  << " to patch " << destPatchi << endl;
1248 
1249  changePatchID
1250  (
1251  mesh,
1252  pp.start()+identity(pp.size()),
1253  destPatchi,
1254  isRepatchedBoundary,
1255  meshMod
1256  );
1257  }
1258  }
1259  else if (sourceType == "set")
1260  {
1261  const label destPatchi = patches.findPatchID(patchName, false);
1262  const word setName(dict.get<word>("set"));
1263 
1264  faceSet set(mesh, setName);
1265 
1266  Info<< "Read " << returnReduce(set.size(), sumOp<label>())
1267  << " faces from faceSet " << set.name() << endl;
1268 
1269  // Sort (since faceSet contains faces in arbitrary order)
1271 
1272  changePatchID
1273  (
1274  mesh,
1275  faceLabels,
1276  destPatchi,
1277  isRepatchedBoundary,
1278  meshMod
1279  );
1280  }
1281  else
1282  {
1284  << "Invalid source type " << sourceType << endl
1285  << "Valid source types are 'patches' 'set'"
1286  << exit(FatalError);
1287  }
1288  }
1289 
1290 
1291  // Change mesh, use inflation to reforce calculation of transformation
1292  // tensors.
1293  Info<< "Doing topology modification to order faces." << nl << endl;
1294  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, true);
1295 
1296  if (map().hasMotionPoints())
1297  {
1298  mesh.movePoints(map().preMotionPoints());
1299  }
1300  else
1301  {
1302  // Force calculation of transformation tensors
1304  }
1305 
1306 
1307  // Update fields
1308  mesh.updateMesh(map());
1309 
1310 
1311  // Update numbering pointing to meshi
1312  const auto& oldToNew = map().reverseFaceMap();
1313  forAll(interfaceMesh0, inti)
1314  {
1315  const labelList& allMeshes0 = interfaceMesh0[inti];
1316  forAll(allMeshes0, sourcei)
1317  {
1318  if (allMeshes0[sourcei] == meshi)
1319  {
1320  inplaceRenumber(oldToNew, interfaceFaces0[inti][sourcei]);
1321  }
1322  }
1323  }
1324  forAll(interfaceMesh1, inti)
1325  {
1326  const labelList& allMeshes1 = interfaceMesh1[inti];
1327  forAll(allMeshes1, sourcei)
1328  {
1329  if (allMeshes1[sourcei] == meshi)
1330  {
1331  inplaceRenumber(oldToNew, interfaceFaces1[inti][sourcei]);
1332  }
1333  }
1334  }
1335 
1336 
1337  if (writeObj)
1338  {
1339  dumpCyclicMatch("coupled_", mesh);
1340  }
1341 
1342  // Synchronise points.
1343  if (!pointSync)
1344  {
1345  Info<< "Not synchronising points." << nl << endl;
1346  }
1347  else
1348  {
1349  Info<< "Synchronising points." << nl << endl;
1350 
1351  // This is a bit tricky. Both normal and position might be out and
1352  // current separation also includes the normal
1353  // ( separation_ = (nf&(Cr - Cf))*nf ).
1354 
1355  // For cyclic patches:
1356  // - for separated ones use user specified offset vector
1357 
1358  forAll(mesh.boundaryMesh(), patchi)
1359  {
1360  const polyPatch& pp = mesh.boundaryMesh()[patchi];
1361 
1362  if (pp.size() && isA<coupledPolyPatch>(pp))
1363  {
1364  const coupledPolyPatch& cpp =
1365  refCast<const coupledPolyPatch>(pp);
1366 
1367  if (cpp.separated())
1368  {
1369  Info<< "On coupled patch " << pp.name()
1370  << " separation[0] was "
1371  << cpp.separation()[0] << endl;
1372 
1373  if (isA<cyclicPolyPatch>(pp) && pp.size())
1374  {
1375  const auto& cycpp =
1376  refCast<const cyclicPolyPatch>(pp);
1377 
1378  if
1379  (
1380  cycpp.transform()
1382  )
1383  {
1384  // Force to wanted separation
1385  Info<< "On cyclic translation patch "
1386  << pp.name()
1387  << " forcing uniform separation of "
1388  << cycpp.separationVector() << endl;
1389  const_cast<vectorField&>(cpp.separation()) =
1390  pointField(1, cycpp.separationVector());
1391  }
1392  else
1393  {
1394  const auto& nbr = cycpp.neighbPatch();
1395  const_cast<vectorField&>(cpp.separation()) =
1396  pointField
1397  (
1398  1,
1399  nbr[0].centre(mesh.points())
1400  - cycpp[0].centre(mesh.points())
1401  );
1402  }
1403  }
1404  Info<< "On coupled patch " << pp.name()
1405  << " forcing uniform separation of "
1406  << cpp.separation() << endl;
1407  }
1408  else if (!cpp.parallel())
1409  {
1410  Info<< "On coupled patch " << pp.name()
1411  << " forcing uniform rotation of "
1412  << cpp.forwardT()[0] << endl;
1413 
1414  const_cast<tensorField&>
1415  (
1416  cpp.forwardT()
1417  ).setSize(1);
1418  const_cast<tensorField&>
1419  (
1420  cpp.reverseT()
1421  ).setSize(1);
1422 
1423  Info<< "On coupled patch " << pp.name()
1424  << " forcing uniform rotation of "
1425  << cpp.forwardT() << endl;
1426  }
1427  }
1428  }
1429 
1430  Info<< "Synchronising points." << endl;
1431 
1432  pointField newPoints(mesh.points());
1433 
1434  syncPoints
1435  (
1436  mesh,
1437  newPoints,
1439  point(GREAT, GREAT, GREAT)
1440  );
1441 
1442  scalarField diff(mag(newPoints-mesh.points()));
1443  Info<< "Points changed by average:" << gAverage(diff)
1444  << " max:" << gMax(diff) << nl << endl;
1445 
1446  mesh.movePoints(newPoints);
1447  }
1448 
1449 
1450  // 3. Remove zeros-sized patches
1451  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1452 
1453  Info<< "Removing patches with no faces in them." << nl << endl;
1454  const wordList oldPatchNames(mesh.boundaryMesh().names());
1455  const wordList oldPatchTypes(mesh.boundaryMesh().types());
1457  forAll(oldPatchNames, patchi)
1458  {
1459  const word& pName = oldPatchNames[patchi];
1460  if (mesh.boundaryMesh().findPatchID(pName) == -1)
1461  {
1462  Info<< "Removed zero-sized patch " << pName
1463  << " type " << oldPatchTypes[patchi]
1464  << " at position " << patchi << endl;
1465  }
1466  }
1467 
1468 
1469  if (writeObj)
1470  {
1471  dumpCyclicMatch("final_", mesh);
1472  }
1473  }
1474 
1475  if (!overwrite)
1476  {
1477  ++runTime;
1478  }
1479  else
1480  {
1481  forAll(meshes, meshi)
1482  {
1483  fvMesh& mesh = meshes[meshi];
1484 
1485  mesh.setInstance(oldInstances[meshi]);
1486  }
1487  }
1488 
1489  // More precision (for points data)
1491 
1492  // Write resulting mesh
1493  forAll(meshes, meshi)
1494  {
1495  fvMesh& mesh = meshes[meshi];
1496  Info<< "\n\nWriting repatched mesh " << mesh.name()
1497  << " to " << runTime.timeName() << nl << endl;
1498  mesh.clearOut(); // remove meshPhi
1499  mesh.write();
1502  }
1503 
1504  Info<< "End\n" << endl;
1505 
1506  return 0;
1507 }
1508 
1509 
1510 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream. ...
Definition: dictionary.C:367
label nPoints() const
Number of points supporting patch faces.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
dictionary dict
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:561
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:476
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
rDeltaTY field()
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:498
A class for handling file names.
Definition: fileName.H:72
A list of face labels.
Definition: faceSet.H:47
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.
virtual bool separated() const
Are the planes separated.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Class describing modification of a face.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:600
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Output to file stream as an OSstream, normally using std::ofstream for the actual output...
Definition: OFstream.H:71
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const word dictName("faMeshDefinition")
const DynamicList< label > & region() const
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:884
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1586
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
Definition: fvMesh.C:227
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:388
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicListI.H:308
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:625
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:397
Field reading functions for post-processing utilities.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
wordList types() const
Return a list of patch types.
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.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:134
labelList faceLabels(nFaceLabels)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
virtual const vectorField & separation() const
If the planes are separated the separation vector.
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:47
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
virtual const tensorField & forwardT() const
Return face transformation tensor.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1063
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
points setSize(newPointi)
virtual bool parallel() const
Are the cyclic planes parallel.
virtual bool owner() const
Does this side own the patch ?
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition: topoSet.C:693
label nGlobalPoints() const
Return number of globally shared points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:838
void setSize(const label n)
Alias for resize()
Definition: List.H:325
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:960
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition: IPstream.H:80
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
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...
Definition: labelLists.C:44
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:38
Cyclic plane patch.
T & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
Definition: DynamicListI.H:541
wordList names() const
Return a list of patch names.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
wordList patchNames(nPatches)
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1101
static const word null
An empty word.
Definition: word.H:84
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:399
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.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1296
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1088
String literal.
Definition: keyType.H:82
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition: IOstream.H:440
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1068
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh...
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:192
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type gMax(const FieldField< Field, Type > &f)
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:418
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary. ...
Definition: dictionary.C:521
Output inter-processor communications stream.
Definition: OPstream.H:49
void push_back(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:558
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:410
const vectorField & faceCentres() const
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:49
void transformList(const tensor &rotTensor, UList< T > &field)
Inplace transform a list of elements.
Definition: transformList.C:48
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
vector point
Point is a vector.
Definition: point.H:37
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
const word & name() const
Return reference to name.
Definition: fvMesh.H:387
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:54
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
static void listCombineReduce(UList< T > &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::listReduce with an in-place cop.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:156
Direct mesh changes based on v1.3 polyTopoChange syntax.
const polyBoundaryMesh & patches
DynamicList< T, SizeMin > & shrink()
Calls shrink_to_fit() and returns a reference to the DynamicList.
Definition: DynamicList.H:598
Required Classes.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:185
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:367
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject dictIO
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:28
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
std::enable_if_t< std::is_same_v< bool, TypeT >, bool > set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:498
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:765
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
"buffered" : (MPI_Bsend, MPI_Recv)
const cyclicPolyPatch & neighbPatch() const
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Foam::argList args(argc, argv)
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points list.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Remove zero sized patches. All but processor patches are.
Definition: fvMeshTools.C:363
bool set(const Key &key, const T &obj)
Copy assign a new entry, overwriting existing entries.
Definition: HashTableI.H:174
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
Required Variables.