createPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  createPatch
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Create patches out of selected boundary faces, which are either
35  from existing patches or from a faceSet.
36 
37  More specifically it:
38  - creates new patches (from selected boundary faces).
39  Synchronise faces on coupled patches.
40  - synchronises points on coupled boundaries
41  - remove patches with 0 faces in them
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #include "cyclicPolyPatch.H"
46 #include "syncTools.H"
47 #include "argList.H"
48 #include "Time.H"
49 #include "OFstream.H"
50 #include "meshTools.H"
51 #include "faceSet.H"
52 #include "IOPtrList.H"
53 #include "polyTopoChange.H"
54 #include "polyModifyFace.H"
55 #include "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 (UPstream::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  // We do not know the number of points on the other side
591  // so cannot use UIPstream::read
592 
593  pointField nbrPatchInfo;
594  IPstream::recv(nbrPatchInfo, procPatch.neighbProcNo());
595 
596  // Null any value which is not on neighbouring processor
597  nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
598 
599  if (!procPatch.parallel())
600  {
601  hasTransformation = true;
602  transformList(procPatch.forwardT(), nbrPatchInfo);
603  }
604  else if (procPatch.separated())
605  {
606  hasTransformation = true;
607  separateList(-procPatch.separation(), nbrPatchInfo);
608  }
609 
610  const labelList& meshPts = procPatch.meshPoints();
611 
612  forAll(meshPts, pointi)
613  {
614  label meshPointi = meshPts[pointi];
615  points[meshPointi] = nbrPatchInfo[pointi];
616  }
617  }
618  }
619  }
620 
621  // Do the cyclics.
622  for (const polyPatch& pp : patches)
623  {
624  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(pp);
625 
626  if (cpp && cpp->owner())
627  {
628  // Owner does all.
629 
630  const auto& cycPatch = *cpp;
631  const auto& nbrPatch = cycPatch.neighbPatch();
632 
633  const edgeList& coupledPoints = cycPatch.coupledPoints();
634  const labelList& meshPts = cycPatch.meshPoints();
635  const labelList& nbrMeshPts = nbrPatch.meshPoints();
636 
637  pointField half0Values(coupledPoints.size());
638 
639  forAll(coupledPoints, i)
640  {
641  const edge& e = coupledPoints[i];
642  label point0 = meshPts[e[0]];
643  half0Values[i] = points[point0];
644  }
645 
646  if (!cycPatch.parallel())
647  {
648  hasTransformation = true;
649  transformList(cycPatch.reverseT(), half0Values);
650  }
651  else if (cycPatch.separated())
652  {
653  hasTransformation = true;
654  separateList(cycPatch.separation(), half0Values);
655  }
656 
657  forAll(coupledPoints, i)
658  {
659  const edge& e = coupledPoints[i];
660  label point1 = nbrMeshPts[e[1]];
661  points[point1] = half0Values[i];
662  }
663  }
664  }
665 
666  //- Note: hasTransformation is only used for warning messages so
667  // reduction not strictly necessary.
668  //Pstream::reduceOr(hasTransformation);
669 
670  // Synchronize multiple shared points.
671  const globalMeshData& pd = mesh.globalData();
672 
673  if (pd.nGlobalPoints() > 0)
674  {
675  if (hasTransformation)
676  {
678  << "There are decomposed cyclics in this mesh with"
679  << " transformations." << endl
680  << "This is not supported. The result will be incorrect"
681  << endl;
682  }
683 
684 
685  // Values on shared points.
686  pointField sharedPts(pd.nGlobalPoints(), nullValue);
687 
688  forAll(pd.sharedPointLabels(), i)
689  {
690  label meshPointi = pd.sharedPointLabels()[i];
691  // Fill my entries in the shared points
692  sharedPts[pd.sharedPointAddr()[i]] = points[meshPointi];
693  }
694 
695  // Combine - globally consistent
696  Pstream::listCombineReduce(sharedPts, cop);
697 
698  // Now we will all have the same information. Merge it back with
699  // my local information.
700  forAll(pd.sharedPointLabels(), i)
701  {
702  label meshPointi = pd.sharedPointLabels()[i];
703  points[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
704  }
705  }
706 }
707 
708 
709 
710 int main(int argc, char *argv[])
711 {
713  (
714  "Create patches out of selected boundary faces, which are either"
715  " from existing patches or from a faceSet"
716  );
717 
718  #include "addOverwriteOption.H"
719  #include "addAllRegionOptions.H"
720 
721  argList::addOption("dict", "file", "Alternative createPatchDict");
723  (
724  "writeObj",
725  "Write obj files showing the cyclic matching process"
726  );
727 
728  argList::noFunctionObjects(); // Never use function objects
729 
730  #include "setRootCase.H"
731  #include "createTime.H"
732  #include "getAllRegionOptions.H"
733 
734  const bool overwrite = args.found("overwrite");
735 
736  #include "createNamedMeshes.H"
737 
738  const bool writeObj = args.found("writeObj");
739 
740 
741  // Read dictionaries and extract various
742  wordList oldInstances(meshes.size());
745  List<labelListList> matchPatchIDs(meshes.size());
746  List<wordList> matchMethods(meshes.size());
747  List<PtrList<dictionary>> patchInfoDicts(meshes.size());
748 
749  List<DynamicList<label>> interRegionSources(meshes.size());
750  forAll(meshes, meshi)
751  {
752  fvMesh& mesh = meshes[meshi];
754 
755  // If running parallel check same patches everywhere
757 
758  oldInstances[meshi] = mesh.pointsInstance();
759 
760  const word dictName("createPatchDict");
761  #include "setSystemMeshDictionaryIO.H"
762  Info<< "Reading " << dictIO.instance()/dictIO.name() << nl << endl;
763 
764  dicts.set(meshi, new IOdictionary(dictIO));
765 
766  const auto& dict = dicts[meshi];
767  PtrList<dictionary> patchSources(dict.lookup("patches"));
768  patchNames[meshi].setSize(patchSources.size());
769  matchPatchIDs[meshi].setSize(patchSources.size());
770  matchMethods[meshi].setSize(patchSources.size());
771  patchInfoDicts[meshi].setSize(patchSources.size());
772 
773  // Read patch construct info from dictionary
774  forAll(patchSources, sourcei)
775  {
776  const auto& pDict = patchSources[sourcei];
777  patchNames[meshi][sourcei] = pDict.getOrDefault<word>
778  (
779  "name",
780  word::null,
782  );
783 
784  patchInfoDicts[meshi].set
785  (
786  sourcei,
787  new dictionary(pDict.subDict("patchInfo"))
788  );
789  dictionary& patchDict = patchInfoDicts[meshi][sourcei];
790  if (patchDict.found("AMIMethod"))
791  {
792  matchMethods[meshi][sourcei] = patchDict.get<word>("AMIMethod");
793  // Disable full matching since we're trying to use AMIMethod to
794  // find out actual overlap
795  patchDict.add("requireMatch", false);
796  }
797 
798  wordRes matchNames;
799  if (pDict.readIfPresent("patches", matchNames, keyType::LITERAL))
800  {
801  matchPatchIDs[meshi][sourcei] =
802  patches.patchSet(matchNames).sortedToc();
803  }
804 
805  if (pDict.get<word>("constructFrom") == "autoPatch")
806  {
807  interRegionSources[meshi].push_back(sourcei);
808  }
809  }
810  }
811 
812 
813 
814  // Determine matching faces. This is an interface between two
815  // sets of faces:
816  // - originating mesh
817  // - originating patch
818  // - originating (mesh)facelabels
819  // It matches all mesh against each other. Lower numbered mesh gets
820  // postfix 0, higher numbered mesh postfix 1.
821 
822  // Per interface, per patchSource:
823  // 1. the lower numbered mesh
824  DynamicList<labelList> interfaceMesh0;
825  // 1b. the source index (i.e. the patch dictionary)
826  DynamicList<label> interfaceSource0;
827  // 2. the patch on the interfaceMesh0
828  DynamicList<labelList> interfacePatch0;
829  // 3. the facelabels on the interfaceMesh0
830  DynamicList<List<DynamicList<label>>> interfaceFaces0;
831  // 4. generated interface name
832  DynamicList<wordList> interfaceNames0;
833 
834  // Same for the higher numbered mesh
835  DynamicList<labelList> interfaceMesh1;
836  DynamicList<label> interfaceSource1;
837  DynamicList<labelList> interfacePatch1;
838  DynamicList<List<DynamicList<label>>> interfaceFaces1;
839  DynamicList<wordList> interfaceNames1;
840  {
841  // Whether to match to patches in own mesh
842  const bool includeOwn = (meshes.size() == 1);
843 
844  matchPatchFaces
845  (
846  includeOwn,
847  meshes,
848  interRegionSources,
849  patchNames,
850  matchPatchIDs,
851  matchMethods,
852  patchInfoDicts,
853 
854  interfaceMesh0,
855  interfaceSource0,
856  interfacePatch0,
857  interfaceFaces0,
858  interfaceNames0,
859 
860  interfaceMesh1,
861  interfaceSource1,
862  interfacePatch1,
863  interfaceFaces1,
864  interfaceNames1
865  );
866  }
867 
868 
869 
870  // Read fields
881 
882  {
883  forAll(meshes, meshi)
884  {
885  const fvMesh& mesh = meshes[meshi];
886 
887  bool noFields = true;
888  for (const auto& d : patchInfoDicts[meshi])
889  {
890  if (d.found("patchFields"))
891  {
892  noFields = false;
893  }
894  }
895 
896  if (!noFields)
897  {
898  // Read objects in time directory
899  IOobjectList objects(mesh, runTime.timeName());
900 
901  // Read volume fields.
902  ReadFields(mesh, objects, vsFlds[meshi]);
903  ReadFields(mesh, objects, vvFlds[meshi]);
904  ReadFields(mesh, objects, vstFlds[meshi]);
905  ReadFields(mesh, objects, vsymtFlds[meshi]);
906  ReadFields(mesh, objects, vtFlds[meshi]);
907 
908  // Read surface fields.
909  ReadFields(mesh, objects, ssFlds[meshi]);
910  ReadFields(mesh, objects, svFlds[meshi]);
911  ReadFields(mesh, objects, sstFlds[meshi]);
912  ReadFields(mesh, objects, ssymtFlds[meshi]);
913  ReadFields(mesh, objects, stFlds[meshi]);
914  }
915  }
916  }
917 
918 
919 
920  // Loop over all regions
921 
922  forAll(meshes, meshi)
923  {
924  fvMesh& mesh = meshes[meshi];
925 
926  Info<< "\n\nAdding patches to mesh " << mesh.name() << nl << endl;
927 
929  const dictionary& dict = dicts[meshi];
930 
931  if (writeObj)
932  {
933  dumpCyclicMatch("initial_", mesh);
934  }
935 
936  // Read patch construct info from dictionary
937  PtrList<dictionary> patchSources(dict.lookup("patches"));
938 
939 
940  // 1. Add all new patches
941  // ~~~~~~~~~~~~~~~~~~~~~~
942 
943  forAll(patchSources, sourcei)
944  {
945  const dictionary& dict = patchSources[sourcei];
946  const word sourceType(dict.get<word>("constructFrom"));
947  const word patchName
948  (
950  (
951  "name",
952  word::null,
954  )
955  );
956 
957  dictionary patchDict(patchInfoDicts[meshi][sourcei]);
958  patchDict.set("nFaces", 0);
959  patchDict.set("startFace", 0); // Gets overwritten
960 
961 
962  if (sourceType == "autoPatch")
963  {
964  // Special : automatically create the necessary inter-region
965  // coupling patches
966 
967  forAll(interfaceMesh0, inti)
968  {
969  const labelList& allMeshes0 = interfaceMesh0[inti];
970  const wordList& allNames0 = interfaceNames0[inti];
971  const labelList& allMeshes1 = interfaceMesh1[inti];
972  const wordList& allNames1 = interfaceNames1[inti];
973 
974  if
975  (
976  interfaceSource0[inti] == sourcei
977  && allMeshes0[sourcei] == meshi
978  )
979  {
980  // Current mesh is mesh0. mesh1 is the remote mesh.
981  const label sourcej = interfaceSource1[inti];
982  const word& patchName = allNames0[sourcei];
983  if (patches.findPatchID(patchName) == -1)
984  {
985  dictionary allDict(patchDict);
986  const auto& mesh1 = meshes[allMeshes1[sourcej]];
987  allDict.set("sampleRegion", mesh1.name());
988  const auto& destPatch = allNames1[sourcej];
989  allDict.set("samplePatch", destPatch);
990  allDict.set("neighbourPatch", destPatch);
991 
992  Info<< "Adding new patch " << patchName
993  << " from " << allDict << endl;
994 
995  autoPtr<polyPatch> ppPtr
996  (
998  (
999  patchName,
1000  allDict,
1001  0, // overwritten
1002  patches
1003  )
1004  );
1006  (
1007  mesh,
1008  ppPtr(),
1009  patchDict.subOrEmptyDict("patchFields"),
1011  true
1012  );
1013  }
1014  }
1015  }
1016 
1017  forAll(interfaceMesh1, inti)
1018  {
1019  const labelList& allMeshes0 = interfaceMesh0[inti];
1020  const wordList& allNames0 = interfaceNames0[inti];
1021  const labelList& allMeshes1 = interfaceMesh1[inti];
1022  const wordList& allNames1 = interfaceNames1[inti];
1023 
1024  if
1025  (
1026  interfaceSource1[inti] == sourcei
1027  && allMeshes1[sourcei] == meshi
1028  )
1029  {
1030  // Current mesh is mesh1. mesh0 is the remote mesh.
1031 
1032  const label sourcej = interfaceSource0[inti];
1033  const word& patchName = allNames1[sourcei];
1034  if (patches.findPatchID(patchName) == -1)
1035  {
1036  dictionary allDict(patchDict);
1037  const auto& destPatch = allNames0[sourcej];
1038  const auto& mesh0 = meshes[allMeshes0[sourcej]];
1039  allDict.set("sampleRegion", mesh0.name());
1040  allDict.set("samplePatch", destPatch);
1041  allDict.set("neighbourPatch", destPatch);
1042 
1043  Info<< "Adding new patch " << patchName
1044  << " from " << allDict << endl;
1045 
1046  autoPtr<polyPatch> ppPtr
1047  (
1049  (
1050  patchName,
1051  allDict,
1052  0, // overwritten
1053  patches
1054  )
1055  );
1057  (
1058  mesh,
1059  ppPtr(),
1060  patchDict.subOrEmptyDict("patchFields"),
1062  true
1063  );
1064  }
1065  }
1066  }
1067  }
1068  else
1069  {
1070  if (patches.findPatchID(patchName) == -1)
1071  {
1072  Info<< "Adding new patch " << patchName
1073  << " from " << patchDict << endl;
1074 
1075  autoPtr<polyPatch> ppPtr
1076  (
1078  (
1079  patchName,
1080  patchDict,
1081  0, // overwritten
1082  patches
1083  )
1084  );
1086  (
1087  mesh,
1088  ppPtr(),
1089  patchDict.subOrEmptyDict("patchFields"),
1091  true
1092  );
1093  }
1094  }
1095  }
1096 
1097  Info<< endl;
1098  }
1099 
1100 
1101  // 2. Repatch faces
1102  // ~~~~~~~~~~~~~~~~
1103 
1104  forAll(meshes, meshi)
1105  {
1106  fvMesh& mesh = meshes[meshi];
1107 
1108  Info<< "\n\nRepatching mesh " << mesh.name() << nl << endl;
1109 
1110 
1112  const dictionary& dict = dicts[meshi];
1113 
1114  // Whether to synchronise points
1115  const bool pointSync(dict.get<bool>("pointSync"));
1116 
1117  // Read patch construct info from dictionary
1118  PtrList<dictionary> patchSources(dict.lookup("patches"));
1119 
1120 
1121  polyTopoChange meshMod(mesh);
1122 
1123  // Mark all repatched faces. This makes sure that the faces to repatch
1124  // do not overlap
1125  bitSet isRepatchedBoundary(mesh.nBoundaryFaces());
1126 
1127  forAll(patchSources, sourcei)
1128  {
1129  const dictionary& dict = patchSources[sourcei];
1130  const word patchName
1131  (
1133  (
1134  "name",
1135  word::null,
1137  )
1138  );
1139  const word sourceType(dict.get<word>("constructFrom"));
1140 
1141  if (sourceType == "autoPatch")
1142  {
1143  forAll(interfaceMesh0, inti)
1144  {
1145  const labelList& allMeshes0 = interfaceMesh0[inti];
1146  const wordList& allNames0 = interfaceNames0[inti];
1147  const auto& allFaces0 = interfaceFaces0[inti];
1148 
1149  const labelList& allMeshes1 = interfaceMesh1[inti];
1150  const wordList& allNames1 = interfaceNames1[inti];
1151  const auto& allFaces1 = interfaceFaces1[inti];
1152 
1153  if
1154  (
1155  interfaceSource0[inti] == sourcei
1156  && allMeshes0[sourcei] == meshi
1157  )
1158  {
1159  // Current mesh is mesh0. mesh1 is the remote mesh.
1160 
1161  const label destPatchi =
1162  patches.findPatchID(allNames0[sourcei], false);
1163 
1164  //const auto& mesh1 =
1165  // meshes[allMeshes1[interfaceSource1[inti]]];
1166  //Pout<< "Matched mesh:" << mesh.name()
1167  // << " to mesh:" << mesh1.name()
1168  // << " through:" << allNames0[sourcei] << endl;
1169 
1170  changePatchID
1171  (
1172  mesh,
1173  allFaces0[sourcei],
1174  destPatchi,
1175  isRepatchedBoundary,
1176  meshMod
1177  );
1178  }
1179  if
1180  (
1181  interfaceSource1[inti] == sourcei
1182  && allMeshes1[sourcei] == meshi
1183  )
1184  {
1185  // Current mesh is mesh1. mesh0 is the remote mesh.
1186 
1187  const label destPatchi =
1188  patches.findPatchID(allNames1[sourcei], false);
1189 
1190  //const auto& mesh0 =
1191  // meshes[allMeshes0[interfaceSource0[inti]]];
1192  //Pout<< "Matched mesh:" << mesh.name()
1193  // << " to mesh:" << mesh0.name()
1194  // << " through:" << allNames1[sourcei] << endl;
1195 
1196  changePatchID
1197  (
1198  mesh,
1199  allFaces1[sourcei],
1200  destPatchi,
1201  isRepatchedBoundary,
1202  meshMod
1203  );
1204  }
1205  }
1206  }
1207  else if (sourceType == "patches")
1208  {
1209  const label destPatchi = patches.findPatchID(patchName, false);
1210  labelHashSet patchSources
1211  (
1212  patches.patchSet(dict.get<wordRes>("patches"))
1213  );
1214 
1215  // Repatch faces of the patches.
1216  for (const label patchi : patchSources)
1217  {
1218  const polyPatch& pp = patches[patchi];
1219 
1220  Info<< "Moving faces from patch " << pp.name()
1221  << " to patch " << destPatchi << endl;
1222 
1223  changePatchID
1224  (
1225  mesh,
1226  pp.start()+identity(pp.size()),
1227  destPatchi,
1228  isRepatchedBoundary,
1229  meshMod
1230  );
1231  }
1232  }
1233  else if (sourceType == "set")
1234  {
1235  const label destPatchi = patches.findPatchID(patchName, false);
1236  const word setName(dict.get<word>("set"));
1237 
1238  faceSet set(mesh, setName);
1239 
1240  Info<< "Read " << returnReduce(set.size(), sumOp<label>())
1241  << " faces from faceSet " << set.name() << endl;
1242 
1243  // Sort (since faceSet contains faces in arbitrary order)
1245 
1246  changePatchID
1247  (
1248  mesh,
1249  faceLabels,
1250  destPatchi,
1251  isRepatchedBoundary,
1252  meshMod
1253  );
1254  }
1255  else
1256  {
1258  << "Invalid source type " << sourceType << endl
1259  << "Valid source types are 'patches' 'set'"
1260  << exit(FatalError);
1261  }
1262  }
1263 
1264 
1265  // Change mesh, use inflation to reforce calculation of transformation
1266  // tensors.
1267  Info<< "Doing topology modification to order faces." << nl << endl;
1268  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, true);
1269 
1270  if (map().hasMotionPoints())
1271  {
1272  mesh.movePoints(map().preMotionPoints());
1273  }
1274  else
1275  {
1276  // Force calculation of transformation tensors
1278  }
1279 
1280 
1281  // Update fields
1282  mesh.updateMesh(map());
1283 
1284 
1285  // Update numbering pointing to meshi
1286  const auto& oldToNew = map().reverseFaceMap();
1287  forAll(interfaceMesh0, inti)
1288  {
1289  const labelList& allMeshes0 = interfaceMesh0[inti];
1290  forAll(allMeshes0, sourcei)
1291  {
1292  if (allMeshes0[sourcei] == meshi)
1293  {
1294  inplaceRenumber(oldToNew, interfaceFaces0[inti][sourcei]);
1295  }
1296  }
1297  }
1298  forAll(interfaceMesh1, inti)
1299  {
1300  const labelList& allMeshes1 = interfaceMesh1[inti];
1301  forAll(allMeshes1, sourcei)
1302  {
1303  if (allMeshes1[sourcei] == meshi)
1304  {
1305  inplaceRenumber(oldToNew, interfaceFaces1[inti][sourcei]);
1306  }
1307  }
1308  }
1309 
1310 
1311  if (writeObj)
1312  {
1313  dumpCyclicMatch("coupled_", mesh);
1314  }
1315 
1316  // Synchronise points.
1317  if (!pointSync)
1318  {
1319  Info<< "Not synchronising points." << nl << endl;
1320  }
1321  else
1322  {
1323  Info<< "Synchronising points." << nl << endl;
1324 
1325  // This is a bit tricky. Both normal and position might be out and
1326  // current separation also includes the normal
1327  // ( separation_ = (nf&(Cr - Cf))*nf ).
1328 
1329  // For cyclic patches:
1330  // - for separated ones use user specified offset vector
1331 
1332  forAll(mesh.boundaryMesh(), patchi)
1333  {
1334  const polyPatch& pp = mesh.boundaryMesh()[patchi];
1335 
1336  if (pp.size() && isA<coupledPolyPatch>(pp))
1337  {
1338  const coupledPolyPatch& cpp =
1339  refCast<const coupledPolyPatch>(pp);
1340 
1341  if (cpp.separated())
1342  {
1343  Info<< "On coupled patch " << pp.name()
1344  << " separation[0] was "
1345  << cpp.separation()[0] << endl;
1346 
1347  if (isA<cyclicPolyPatch>(pp) && pp.size())
1348  {
1349  const auto& cycpp =
1350  refCast<const cyclicPolyPatch>(pp);
1351 
1352  if
1353  (
1354  cycpp.transform()
1356  )
1357  {
1358  // Force to wanted separation
1359  Info<< "On cyclic translation patch "
1360  << pp.name()
1361  << " forcing uniform separation of "
1362  << cycpp.separationVector() << endl;
1363  const_cast<vectorField&>(cpp.separation()) =
1364  pointField(1, cycpp.separationVector());
1365  }
1366  else
1367  {
1368  const auto& nbr = cycpp.neighbPatch();
1369  const_cast<vectorField&>(cpp.separation()) =
1370  pointField
1371  (
1372  1,
1373  nbr[0].centre(mesh.points())
1374  - cycpp[0].centre(mesh.points())
1375  );
1376  }
1377  }
1378  Info<< "On coupled patch " << pp.name()
1379  << " forcing uniform separation of "
1380  << cpp.separation() << endl;
1381  }
1382  else if (!cpp.parallel())
1383  {
1384  Info<< "On coupled patch " << pp.name()
1385  << " forcing uniform rotation of "
1386  << cpp.forwardT()[0] << endl;
1387 
1388  const_cast<tensorField&>
1389  (
1390  cpp.forwardT()
1391  ).setSize(1);
1392  const_cast<tensorField&>
1393  (
1394  cpp.reverseT()
1395  ).setSize(1);
1396 
1397  Info<< "On coupled patch " << pp.name()
1398  << " forcing uniform rotation of "
1399  << cpp.forwardT() << endl;
1400  }
1401  }
1402  }
1403 
1404  Info<< "Synchronising points." << endl;
1405 
1406  pointField newPoints(mesh.points());
1407 
1408  syncPoints
1409  (
1410  mesh,
1411  newPoints,
1413  point(GREAT, GREAT, GREAT)
1414  );
1415 
1416  scalarField diff(mag(newPoints-mesh.points()));
1417  Info<< "Points changed by average:" << gAverage(diff)
1418  << " max:" << gMax(diff) << nl << endl;
1419 
1420  mesh.movePoints(newPoints);
1421  }
1422 
1423 
1424  // 3. Remove zeros-sized patches
1425  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1426 
1427  Info<< "Removing patches with no faces in them." << nl << endl;
1428  const wordList oldPatchNames(mesh.boundaryMesh().names());
1429  const wordList oldPatchTypes(mesh.boundaryMesh().types());
1431  forAll(oldPatchNames, patchi)
1432  {
1433  const word& pName = oldPatchNames[patchi];
1434  if (mesh.boundaryMesh().findPatchID(pName) == -1)
1435  {
1436  Info<< "Removed zero-sized patch " << pName
1437  << " type " << oldPatchTypes[patchi]
1438  << " at position " << patchi << endl;
1439  }
1440  }
1441 
1442 
1443  if (writeObj)
1444  {
1445  dumpCyclicMatch("final_", mesh);
1446  }
1447  }
1448 
1449  if (!overwrite)
1450  {
1451  ++runTime;
1452  }
1453  else
1454  {
1455  forAll(meshes, meshi)
1456  {
1457  fvMesh& mesh = meshes[meshi];
1458 
1459  mesh.setInstance(oldInstances[meshi]);
1460  }
1461  }
1462 
1463  // More precision (for points data)
1465 
1466  // Write resulting mesh
1467  forAll(meshes, meshi)
1468  {
1469  fvMesh& mesh = meshes[meshi];
1470  Info<< "\n\nWriting repatched mesh " << mesh.name()
1471  << " to " << runTime.timeName() << nl << endl;
1472  mesh.clearOut(); // remove meshPhi
1473  mesh.write();
1476  }
1477 
1478  Info<< "End\n" << endl;
1479 
1480  return 0;
1481 }
1482 
1483 
1484 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream. ...
Definition: dictionary.C:367
label nPoints() const
Number of points supporting patch faces.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
dictionary dict
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C: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)
rDeltaTY field()
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:498
A class for handling file names.
Definition: fileName.H:72
A list of face labels.
Definition: faceSet.H:47
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
virtual bool separated() const
Are the planes separated.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Class describing modification of a face.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
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.
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:493
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Output to file stream as an OSstream, normally using std::ofstream for the actual output...
Definition: OFstream.H:71
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const word dictName("faMeshDefinition")
const DynamicList< label > & region() const
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:929
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
Definition: fvMesh.C:227
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C: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:134
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
points setSize(newPointi)
virtual bool parallel() const
Are the cyclic planes parallel.
virtual bool owner() const
Does this side own the patch ?
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition: topoSet.C:693
label nGlobalPoints() const
Return number of globally shared points.
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition: IPstream.H:81
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:320
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1005
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:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:38
Cyclic plane patch.
T & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
Definition: DynamicListI.H: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
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition: IOstream.H:440
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1113
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:671
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:406
const vectorField & faceCentres() const
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:49
void transformList(const tensor &rotTensor, UList< T > &field)
Inplace transform a list of elements.
Definition: transformList.C:48
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
vector point
Point is a vector.
Definition: point.H:37
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
const word & name() const
Return reference to name.
Definition: fvMesh.H:387
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:163
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:75
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:765
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
"buffered" : (MPI_Bsend, MPI_Recv)
const cyclicPolyPatch & neighbPatch() const
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Foam::argList args(argc, argv)
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points list.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Remove zero sized patches. All but processor patches are.
Definition: fvMeshTools.C:363
bool set(const Key &key, const T &obj)
Copy assign a new entry, overwriting existing entries.
Definition: HashTableI.H:174
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Combines List elements. 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.