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-2022 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 "wordRes.H"
56 #include "processorMeshes.H"
57 #include "IOdictionary.H"
58 #include "regionProperties.H"
59 #include "faceAreaWeightAMI2D.H"
60 #include "fvMeshTools.H"
61 #include "ReadFields.H"
62 
63 using namespace Foam;
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
70 }
71 
72 
73 word patchName(const word& name, const fvMesh& mesh0, const fvMesh& mesh1)
74 {
75  word pName(name != "none" ? name : word::null);
76  pName += mesh0.name();
77  pName += "_to_";
78  pName += mesh1.name();
79  return pName;
80 }
81 
82 
83 void matchPatchFaces
84 (
85  const word& entryName,
86  const word& AMIMethod,
87  const dictionary& AMIDict,
88  const PtrList<fvMesh>& meshes,
89 
90  const label meshi,
91  const label nSourcei,
92  const label sourcei,
93  const labelList& patchesi,
94 
95  const label meshj,
96  const label nSourcej,
97  const label sourcej,
98  const labelList& patchesj,
99 
100  DynamicList<labelList>& interfaceMesh0,
101  DynamicList<label>& interfaceSource0,
102  DynamicList<labelList>& interfacePatch0,
103  DynamicList<wordList>& interfaceNames0,
104  DynamicList<List<DynamicList<label>>>& interfaceFaces0,
105 
106  DynamicList<labelList>& interfaceMesh1,
107  DynamicList<label>& interfaceSource1,
108  DynamicList<labelList>& interfacePatch1,
109  DynamicList<wordList>& interfaceNames1,
110  DynamicList<List<DynamicList<label>>>& interfaceFaces1
111 )
112 {
113  // Now we have:
114  // - meshi, sourcei, patchesi
115  // - meshj, sourcej, patchesj
116 
117 
118  // Attempt to match patches
119  forAll(patchesi, i)
120  {
121  const auto& ppi = meshes[meshi].boundaryMesh()[patchesi[i]];
122 
123  forAll(patchesj, j)
124  {
125  const auto& ppj = meshes[meshj].boundaryMesh()[patchesj[j]];
126 
127  // Use AMI to try and find matches
128  auto AMPtr(AMIInterpolation::New(AMIMethod, AMIDict));
129 
130  AMPtr->calculate(ppi, ppj, nullptr);
131  if
132  (
133  gAverage(AMPtr->tgtWeightsSum()) > SMALL
134  || gAverage(AMPtr->srcWeightsSum()) > SMALL
135  )
136  {
137  const label inti = interfaceMesh0.size();
138 
139  Info<< "Introducing interface " << inti << " between"
140  << " mesh " << meshes[meshi].name()
141  //<< " source:" << sourcei
142  << " patch " << ppi.name()
143  << " and mesh " << meshes[meshj].name()
144  //<< " source:" << sourcej
145  << " patch " << ppj.name()
146  << endl;
147 
148  // Mesh 0
149  //~~~~~~~
150 
151  auto& intMesh0 = interfaceMesh0.emplace_back(nSourcei, -1);
152  intMesh0[sourcei] = meshi;
153 
154  interfaceSource0.push_back(sourcei);
155 
156  auto& intPatch0 = interfacePatch0.emplace_back(nSourcei, -1);
157  intPatch0[sourcei] = ppi.index();
158 
159  auto& intNames0 = interfaceNames0.emplace_back(nSourcei);
160  intNames0[sourcei] =
161  patchName(entryName, meshes[meshi], meshes[meshj]);
162 
163 
164  // Mesh 1
165  //~~~~~~~
166 
167  auto& intMesh1 = interfaceMesh1.emplace_back(nSourcej, -1);
168  intMesh1[sourcej] = meshj;
169 
170  interfaceSource1.push_back(sourcej);
171 
172  auto& intPatch1 = interfacePatch1.emplace_back(nSourcej, -1);
173  intPatch1[sourcej] = ppj.index();
174 
175  auto& intNames1 = interfaceNames1.emplace_back(nSourcej);
176  intNames1[sourcej] =
177  patchName(entryName, meshes[meshj], meshes[meshi]);
178 
179  auto& intFaces0 = interfaceFaces0.emplace_back(nSourcei);
180  DynamicList<label>& faces0 = intFaces0[sourcei];
181  faces0.setCapacity(ppi.size());
182 
183  auto& intFaces1 = interfaceFaces1.emplace_back(nSourcej);
184  DynamicList<label>& faces1 = intFaces1[sourcej];
185  faces1.setCapacity(ppj.size());
186 
187 
188  // Mark any mesh faces:
189  // - that have a weight > 0.5.
190  // - contribute to these faces (even if < 0.5)
191 
192  faces0.clear();
193  faces1.clear();
194 
195  // Marked as donor of face with high weight
196  scalarField targetMask;
197  {
198  const scalarField& weights = AMPtr->srcWeightsSum();
199  scalarField mask(weights.size(), Zero);
200  forAll(weights, facei)
201  {
202  const scalar sum = weights[facei];
203  if (sum > 0.5)
204  {
205  mask[facei] = sum;
206  }
207  }
208 
209  // Push source field mask to target
210  targetMask = AMPtr->interpolateToTarget(mask);
211  }
212 
213  scalarField sourceMask;
214  {
215  const scalarField& weights = AMPtr->tgtWeightsSum();
216  scalarField mask(weights.size(), Zero);
217  forAll(weights, facei)
218  {
219  const scalar sum = weights[facei];
220  if (sum > 0.5)
221  {
222  mask[facei] = sum;
223  }
224  }
225 
226  // Push target mask back to source
227  sourceMask = AMPtr->interpolateToSource(mask);
228  }
229 
230  {
231  const scalarField& weights = AMPtr->srcWeightsSum();
232  forAll(weights, facei)
233  {
234  if (weights[facei] > 0.5 || sourceMask[facei] > SMALL)
235  {
236  faces0.push_back(ppi.start()+facei);
237  }
238  }
239  }
240  {
241  const scalarField& weights = AMPtr->tgtWeightsSum();
242  forAll(weights, facei)
243  {
244  if (weights[facei] > 0.5 || targetMask[facei] > SMALL)
245  {
246  faces1.push_back(ppj.start()+facei);
247  }
248  }
249  }
250 
251  faces0.shrink();
252  faces1.shrink();
253  }
254  }
255  }
256 }
257 
258 
259 void matchPatchFaces
260 (
261  const bool includeOwn,
262  const PtrList<fvMesh>& meshes,
263  const List<DynamicList<label>>& interRegionSources,
264  const List<wordList>& patchNames,
265  const labelListListList& matchPatchIDs,
266  const List<wordList>& AMIMethods,
267  List<PtrList<dictionary>> patchInfoDicts,
268 
269  DynamicList<labelList>& interfaceMesh0,
270  DynamicList<label>& interfaceSource0,
271  DynamicList<labelList>& interfacePatch0,
272  DynamicList<List<DynamicList<label>>>& interfaceFaces0,
273  DynamicList<wordList>& interfaceNames0,
274 
275  DynamicList<labelList>& interfaceMesh1,
276  DynamicList<label>& interfaceSource1,
277  DynamicList<labelList>& interfacePatch1,
278  DynamicList<List<DynamicList<label>>>& interfaceFaces1,
279  DynamicList<wordList>& interfaceNames1
280 )
281 {
282  // Add approximate matches
283  forAll(meshes, meshi)
284  {
285  const labelListList& allPatchesi = matchPatchIDs[meshi];
286 
287  const label meshjStart(includeOwn ? meshi : meshi+1);
288 
289  for (label meshj = meshjStart; meshj < meshes.size(); meshj++)
290  {
291  const labelListList& allPatchesj = matchPatchIDs[meshj];
292 
293  for (const label sourcei : interRegionSources[meshi])
294  {
295  const labelList& patchesi = allPatchesi[sourcei];
296 
297  for (const label sourcej : interRegionSources[meshj])
298  {
299  const labelList& patchesj = allPatchesj[sourcej];
300 
301  // Now we have:
302  // - meshi, sourcei, patchesi
303  // - meshj, sourcej, patchesj
304 
305  matchPatchFaces
306  (
307  patchNames[meshi][sourcei],
308  (
309  AMIMethods[meshi][sourcei].size()
310  ? AMIMethods[meshi][sourcei]
311  : faceAreaWeightAMI2D::typeName
312  ),
313  patchInfoDicts[meshi][sourcei],
314  meshes,
315 
316  meshi,
317  allPatchesi.size(), // nSourcei,
318  sourcei,
319  patchesi,
320 
321  meshj,
322  allPatchesj.size(), // nSourcej,
323  sourcej,
324  patchesj,
325 
326  interfaceMesh0,
327  interfaceSource0,
328  interfacePatch0,
329  interfaceNames0,
330  interfaceFaces0,
331 
332  interfaceMesh1,
333  interfaceSource1,
334  interfacePatch1,
335  interfaceNames1,
336  interfaceFaces1
337  );
338  }
339  }
340  }
341  }
342 }
343 
344 
345 void changePatchID
346 (
347  const fvMesh& mesh,
348  const label faceID,
349  const label patchID,
350  polyTopoChange& meshMod
351 )
352 {
353  const label zoneID = mesh.faceZones().whichZone(faceID);
354 
355  bool zoneFlip = false;
356 
357  if (zoneID >= 0)
358  {
359  const faceZone& fZone = mesh.faceZones()[zoneID];
360 
361  zoneFlip = fZone.flipMap()[fZone.whichFace(faceID)];
362  }
363 
364  meshMod.setAction
365  (
367  (
368  mesh.faces()[faceID], // face
369  faceID, // face ID
370  mesh.faceOwner()[faceID], // owner
371  -1, // neighbour
372  false, // flip flux
373  patchID, // patch ID
374  false, // remove from zone
375  zoneID, // zone ID
376  zoneFlip // zone flip
377  )
378  );
379 }
380 
381 
382 void changePatchID
383 (
384  const fvMesh& mesh,
385  const labelList& faceLabels,
386  const label patchID,
387  bitSet& isRepatchedBoundary,
388  polyTopoChange& meshMod
389 )
390 {
391  for (const label facei : faceLabels)
392  {
393  if (mesh.isInternalFace(facei))
394  {
396  << "Face " << facei
397  << " is not an external face of the mesh." << endl
398  << "This application can only repatch"
399  << " existing boundary faces." << exit(FatalError);
400  }
401 
402  if (!isRepatchedBoundary.set(facei-mesh.nInternalFaces()))
403  {
404  static label nWarnings = 0;
405  if (nWarnings == 0)
406  {
407  const label newPatchi = meshMod.region()[facei];
408  //FatalErrorInFunction
410  << "Face " << facei
411  << " at " << mesh.faceCentres()[facei]
412  << " marked for patch " << patchID
413  << " name " << mesh.boundaryMesh()[patchID].name()
414  << " is already marked for patch " << newPatchi
415  << " name " << mesh.boundaryMesh()[newPatchi].name()
416  << ". Suppressing further warnings"
417  //<< exit(FatalError);
418  << endl;
419  }
420  nWarnings++;
421  }
422 
423  changePatchID(mesh, facei, patchID, meshMod);
424  }
425 }
426 
427 
428 // Dump for all patches the current match
429 void dumpCyclicMatch(const fileName& prefix, const polyMesh& mesh)
430 {
431  for (const polyPatch& pp : mesh.boundaryMesh())
432  {
433  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(pp);
434 
435  if (cpp && cpp->owner())
436  {
437  const auto& cycPatch = *cpp;
438  const auto& nbrPatch = cycPatch.neighbPatch();
439 
440  // Dump patches
441  {
442  OFstream str(prefix+cycPatch.name()+".obj");
443  Pout<< "Dumping " << cycPatch.name()
444  << " faces to " << str.name() << endl;
446  (
447  str,
448  cycPatch,
449  cycPatch.points()
450  );
451  }
452 
453  {
454  OFstream str(prefix+nbrPatch.name()+".obj");
455  Pout<< "Dumping " << nbrPatch.name()
456  << " faces to " << str.name() << endl;
458  (
459  str,
460  nbrPatch,
461  nbrPatch.points()
462  );
463  }
464 
465 
466  // Lines between corresponding face centres
467  OFstream str(prefix+cycPatch.name()+nbrPatch.name()+"_match.obj");
468  label vertI = 0;
469 
470  Pout<< "Dumping cyclic match as lines between face centres to "
471  << str.name() << endl;
472 
473  forAll(cycPatch, facei)
474  {
475  const point& fc0 = mesh.faceCentres()[cycPatch.start()+facei];
476  meshTools::writeOBJ(str, fc0);
477  vertI++;
478  const point& fc1 = mesh.faceCentres()[nbrPatch.start()+facei];
479  meshTools::writeOBJ(str, fc1);
480  vertI++;
481 
482  str<< "l " << vertI-1 << ' ' << vertI << nl;
483  }
484  }
485  }
486 }
487 
488 
489 void separateList
490 (
491  const vectorField& separation,
493 )
494 {
495  if (separation.size() == 1)
496  {
497  // Single value for all.
498 
499  forAll(field, i)
500  {
501  field[i] += separation[0];
502  }
503  }
504  else if (separation.size() == field.size())
505  {
506  forAll(field, i)
507  {
508  field[i] += separation[i];
509  }
510  }
511  else
512  {
514  << "Sizes of field and transformation not equal. field:"
515  << field.size() << " transformation:" << separation.size()
516  << abort(FatalError);
517  }
518 }
519 
520 
521 // Synchronise points on both sides of coupled boundaries.
522 template<class CombineOp>
523 void syncPoints
524 (
525  const polyMesh& mesh,
527  const CombineOp& cop,
528  const point& nullValue
529 )
530 {
531  if (points.size() != mesh.nPoints())
532  {
534  << "Number of values " << points.size()
535  << " is not equal to the number of points in the mesh "
536  << mesh.nPoints() << abort(FatalError);
537  }
538 
540 
541  // Is there any coupled patch with transformation?
542  bool hasTransformation = false;
543 
544  if (Pstream::parRun())
545  {
546  const labelList& procPatches = mesh.globalData().processorPatches();
547 
548  // Send
549  for (const label patchi : procPatches)
550  {
551  const polyPatch& pp = patches[patchi];
552  const auto& procPatch = refCast<const processorPolyPatch>(pp);
553 
554  if (pp.nPoints() && procPatch.owner())
555  {
556  // Get data per patchPoint in neighbouring point numbers.
557  pointField patchInfo(procPatch.nPoints(), nullValue);
558 
559  const labelList& meshPts = procPatch.meshPoints();
560  const labelList& nbrPts = procPatch.neighbPoints();
561 
562  forAll(nbrPts, pointi)
563  {
564  label nbrPointi = nbrPts[pointi];
565  if (nbrPointi >= 0 && nbrPointi < patchInfo.size())
566  {
567  patchInfo[nbrPointi] = points[meshPts[pointi]];
568  }
569  }
570 
571  OPstream toNbr
572  (
574  procPatch.neighbProcNo()
575  );
576  toNbr << patchInfo;
577  }
578  }
579 
580 
581  // Receive and set.
582 
583  for (const label patchi : procPatches)
584  {
585  const polyPatch& pp = patches[patchi];
586  const auto& procPatch = refCast<const processorPolyPatch>(pp);
587 
588  if (pp.nPoints() && !procPatch.owner())
589  {
590  pointField nbrPatchInfo(procPatch.nPoints());
591  {
592  // We do not know the number of points on the other side
593  // so cannot use UIPstream::read
594  IPstream fromNbr
595  (
597  procPatch.neighbProcNo()
598  );
599  fromNbr >> nbrPatchInfo;
600  }
601  // Null any value which is not on neighbouring processor
602  nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
603 
604  if (!procPatch.parallel())
605  {
606  hasTransformation = true;
607  transformList(procPatch.forwardT(), nbrPatchInfo);
608  }
609  else if (procPatch.separated())
610  {
611  hasTransformation = true;
612  separateList(-procPatch.separation(), nbrPatchInfo);
613  }
614 
615  const labelList& meshPts = procPatch.meshPoints();
616 
617  forAll(meshPts, pointi)
618  {
619  label meshPointi = meshPts[pointi];
620  points[meshPointi] = nbrPatchInfo[pointi];
621  }
622  }
623  }
624  }
625 
626  // Do the cyclics.
627  for (const polyPatch& pp : patches)
628  {
629  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(pp);
630 
631  if (cpp && cpp->owner())
632  {
633  // Owner does all.
634 
635  const auto& cycPatch = *cpp;
636  const auto& nbrPatch = cycPatch.neighbPatch();
637 
638  const edgeList& coupledPoints = cycPatch.coupledPoints();
639  const labelList& meshPts = cycPatch.meshPoints();
640  const labelList& nbrMeshPts = nbrPatch.meshPoints();
641 
642  pointField half0Values(coupledPoints.size());
643 
644  forAll(coupledPoints, i)
645  {
646  const edge& e = coupledPoints[i];
647  label point0 = meshPts[e[0]];
648  half0Values[i] = points[point0];
649  }
650 
651  if (!cycPatch.parallel())
652  {
653  hasTransformation = true;
654  transformList(cycPatch.reverseT(), half0Values);
655  }
656  else if (cycPatch.separated())
657  {
658  hasTransformation = true;
659  separateList(cycPatch.separation(), half0Values);
660  }
661 
662  forAll(coupledPoints, i)
663  {
664  const edge& e = coupledPoints[i];
665  label point1 = nbrMeshPts[e[1]];
666  points[point1] = half0Values[i];
667  }
668  }
669  }
670 
671  //- Note: hasTransformation is only used for warning messages so
672  // reduction not strictly necessary.
673  //Pstream::reduceOr(hasTransformation);
674 
675  // Synchronize multiple shared points.
676  const globalMeshData& pd = mesh.globalData();
677 
678  if (pd.nGlobalPoints() > 0)
679  {
680  if (hasTransformation)
681  {
683  << "There are decomposed cyclics in this mesh with"
684  << " transformations." << endl
685  << "This is not supported. The result will be incorrect"
686  << endl;
687  }
688 
689 
690  // Values on shared points.
691  pointField sharedPts(pd.nGlobalPoints(), nullValue);
692 
693  forAll(pd.sharedPointLabels(), i)
694  {
695  label meshPointi = pd.sharedPointLabels()[i];
696  // Fill my entries in the shared points
697  sharedPts[pd.sharedPointAddr()[i]] = points[meshPointi];
698  }
699 
700  // Combine - globally consistent
701  Pstream::listCombineReduce(sharedPts, cop);
702 
703  // Now we will all have the same information. Merge it back with
704  // my local information.
705  forAll(pd.sharedPointLabels(), i)
706  {
707  label meshPointi = pd.sharedPointLabels()[i];
708  points[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
709  }
710  }
711 }
712 
713 
714 
715 int main(int argc, char *argv[])
716 {
718  (
719  "Create patches out of selected boundary faces, which are either"
720  " from existing patches or from a faceSet"
721  );
722 
723  #include "addOverwriteOption.H"
724  #include "addAllRegionOptions.H"
725 
726  argList::addOption("dict", "file", "Alternative createPatchDict");
728  (
729  "writeObj",
730  "Write obj files showing the cyclic matching process"
731  );
732 
733  argList::noFunctionObjects(); // Never use function objects
734 
735  #include "setRootCase.H"
736  #include "createTime.H"
737  #include "getAllRegionOptions.H"
738 
739  const bool overwrite = args.found("overwrite");
740 
741  #include "createNamedMeshes.H"
742 
743  const bool writeObj = args.found("writeObj");
744 
745 
746  // Read dictionaries and extract various
747  wordList oldInstances(meshes.size());
750  List<labelListList> matchPatchIDs(meshes.size());
751  List<wordList> matchMethods(meshes.size());
752  List<PtrList<dictionary>> patchInfoDicts(meshes.size());
753 
754  List<DynamicList<label>> interRegionSources(meshes.size());
755  forAll(meshes, meshi)
756  {
757  fvMesh& mesh = meshes[meshi];
759 
760  // If running parallel check same patches everywhere
762 
763  oldInstances[meshi] = mesh.pointsInstance();
764 
765  const word dictName("createPatchDict");
766  #include "setSystemMeshDictionaryIO.H"
767  Info<< "Reading " << dictIO.instance()/dictIO.name() << nl << endl;
768 
769  dicts.set(meshi, new IOdictionary(dictIO));
770 
771  const auto& dict = dicts[meshi];
772  PtrList<dictionary> patchSources(dict.lookup("patches"));
773  patchNames[meshi].setSize(patchSources.size());
774  matchPatchIDs[meshi].setSize(patchSources.size());
775  matchMethods[meshi].setSize(patchSources.size());
776  patchInfoDicts[meshi].setSize(patchSources.size());
777 
778  // Read patch construct info from dictionary
779  forAll(patchSources, sourcei)
780  {
781  const auto& pDict = patchSources[sourcei];
782  patchNames[meshi][sourcei] = pDict.getOrDefault<word>
783  (
784  "name",
785  word::null,
787  );
788 
789  patchInfoDicts[meshi].set
790  (
791  sourcei,
792  new dictionary(pDict.subDict("patchInfo"))
793  );
794  dictionary& patchDict = patchInfoDicts[meshi][sourcei];
795  if (patchDict.found("AMIMethod"))
796  {
797  matchMethods[meshi][sourcei] = patchDict.get<word>("AMIMethod");
798  // Disable full matching since we're trying to use AMIMethod to
799  // find out actual overlap
800  patchDict.add("requireMatch", false);
801  }
802 
803  wordRes matchNames;
804  if (pDict.readIfPresent("patches", matchNames, keyType::LITERAL))
805  {
806  matchPatchIDs[meshi][sourcei] =
807  patches.patchSet(matchNames).sortedToc();
808  }
809 
810  if (pDict.get<word>("constructFrom") == "autoPatch")
811  {
812  interRegionSources[meshi].push_back(sourcei);
813  }
814  }
815  }
816 
817 
818 
819  // Determine matching faces. This is an interface between two
820  // sets of faces:
821  // - originating mesh
822  // - originating patch
823  // - originating (mesh)facelabels
824  // It matches all mesh against each other. Lower numbered mesh gets
825  // postfix 0, higher numbered mesh postfix 1.
826 
827  // Per interface, per patchSource:
828  // 1. the lower numbered mesh
829  DynamicList<labelList> interfaceMesh0;
830  // 1b. the source index (i.e. the patch dictionary)
831  DynamicList<label> interfaceSource0;
832  // 2. the patch on the interfaceMesh0
833  DynamicList<labelList> interfacePatch0;
834  // 3. the facelabels on the interfaceMesh0
835  DynamicList<List<DynamicList<label>>> interfaceFaces0;
836  // 4. generated interface name
837  DynamicList<wordList> interfaceNames0;
838 
839  // Same for the higher numbered mesh
840  DynamicList<labelList> interfaceMesh1;
841  DynamicList<label> interfaceSource1;
842  DynamicList<labelList> interfacePatch1;
843  DynamicList<List<DynamicList<label>>> interfaceFaces1;
844  DynamicList<wordList> interfaceNames1;
845  {
846  // Whether to match to patches in own mesh
847  const bool includeOwn = (meshes.size() == 1);
848 
849  matchPatchFaces
850  (
851  includeOwn,
852  meshes,
853  interRegionSources,
854  patchNames,
855  matchPatchIDs,
856  matchMethods,
857  patchInfoDicts,
858 
859  interfaceMesh0,
860  interfaceSource0,
861  interfacePatch0,
862  interfaceFaces0,
863  interfaceNames0,
864 
865  interfaceMesh1,
866  interfaceSource1,
867  interfacePatch1,
868  interfaceFaces1,
869  interfaceNames1
870  );
871  }
872 
873 
874 
875  // Read fields
886 
887  {
888  forAll(meshes, meshi)
889  {
890  const fvMesh& mesh = meshes[meshi];
891 
892  bool noFields = true;
893  for (const auto& d : patchInfoDicts[meshi])
894  {
895  if (d.found("patchFields"))
896  {
897  noFields = false;
898  }
899  }
900 
901  if (!noFields)
902  {
903  // Read objects in time directory
904  IOobjectList objects(mesh, runTime.timeName());
905 
906  // Read volume fields.
907  ReadFields(mesh, objects, vsFlds[meshi]);
908  ReadFields(mesh, objects, vvFlds[meshi]);
909  ReadFields(mesh, objects, vstFlds[meshi]);
910  ReadFields(mesh, objects, vsymtFlds[meshi]);
911  ReadFields(mesh, objects, vtFlds[meshi]);
912 
913  // Read surface fields.
914  ReadFields(mesh, objects, ssFlds[meshi]);
915  ReadFields(mesh, objects, svFlds[meshi]);
916  ReadFields(mesh, objects, sstFlds[meshi]);
917  ReadFields(mesh, objects, ssymtFlds[meshi]);
918  ReadFields(mesh, objects, stFlds[meshi]);
919  }
920  }
921  }
922 
923 
924 
925  // Loop over all regions
926 
927  forAll(meshes, meshi)
928  {
929  fvMesh& mesh = meshes[meshi];
930 
931  Info<< "\n\nAdding patches to mesh " << mesh.name() << nl << endl;
932 
934  const dictionary& dict = dicts[meshi];
935 
936  if (writeObj)
937  {
938  dumpCyclicMatch("initial_", mesh);
939  }
940 
941  // Read patch construct info from dictionary
942  PtrList<dictionary> patchSources(dict.lookup("patches"));
943 
944 
945  // 1. Add all new patches
946  // ~~~~~~~~~~~~~~~~~~~~~~
947 
948  forAll(patchSources, sourcei)
949  {
950  const dictionary& dict = patchSources[sourcei];
951  const word sourceType(dict.get<word>("constructFrom"));
952  const word patchName
953  (
955  (
956  "name",
957  word::null,
959  )
960  );
961 
962  dictionary patchDict(patchInfoDicts[meshi][sourcei]);
963  patchDict.set("nFaces", 0);
964  patchDict.set("startFace", 0); // Gets overwritten
965 
966 
967  if (sourceType == "autoPatch")
968  {
969  // Special : automatically create the necessary inter-region
970  // coupling patches
971 
972  forAll(interfaceMesh0, inti)
973  {
974  const labelList& allMeshes0 = interfaceMesh0[inti];
975  const wordList& allNames0 = interfaceNames0[inti];
976  const labelList& allMeshes1 = interfaceMesh1[inti];
977  const wordList& allNames1 = interfaceNames1[inti];
978 
979  if
980  (
981  interfaceSource0[inti] == sourcei
982  && allMeshes0[sourcei] == meshi
983  )
984  {
985  // Current mesh is mesh0. mesh1 is the remote mesh.
986  const label sourcej = interfaceSource1[inti];
987  const word& patchName = allNames0[sourcei];
988  if (patches.findPatchID(patchName) == -1)
989  {
990  dictionary allDict(patchDict);
991  const auto& mesh1 = meshes[allMeshes1[sourcej]];
992  allDict.set("sampleRegion", mesh1.name());
993  const auto& destPatch = allNames1[sourcej];
994  allDict.set("samplePatch", destPatch);
995  allDict.set("neighbourPatch", destPatch);
996 
997  Info<< "Adding new patch " << patchName
998  << " from " << allDict << endl;
999 
1000  autoPtr<polyPatch> ppPtr
1001  (
1003  (
1004  patchName,
1005  allDict,
1006  0, // overwritten
1007  patches
1008  )
1009  );
1011  (
1012  mesh,
1013  ppPtr(),
1014  patchDict.subOrEmptyDict("patchFields"),
1016  true
1017  );
1018  }
1019  }
1020  }
1021 
1022  forAll(interfaceMesh1, inti)
1023  {
1024  const labelList& allMeshes0 = interfaceMesh0[inti];
1025  const wordList& allNames0 = interfaceNames0[inti];
1026  const labelList& allMeshes1 = interfaceMesh1[inti];
1027  const wordList& allNames1 = interfaceNames1[inti];
1028 
1029  if
1030  (
1031  interfaceSource1[inti] == sourcei
1032  && allMeshes1[sourcei] == meshi
1033  )
1034  {
1035  // Current mesh is mesh1. mesh0 is the remote mesh.
1036 
1037  const label sourcej = interfaceSource0[inti];
1038  const word& patchName = allNames1[sourcei];
1039  if (patches.findPatchID(patchName) == -1)
1040  {
1041  dictionary allDict(patchDict);
1042  const auto& destPatch = allNames0[sourcej];
1043  const auto& mesh0 = meshes[allMeshes0[sourcej]];
1044  allDict.set("sampleRegion", mesh0.name());
1045  allDict.set("samplePatch", destPatch);
1046  allDict.set("neighbourPatch", destPatch);
1047 
1048  Info<< "Adding new patch " << patchName
1049  << " from " << allDict << endl;
1050 
1051  autoPtr<polyPatch> ppPtr
1052  (
1054  (
1055  patchName,
1056  allDict,
1057  0, // overwritten
1058  patches
1059  )
1060  );
1062  (
1063  mesh,
1064  ppPtr(),
1065  patchDict.subOrEmptyDict("patchFields"),
1067  true
1068  );
1069  }
1070  }
1071  }
1072  }
1073  else
1074  {
1075  if (patches.findPatchID(patchName) == -1)
1076  {
1077  Info<< "Adding new patch " << patchName
1078  << " from " << patchDict << endl;
1079 
1080  autoPtr<polyPatch> ppPtr
1081  (
1083  (
1084  patchName,
1085  patchDict,
1086  0, // overwritten
1087  patches
1088  )
1089  );
1091  (
1092  mesh,
1093  ppPtr(),
1094  patchDict.subOrEmptyDict("patchFields"),
1096  true
1097  );
1098  }
1099  }
1100  }
1101 
1102  Info<< endl;
1103  }
1104 
1105 
1106  // 2. Repatch faces
1107  // ~~~~~~~~~~~~~~~~
1108 
1109  forAll(meshes, meshi)
1110  {
1111  fvMesh& mesh = meshes[meshi];
1112 
1113  Info<< "\n\nRepatching mesh " << mesh.name() << nl << endl;
1114 
1115 
1117  const dictionary& dict = dicts[meshi];
1118 
1119  // Whether to synchronise points
1120  const bool pointSync(dict.get<bool>("pointSync"));
1121 
1122  // Read patch construct info from dictionary
1123  PtrList<dictionary> patchSources(dict.lookup("patches"));
1124 
1125 
1126  polyTopoChange meshMod(mesh);
1127 
1128  // Mark all repatched faces. This makes sure that the faces to repatch
1129  // do not overlap
1130  bitSet isRepatchedBoundary(mesh.nBoundaryFaces());
1131 
1132  forAll(patchSources, sourcei)
1133  {
1134  const dictionary& dict = patchSources[sourcei];
1135  const word patchName
1136  (
1138  (
1139  "name",
1140  word::null,
1142  )
1143  );
1144  const word sourceType(dict.get<word>("constructFrom"));
1145 
1146  if (sourceType == "autoPatch")
1147  {
1148  forAll(interfaceMesh0, inti)
1149  {
1150  const labelList& allMeshes0 = interfaceMesh0[inti];
1151  const wordList& allNames0 = interfaceNames0[inti];
1152  const auto& allFaces0 = interfaceFaces0[inti];
1153 
1154  const labelList& allMeshes1 = interfaceMesh1[inti];
1155  const wordList& allNames1 = interfaceNames1[inti];
1156  const auto& allFaces1 = interfaceFaces1[inti];
1157 
1158  if
1159  (
1160  interfaceSource0[inti] == sourcei
1161  && allMeshes0[sourcei] == meshi
1162  )
1163  {
1164  // Current mesh is mesh0. mesh1 is the remote mesh.
1165 
1166  const label destPatchi =
1167  patches.findPatchID(allNames0[sourcei], false);
1168 
1169  //const auto& mesh1 =
1170  // meshes[allMeshes1[interfaceSource1[inti]]];
1171  //Pout<< "Matched mesh:" << mesh.name()
1172  // << " to mesh:" << mesh1.name()
1173  // << " through:" << allNames0[sourcei] << endl;
1174 
1175  changePatchID
1176  (
1177  mesh,
1178  allFaces0[sourcei],
1179  destPatchi,
1180  isRepatchedBoundary,
1181  meshMod
1182  );
1183  }
1184  if
1185  (
1186  interfaceSource1[inti] == sourcei
1187  && allMeshes1[sourcei] == meshi
1188  )
1189  {
1190  // Current mesh is mesh1. mesh0 is the remote mesh.
1191 
1192  const label destPatchi =
1193  patches.findPatchID(allNames1[sourcei], false);
1194 
1195  //const auto& mesh0 =
1196  // meshes[allMeshes0[interfaceSource0[inti]]];
1197  //Pout<< "Matched mesh:" << mesh.name()
1198  // << " to mesh:" << mesh0.name()
1199  // << " through:" << allNames1[sourcei] << endl;
1200 
1201  changePatchID
1202  (
1203  mesh,
1204  allFaces1[sourcei],
1205  destPatchi,
1206  isRepatchedBoundary,
1207  meshMod
1208  );
1209  }
1210  }
1211  }
1212  else if (sourceType == "patches")
1213  {
1214  const label destPatchi = patches.findPatchID(patchName, false);
1215  labelHashSet patchSources
1216  (
1217  patches.patchSet(dict.get<wordRes>("patches"))
1218  );
1219 
1220  // Repatch faces of the patches.
1221  for (const label patchi : patchSources)
1222  {
1223  const polyPatch& pp = patches[patchi];
1224 
1225  Info<< "Moving faces from patch " << pp.name()
1226  << " to patch " << destPatchi << endl;
1227 
1228  changePatchID
1229  (
1230  mesh,
1231  pp.start()+identity(pp.size()),
1232  destPatchi,
1233  isRepatchedBoundary,
1234  meshMod
1235  );
1236  }
1237  }
1238  else if (sourceType == "set")
1239  {
1240  const label destPatchi = patches.findPatchID(patchName, false);
1241  const word setName(dict.get<word>("set"));
1242 
1243  faceSet set(mesh, setName);
1244 
1245  Info<< "Read " << returnReduce(set.size(), sumOp<label>())
1246  << " faces from faceSet " << set.name() << endl;
1247 
1248  // Sort (since faceSet contains faces in arbitrary order)
1250 
1251  changePatchID
1252  (
1253  mesh,
1254  faceLabels,
1255  destPatchi,
1256  isRepatchedBoundary,
1257  meshMod
1258  );
1259  }
1260  else
1261  {
1263  << "Invalid source type " << sourceType << endl
1264  << "Valid source types are 'patches' 'set'"
1265  << exit(FatalError);
1266  }
1267  }
1268 
1269 
1270  // Change mesh, use inflation to reforce calculation of transformation
1271  // tensors.
1272  Info<< "Doing topology modification to order faces." << nl << endl;
1273  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, true);
1274 
1275  if (map().hasMotionPoints())
1276  {
1277  mesh.movePoints(map().preMotionPoints());
1278  }
1279  else
1280  {
1281  // Force calculation of transformation tensors
1283  }
1284 
1285 
1286  // Update fields
1287  mesh.updateMesh(map());
1288 
1289 
1290  // Update numbering pointing to meshi
1291  const auto& oldToNew = map().reverseFaceMap();
1292  forAll(interfaceMesh0, inti)
1293  {
1294  const labelList& allMeshes0 = interfaceMesh0[inti];
1295  forAll(allMeshes0, sourcei)
1296  {
1297  if (allMeshes0[sourcei] == meshi)
1298  {
1299  inplaceRenumber(oldToNew, interfaceFaces0[inti][sourcei]);
1300  }
1301  }
1302  }
1303  forAll(interfaceMesh1, inti)
1304  {
1305  const labelList& allMeshes1 = interfaceMesh1[inti];
1306  forAll(allMeshes1, sourcei)
1307  {
1308  if (allMeshes1[sourcei] == meshi)
1309  {
1310  inplaceRenumber(oldToNew, interfaceFaces1[inti][sourcei]);
1311  }
1312  }
1313  }
1314 
1315 
1316  if (writeObj)
1317  {
1318  dumpCyclicMatch("coupled_", mesh);
1319  }
1320 
1321  // Synchronise points.
1322  if (!pointSync)
1323  {
1324  Info<< "Not synchronising points." << nl << endl;
1325  }
1326  else
1327  {
1328  Info<< "Synchronising points." << nl << endl;
1329 
1330  // This is a bit tricky. Both normal and position might be out and
1331  // current separation also includes the normal
1332  // ( separation_ = (nf&(Cr - Cf))*nf ).
1333 
1334  // For cyclic patches:
1335  // - for separated ones use user specified offset vector
1336 
1337  forAll(mesh.boundaryMesh(), patchi)
1338  {
1339  const polyPatch& pp = mesh.boundaryMesh()[patchi];
1340 
1341  if (pp.size() && isA<coupledPolyPatch>(pp))
1342  {
1343  const coupledPolyPatch& cpp =
1344  refCast<const coupledPolyPatch>(pp);
1345 
1346  if (cpp.separated())
1347  {
1348  Info<< "On coupled patch " << pp.name()
1349  << " separation[0] was "
1350  << cpp.separation()[0] << endl;
1351 
1352  if (isA<cyclicPolyPatch>(pp) && pp.size())
1353  {
1354  const auto& cycpp =
1355  refCast<const cyclicPolyPatch>(pp);
1356 
1357  if
1358  (
1359  cycpp.transform()
1361  )
1362  {
1363  // Force to wanted separation
1364  Info<< "On cyclic translation patch "
1365  << pp.name()
1366  << " forcing uniform separation of "
1367  << cycpp.separationVector() << endl;
1368  const_cast<vectorField&>(cpp.separation()) =
1369  pointField(1, cycpp.separationVector());
1370  }
1371  else
1372  {
1373  const auto& nbr = cycpp.neighbPatch();
1374  const_cast<vectorField&>(cpp.separation()) =
1375  pointField
1376  (
1377  1,
1378  nbr[0].centre(mesh.points())
1379  - cycpp[0].centre(mesh.points())
1380  );
1381  }
1382  }
1383  Info<< "On coupled patch " << pp.name()
1384  << " forcing uniform separation of "
1385  << cpp.separation() << endl;
1386  }
1387  else if (!cpp.parallel())
1388  {
1389  Info<< "On coupled patch " << pp.name()
1390  << " forcing uniform rotation of "
1391  << cpp.forwardT()[0] << endl;
1392 
1393  const_cast<tensorField&>
1394  (
1395  cpp.forwardT()
1396  ).setSize(1);
1397  const_cast<tensorField&>
1398  (
1399  cpp.reverseT()
1400  ).setSize(1);
1401 
1402  Info<< "On coupled patch " << pp.name()
1403  << " forcing uniform rotation of "
1404  << cpp.forwardT() << endl;
1405  }
1406  }
1407  }
1408 
1409  Info<< "Synchronising points." << endl;
1410 
1411  pointField newPoints(mesh.points());
1412 
1413  syncPoints
1414  (
1415  mesh,
1416  newPoints,
1418  point(GREAT, GREAT, GREAT)
1419  );
1420 
1421  scalarField diff(mag(newPoints-mesh.points()));
1422  Info<< "Points changed by average:" << gAverage(diff)
1423  << " max:" << gMax(diff) << nl << endl;
1424 
1425  mesh.movePoints(newPoints);
1426  }
1427 
1428 
1429  // 3. Remove zeros-sized patches
1430  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1431 
1432  Info<< "Removing patches with no faces in them." << nl << endl;
1433  const wordList oldPatchNames(mesh.boundaryMesh().names());
1434  const wordList oldPatchTypes(mesh.boundaryMesh().types());
1436  forAll(oldPatchNames, patchi)
1437  {
1438  const word& pName = oldPatchNames[patchi];
1439  if (mesh.boundaryMesh().findPatchID(pName) == -1)
1440  {
1441  Info<< "Removed zero-sized patch " << pName
1442  << " type " << oldPatchTypes[patchi]
1443  << " at position " << patchi << endl;
1444  }
1445  }
1446 
1447 
1448  if (writeObj)
1449  {
1450  dumpCyclicMatch("final_", mesh);
1451  }
1452  }
1453 
1454  if (!overwrite)
1455  {
1456  ++runTime;
1457  }
1458  else
1459  {
1460  forAll(meshes, meshi)
1461  {
1462  fvMesh& mesh = meshes[meshi];
1463 
1464  mesh.setInstance(oldInstances[meshi]);
1465  }
1466  }
1467 
1468  // Set the precision of the points data to 10
1470 
1471  // Write resulting mesh
1472  forAll(meshes, meshi)
1473  {
1474  fvMesh& mesh = meshes[meshi];
1475  Info<< "\n\nWriting repatched mesh " << mesh.name()
1476  << " to " << runTime.timeName() << nl << endl;
1477  mesh.clearOut(); // remove meshPhi
1478  mesh.write();
1481  }
1482 
1483  Info<< "End\n" << endl;
1484 
1485  return 0;
1486 }
1487 
1488 
1489 // ************************************************************************* //
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
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:227
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:547
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
"blocking" : (MPI_Bsend, MPI_Recv)
rDeltaTY field()
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:504
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:598
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.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Output to file stream, using an OSstream.
Definition: OFstream.H:49
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:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:918
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:423
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
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:303
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:128
labelList faceLabels(nFaceLabels)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
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:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Input inter-processor communications stream.
Definition: IPstream.H:49
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:618
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:853
void setSize(const label n)
Alias for resize()
Definition: List.H:316
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:994
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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:608
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:538
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:1116
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:385
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
String literal.
Definition: keyType.H:82
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:1102
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
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 gMax(const FieldField< Field, Type > &f)
DynamicList< T, SizeMin > & shrink()
Calls shrink_to_fit() and returns a reference to the DynamicList.
Definition: DynamicListI.H:447
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:405
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:555
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:346
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 (not the indices) of 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
Type gAverage(const FieldField< Field, Type > &f)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
Direct mesh changes based on v1.3 polyTopoChange syntax.
const polyBoundaryMesh & patches
Required Classes.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
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
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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...
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.
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.
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
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.