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