createNonConformalCouples.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  createNonConformalCouples
26 
27 Description
28  Utility to create non-conformal couples between non-coupled patches.
29 
30 Usage
31  \b createNonConformalCouples <patch1> <patch2>
32 
33  Options:
34  - \par -noOverwrite \n
35  Do not replace the old mesh with the new one, writing the new one
36  into a separate time directory
37 
38  - \par -fields \n
39  Add non-conformal boundary conditions to the fields.
40 
41 Note
42  If run with two arguments, these arguments specify the patches between
43  which a single couple is to be created. The resulting couple will not have
44  a transformation.
45 
46 Usage
47  \b createNonConformalCouples
48 
49 Note
50  If run without arguments then settings are read from a \b
51  system/createNonConformalCouplesDict dictionary (or from a different
52  dictionary specified by the \b -dict option). This dictionary can specify
53  the creation of multiple couples and/or couples with transformations.
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #include "argList.H"
58 #include "fvMeshTools.H"
59 #include "hashedWordList.H"
60 #include "IOobjectList.H"
61 #include "MultiRegionList.H"
66 #include "polyMesh.H"
67 #include "processorPolyPatch.H"
69 #include "systemDict.H"
70 #include "Time.H"
71 
72 #include "ReadFields.H"
73 #include "volFields.H"
74 #include "surfaceFields.H"
75 #include "pointFields.H"
76 
77 using namespace Foam;
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 struct nonConformalCouple
82 {
83  // Public Data
84 
86 
87  Pair<word> origPatchNames;
88 
89  word ncPatchType;
90 
91  Pair<word> ncPatchNames;
92 
93  Pair<dictionary> ncPatchDicts;
94 
95  Pair<dictionary> ncPatchFieldDicts;
96 
98 
99 
100  // Constructors
101 
102  // Construct null (for List)
103  nonConformalCouple()
104  {}
105 
106  // Construct from arguments
107  nonConformalCouple
108  (
109  const word& primaryRegionName,
110  const argList& args
111  )
112  :
113  regionNames(primaryRegionName, primaryRegionName),
114  origPatchNames(args[1], args[2]),
115  ncPatchType(nonConformalCyclicPolyPatch::typeName),
116  ncPatchNames
117  (
118  ncPatchType + "_on_" + args[1],
119  ncPatchType + "_on_" + args[2]
120  ),
121  ncPatchDicts(),
122  ncPatchFieldDicts(),
123  transform(true)
124  {}
125 
126  //- Construct from dictionary
127  nonConformalCouple
128  (
129  const word& primaryRegionName,
130  const dictionary& dict
131  )
132  {
133  const bool haveRegion = dict.found("region");
134  const bool haveRegions = dict.found("regions");
135 
136  if (haveRegion && haveRegions)
137  {
139  << "Regions should be specified either by a \"region\" "
140  << "entry with a single region name (for in-region "
141  << "cyclics), or by a \"regions\" entry with a pair of "
142  << "names (for inter-region mapped walls)"
143  << exit(FatalIOError);
144  }
145 
146  const bool havePatches = dict.found("patches");
147  const bool haveOwnerNeighbour =
148  dict.found("owner") || dict.found("neighbour");
149 
150  if (havePatches == haveOwnerNeighbour)
151  {
153  << "Patches should be specified with either a single "
154  << "\"patches\" entry with a pair of patch names, "
155  << "or with two sub-dictionaries named \"owner\" and "
156  << "\"neighbour\"" << exit(FatalIOError);
157  }
158 
159  if (havePatches)
160  {
161  const word thisRegionName =
162  haveRegion ? dict.lookup<word>("region") : word::null;
163 
164  regionNames =
165  haveRegion ? Pair<word>(thisRegionName, thisRegionName)
166  : haveRegions ? dict.lookup<Pair<word>>("regions")
167  : Pair<word>(primaryRegionName, primaryRegionName);
168  origPatchNames =
169  dict.lookup<Pair<word>>("patches");
170  ncPatchType =
172  (
173  "type",
174  regionNames.first() == regionNames.second()
175  ? nonConformalCyclicPolyPatch::typeName
176  : nonConformalMappedWallPolyPatch::typeName
177  );
178  ncPatchNames =
180  (
181  "names",
182  Pair<word>
183  (
184  dict.dictName() + "_on_" + origPatchNames[0],
185  dict.dictName() + "_on_" + origPatchNames[1]
186  )
187  );
188  forAll(ncPatchDicts, i)
189  {
190  ncPatchDicts[i] = dict;
191  ncPatchDicts[i].remove("region");
192  ncPatchDicts[i].remove("regions");
193  ncPatchDicts[i].remove("patches");
194  ncPatchDicts[i].remove("type");
195  ncPatchDicts[i].remove("names");
196  ncPatchDicts[i].remove(cyclicTransform::keywords);
197  }
198  ncPatchFieldDicts = Pair<dictionary>();
199  transform = cyclicTransform(dict, true);
200  }
201  else
202  {
203  const dictionary& ownerDict = dict.subDict("owner");
204  const dictionary& neighbourDict = dict.subDict("neighbour");
205 
206  regionNames =
207  Pair<word>
208  (
209  ownerDict.lookupOrDefault<word>
210  (
211  "region",
212  primaryRegionName
213  ),
214  neighbourDict.lookupOrDefault<word>
215  (
216  "region",
217  primaryRegionName
218  )
219  );
220  origPatchNames =
221  Pair<word>
222  (
223  ownerDict.lookup<word>("patch"),
224  neighbourDict.lookup<word>("patch")
225  );
226  ncPatchType =
228  (
229  "type",
230  regionNames.first() == regionNames.second()
231  ? nonConformalCyclicPolyPatch::typeName
232  : nonConformalMappedWallPolyPatch::typeName
233  );
234  ncPatchNames =
235  Pair<word>
236  (
237  ownerDict.lookupOrDefault<word>
238  (
239  "name",
240  dict.dictName() + "_on_" + origPatchNames[0]
241  ),
242  neighbourDict.lookupOrDefault<word>
243  (
244  "name",
245  dict.dictName() + "_on_" + origPatchNames[1]
246  )
247  );
248  ncPatchDicts[0] = ownerDict;
249  ncPatchDicts[1] = neighbourDict;
250  forAll(ncPatchDicts, i)
251  {
252  ncPatchDicts[i].remove("region");
253  ncPatchDicts[i].remove("patch");
254  ncPatchDicts[i].remove("name");
255  ncPatchDicts[i].remove("patchFields");
256  }
257  ncPatchFieldDicts =
259  (
260  ownerDict.subOrEmptyDict("patchFields"),
261  neighbourDict.subOrEmptyDict("patchFields")
262  );
263  transform = cyclicTransform(dict, true);
264  }
265  }
266 };
267 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 template<class Type>
272 void evaluateNonConformalProcessorCyclics(const fvMesh& mesh)
273 {
275 
276  forAll(fields, i)
277  {
278  const label nReq = Pstream::nRequests();
279 
281  {
282  typename VolField<Type>::Patch& pf =
283  fields[i].boundaryFieldRef()[patchi];
284 
285  if (isA<nonConformalProcessorCyclicPolyPatch>(pf.patch().patch()))
286  {
287  pf.initEvaluate(Pstream::defaultCommsType);
288  }
289  }
290 
291  if
292  (
295  )
296  {
297  Pstream::waitRequests(nReq);
298  }
299 
301  {
302  typename VolField<Type>::Patch& pf =
303  fields[i].boundaryFieldRef()[patchi];
304 
305  if (isA<nonConformalProcessorCyclicPolyPatch>(pf.patch().patch()))
306  {
307  pf.evaluate(Pstream::defaultCommsType);
308  }
309  }
310  }
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 int main(int argc, char *argv[])
317 {
318  #include "addNoOverwriteOption.H"
319  #include "addMeshOption.H"
320  #include "addRegionOption.H"
321  #include "addDictOption.H"
322 
323  const bool haveArgs = argList::hasArgs(argc, argv);
324  if (haveArgs)
325  {
326  argList::validArgs.append("patch1");
327  argList::validArgs.append("patch2");
329  (
330  "fields",
331  "add non-conformal boundary conditions to the fields"
332  );
333  }
334 
335  #include "setRootCase.H"
336  #include "setMeshPath.H"
338 
339  const Foam::word primaryRegionName =
341 
342  // Flag to determine whether or not patches are added to fields
343  bool fields;
344 
345  // Read the couples to be created from arguments or a system dictionary
346  List<nonConformalCouple> couples;
347  if (haveArgs)
348  {
349  fields = args.optionFound("fields");
350 
351  couples.append(nonConformalCouple(primaryRegionName, args));
352  }
353  else
354  {
355  const dictionary dict
356  (
357  systemDict
358  (
359  "createNonConformalCouplesDict",
360  args,
361  runTime,
362  primaryRegionName
363  )
364  );
365 
366  fields = dict.lookupOrDefault<bool>("fields", false);
367 
368  const dictionary& couplesDict =
369  dict.optionalSubDict("nonConformalCouples");
370 
371  forAllConstIter(dictionary, couplesDict, iter)
372  {
373  if (!iter().isDict()) continue;
374 
375  const dictionary& subDict = iter().dict();
376 
377  const bool havePatches = subDict.found("patches");
378  const bool haveOwnerNeighbour =
379  subDict.found("owner") || subDict.found("neighbour");
380 
381  if (havePatches == haveOwnerNeighbour)
382  {
383  FatalIOErrorInFunction(subDict)
384  << "Patches should be specified with either a single "
385  << "\"patches\" entry with a pair of patch names, "
386  << "or with two sub-dictionaries named \"owner\" and "
387  << "\"neighbour\"." << exit(FatalIOError);
388  }
389 
390  couples.append(nonConformalCouple(primaryRegionName, subDict));
391  }
392  }
393 
394  // Determine if operating on meshes within a sub-path
395  const bool haveMeshPath = meshPath != word::null;
396 
397  // Create all the meshes needed
399  PtrList<fvMesh> regionMeshes;
400  forAll(couples, i)
401  {
402  forAll(couples[i].regionNames, sidei)
403  {
404  const word& regionName = couples[i].regionNames[sidei];
405 
406  if (regionNames.found(regionName)) continue;
407 
408  const IOobject regionMeshIo
409  (
410  regionName,
411  runTime.name(),
412  meshPath,
413  runTime,
415  );
416 
417  // If we are creating couples on a mesh within a mesh path
418  // sub-directory then these couples will not be stitched so loading
419  // the neighbouring regions is optional
420  if (haveMeshPath && !polyMesh::found(regionMeshIo)) continue;
421 
423  regionMeshes.append(new fvMesh(regionMeshIo, false));
424  }
425  }
426 
427  // Early exit if there is nothing to do
428  if (regionMeshes.empty())
429  {
430  Info<< "Nothing to be done" << nl << endl
431  << "End" << nl << endl;
432 
433  return 0;
434  }
435 
436  #include "setNoOverwrite.H"
437 
438  const word oldInstance = regionMeshes[0].pointsInstance();
439 
440  // Read the fields
441  if (fields) Info<< "Reading geometric fields" << nl << endl;
442 
443  #define DeclareRegionGeoTypeFields(Type, Geo) \
444  PtrList<PtrList<Geo##Field<Type>>> \
445  CAT4(region, Geo, CAPITALIZE(Type), Fields)(regionMeshes.size()); \
446  forAll(regionMeshes, regioni) \
447  { \
448  CAT4(region, Geo, CAPITALIZE(Type), Fields).set \
449  ( \
450  regioni, \
451  new PtrList<Geo##Field<Type>> \
452  ); \
453  }
454  FOR_ALL_FIELD_TYPES(DeclareRegionGeoTypeFields, Vol);
455  FOR_ALL_FIELD_TYPES(DeclareRegionGeoTypeFields, Surface);
456  FOR_ALL_FIELD_TYPES(DeclareRegionGeoTypeFields, Point);
457  #undef DeclareRegionGeoTypeFields
458 
459  if (fields)
460  {
461  MultiRegionUList<fvMesh> multiRegionMeshes(regionMeshes, false);
462 
463  forAll(regionMeshes, regioni)
464  {
465  RegionRef<fvMesh> mesh = multiRegionMeshes[regioni];
466 
467  IOobjectList objects(mesh, runTime.name());
468 
469  #define ReadRegionGeoTypeFields(Type, Geo, mesh) \
470  ReadFields \
471  ( \
472  mesh, \
473  objects, \
474  CAT4(region, Geo, CAPITALIZE(Type), Fields)[regioni] \
475  );
476  FOR_ALL_FIELD_TYPES(ReadRegionGeoTypeFields, Vol, mesh);
477  FOR_ALL_FIELD_TYPES(ReadRegionGeoTypeFields, Surface, mesh);
478  const pointMesh& pMesh = pointMesh::New(mesh);
479  FOR_ALL_FIELD_TYPES(ReadRegionGeoTypeFields, Point, pMesh);
480  #undef ReadRegionGeoTypeFields
481  }
482 
483  Info<< endl;
484  }
485 
486  if (!overwrite)
487  {
488  runTime++;
489  }
490 
491  // Make sure the meshes are not connected before couples are added
492  forAll(regionMeshes, regioni)
493  {
494  regionMeshes[regioni].conform();
495  }
496 
497  // Find the first processor patch and face
498  labelList regionFirstProcPatchis(regionMeshes.size());
499  labelList regionFirstProcFaceis(regionMeshes.size());;
500  forAll(regionMeshes, regioni)
501  {
502  const fvMesh& mesh = regionMeshes[regioni];
504  label& firstProcPatchi = regionFirstProcPatchis[regioni];
505  label& firstProcFacei = regionFirstProcFaceis[regioni];
506 
507  firstProcPatchi = patches.size();
508  firstProcFacei = mesh.nFaces();
509 
511  {
512  const polyPatch& pp = patches[patchi];
513 
514  const bool isProcPp = isA<processorPolyPatch>(pp);
515 
516  if (isProcPp && firstProcPatchi == patches.size())
517  {
518  firstProcPatchi = patchi;
519  firstProcFacei = pp.start();
520  }
521 
522  if (!isProcPp && firstProcPatchi != patches.size())
523  {
525  << "Processor patches of region " << regionNames[regioni]
526  << " do not follow boundary patches"
527  << exit(FatalError);
528  }
529  }
530  }
531 
532  // Start building lists of patches and patch-fields to add
533  List<List<polyPatch*>> newPatches(regionMeshes.size());
534  List<boolList> newPatchIsCouple(regionMeshes.size());
535  List<List<dictionary>> newPatchFieldDicts(regionMeshes.size());
536 
537  // Clone the non-processor patches
538  forAll(regionMeshes, regioni)
539  {
540  const fvMesh& mesh = regionMeshes[regioni];
542  const label firstProcPatchi = regionFirstProcPatchis[regioni];
543 
544  for (label patchi = 0; patchi < firstProcPatchi; ++ patchi)
545  {
546  const polyPatch& pp = patches[patchi];
547 
548  newPatches[regioni].append
549  (
550  pp.clone
551  (
552  patches,
553  patchi,
554  pp.size(),
555  pp.start()
556  ).ptr()
557  );
558  newPatchIsCouple[regioni].append(false);
559  newPatchFieldDicts[regioni].append(dictionary());
560  }
561  }
562 
563  // Add the non-processor coupled patches
564  forAll(couples, couplei)
565  {
566  const nonConformalCouple& couple = couples[couplei];
567 
568  Info<< indent << "Adding " << couple.ncPatchType
569  << " interfaces between patches: " << incrIndent << nl
570  << indent << couple.origPatchNames << decrIndent << nl
571  << "In regions: " << incrIndent << nl
572  << indent << couple.regionNames << decrIndent << nl
573  << indent << "Named:" << incrIndent << nl
574  << indent << couple.ncPatchNames << decrIndent << nl
575  << indent << "With transform: " << incrIndent << nl;
576  couple.transform.write(Info);
577  Info<< decrIndent << nl;
578 
579  auto appendPatch = [&](const bool owner)
580  {
581  if (!regionNames.found(couple.regionNames[!owner])) return;
582 
583  const label regioni =
584  regionNames[couple.regionNames[!owner]];
585 
586  dictionary patchDict = dictionary::entries
587  (
588  "type", couple.ncPatchType,
589  "nFaces", 0,
590  "startFace", regionFirstProcFaceis[regioni],
591  "originalPatch", couple.origPatchNames[!owner],
592  "neighbourRegion", couple.regionNames[owner],
593  "neighbourPatch", couple.ncPatchNames[owner],
594  "owner", owner
595  );
596 
597  {
598  OStringStream oss;
599  (owner ? couple.transform : inv(couple.transform)).write(oss);
600  patchDict.merge(IStringStream(oss.str())());
601  }
602 
603  patchDict.merge(couple.ncPatchDicts[!owner]);
604 
605  newPatches[regioni].append
606  (
608  (
609  couple.ncPatchNames[!owner],
610  patchDict,
611  newPatches[regioni].size(),
612  regionMeshes[regioni].boundaryMesh()
613  ).ptr()
614  );
615  newPatchIsCouple[regioni].append(true);
616  newPatchFieldDicts[regioni].append
617  (
618  couple.ncPatchFieldDicts[!owner]
619  );
620  };
621  appendPatch(true);
622  appendPatch(false);
623  }
624 
625  // Add the error patches. Note there is only one for each original patch,
626  // regardless of how many couplings are attached to that patch.
627  {
628  // Create a table of unique original patches
629  HashTable<label, Pair<word>, Hash<Pair<word>>> regionOrigPatchToIndex;
630  DynamicList<Pair<word>> regionOrigPatches;
631  forAll(couples, couplei)
632  {
633  const nonConformalCouple& couple = couples[couplei];
634  forAll(couple.regionNames, i)
635  {
636  const Pair<word> regionOrigPatch
637  (
638  couple.regionNames[i],
639  couple.origPatchNames[i]
640  );
641 
642  if (!regionOrigPatchToIndex.found(regionOrigPatch))
643  {
644  regionOrigPatchToIndex.insert
645  (
646  regionOrigPatch,
647  regionOrigPatches.size()
648  );
649  regionOrigPatches.append(regionOrigPatch);
650  }
651  }
652  }
653 
654  // Add an error patch for each unique original patch
655  forAll(regionOrigPatches, i)
656  {
657  const word& regionName = regionOrigPatches[i].first();
658 
659  if (!regionNames.found(regionName)) continue;
660 
661  const label regioni = regionNames[regionName];
662  const word& origPatchName = regionOrigPatches[i].second();
663 
664  newPatches[regioni].append
665  (
667  (
668  nonConformalErrorPolyPatch::typeName
669  + "_on_"
670  + origPatchName,
671  0,
672  regionFirstProcFaceis[regioni],
673  newPatches[regioni].size(),
674  regionMeshes[regioni].boundaryMesh(),
675  nonConformalErrorPolyPatch::typeName,
676  origPatchName
677  )
678  );
679  newPatchIsCouple[regioni].append(false);
680  newPatchFieldDicts[regioni].append(dictionary());
681  }
682  }
683 
684  // Clone the processor patches
685  forAll(regionMeshes, regioni)
686  {
687  const fvMesh& mesh = regionMeshes[regioni];
689  const label firstProcPatchi = regionFirstProcPatchis[regioni];
690 
691  for (label patchi = firstProcPatchi; patchi < patches.size(); ++ patchi)
692  {
693  const polyPatch& pp = patches[patchi];
694 
695  newPatches[regioni].append
696  (
697  pp.clone
698  (
699  patches,
700  newPatches[regioni].size(),
701  pp.size(),
702  pp.start()
703  ).ptr()
704  );
705  newPatchIsCouple[regioni].append(false);
706  newPatchFieldDicts[regioni].append(dictionary());
707  }
708  }
709 
710  // Add the processor cyclic patches
711  if (Pstream::parRun())
712  {
713  forAll(couples, couplei)
714  {
715  const nonConformalCouple& couple = couples[couplei];
716 
717  if (couple.ncPatchType != nonConformalCyclicPolyPatch::typeName)
718  continue;
719 
720  if (!regionNames.found(couples[couplei].regionNames[0]))
721  continue;
722 
723  const label regioni = regionNames[couples[couplei].regionNames[0]];
724  const fvMesh& mesh = regionMeshes[regioni];
726 
727  const polyPatch& patch1 = patches[couple.origPatchNames[0]];
728  const polyPatch& patch2 = patches[couple.origPatchNames[1]];
729 
730  boolList procHasPatch1(Pstream::nProcs(), false);
731  procHasPatch1[Pstream::myProcNo()] = !patch1.empty();
732  Pstream::gatherList(procHasPatch1);
733  Pstream::scatterList(procHasPatch1);
734 
735  boolList procHasPatch2(Pstream::nProcs(), false);
736  procHasPatch2[Pstream::myProcNo()] = !patch2.empty();
737  Pstream::gatherList(procHasPatch2);
738  Pstream::scatterList(procHasPatch2);
739 
740  // Multiple cyclic interfaces must be ordered in a specific way for
741  // processor communication to function correctly.
742  //
743  // A communication that is sent from the cyclic owner is received
744  // on the cyclic neighbour and vice versa. Therefore, in a coupled
745  // pair of processors if one sends the owner first the other must
746  // receive the neighbour first.
747  //
748  // We ensure the above by ordering the patches so that for the
749  // lower indexed processor the owner interface comes first, and for
750  // the higher indexed processor the neighbour comes first.
751 
752  auto appendProcPatches = [&](const bool owner, const bool first)
753  {
754  const boolList& procHasPatchA =
755  owner ? procHasPatch1 : procHasPatch2;
756  const boolList& procHasPatchB =
757  owner ? procHasPatch2 : procHasPatch1;
758 
759  if (procHasPatchA[Pstream::myProcNo()])
760  {
761  forAll(procHasPatchB, proci)
762  {
763  if
764  (
765  (
766  (first && proci > Pstream::myProcNo())
767  || (!first && proci < Pstream::myProcNo())
768  )
769  && procHasPatchB[proci]
770  )
771  {
772  newPatches[regioni].append
773  (
775  (
776  0,
777  mesh.nFaces(),
778  newPatches.size(),
779  patches,
781  proci,
782  couple.ncPatchNames[!owner],
783  couple.origPatchNames[!owner]
784  )
785  );
786  newPatchIsCouple[regioni].append(true);
787  newPatchFieldDicts[regioni].append
788  (
789  couple.ncPatchFieldDicts[!owner]
790  );
791  }
792  }
793  }
794  };
795 
796  appendProcPatches(true, true);
797  appendProcPatches(false, true);
798  appendProcPatches(false, false);
799  appendProcPatches(true, false);
800  }
801  }
802 
803  // Re-patch the meshes. Create constraint or calculated patch fields at
804  // this stage; don't apply the patch field dictionaries. The patch fields
805  // will be specified properly later after all patches have been added and
806  // the meshes have been stitched. That way the geometry is fully available
807  // for construction of the patch fields.
808  forAll(regionMeshes, regioni)
809  {
810  forAll(newPatches[regioni], newPatchi)
811  {
813  (
814  regionMeshes[regioni],
815  *newPatches[regioni][newPatchi]
816  );
817  }
818  }
819 
820  // Connect the meshes
821  {
822  MultiRegionUList<fvMesh> multiRegionMeshes(regionMeshes, false);
823  forAll(regionMeshes, regioni)
824  {
825  RegionRef<fvMesh> mesh = multiRegionMeshes[regioni];
826  fvMeshStitchers::stationary(mesh).connect(false, false, false);
827  }
828 
829  if (!haveMeshPath) Info<< endl;
830  }
831 
832  // Set the fields on the new patches. This allows constraint types (e.g.,
833  // nonConformalCyclic) to be set by default, but requires a dictionary
834  // setting for non-constrained types (e.g., nonConformalMappedWall). It
835  // also allows constraint types to be overridden (e.g., with jumpCyclic) if
836  // there is a dictionary present.
837  forAll(regionMeshes, regioni)
838  {
839  forAll(newPatches[regioni], newPatchi)
840  {
841  if (newPatchIsCouple[regioni][newPatchi])
842  {
843  fvMeshTools::setPatchFields
844  (
845  regionMeshes[regioni],
846  newPatchi,
847  newPatchFieldDicts[regioni][newPatchi]
848  );
849  }
850  }
851  }
852 
853  // Communicate values across non-conformal processor cyclics so that they
854  // contain valid values that can be written to disk
855  if (Pstream::parRun())
856  {
857  forAll(regionMeshes, regioni)
858  {
859  const fvMesh& mesh = regionMeshes[regioni];
860 
861  #define EVALUATE_NON_CONFORMAL_PROCESSOR_CYCLICS(Type, nullArg) \
862  evaluateNonConformalProcessorCyclics<Type>(mesh);
863  FOR_ALL_FIELD_TYPES(EVALUATE_NON_CONFORMAL_PROCESSOR_CYCLICS)
864  #undef EVALUATE_NON_CONFORMAL_PROCESSOR_CYCLICS
865  }
866  }
867 
868  // Set the precision of the points data to 10
870 
871  // Set the instance so that the mesh writes
872  const word newInstance = overwrite ? oldInstance : runTime.name();
873  forAll(regionMeshes, regioni)
874  {
875  regionMeshes[regioni].setInstance(newInstance);
876  regionMeshes[regioni].setPolyFacesBfInstance(newInstance);
877  }
878 
879  // Write resulting mesh
880  Info<< "Writing mesh to " << runTime.name() << nl << endl;
881  forAll(regionMeshes, regioni)
882  {
883  regionMeshes[regioni].write();
884  }
885 
886  Info<< "End" << nl << endl;
887 
888  return 0;
889 }
890 
891 
892 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Generic GeometricField class.
GeoMesh::template PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
An STL-conforming hash table.
Definition: HashTable.H:127
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:53
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:473
Input from memory buffer stream.
Definition: IStringStream.H:52
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Output to memory buffer stream.
Definition: OStringStream.H:52
string str() const
Return the string.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
T & first()
Return the first element of the list.
Definition: UListI.H:114
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
bool empty() const
Return true if the UPtrList is empty (ie, size() is zero)
Definition: UPtrListI.H:36
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:103
static bool hasArgs(int argc, char *argv[])
Return true if there are arguments.
Definition: argList.C:298
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:294
Cyclic plane transformation.
static const wordList keywords
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary.
Definition: dictionary.C:900
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:927
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
bool merge(const dictionary &)
Merge entries from the given dictionary.
Definition: dictionary.C:1314
static std::tuple< const Entries &... > entries(const Entries &...)
Construct an entries tuple from which to make a dictionary.
bool connect(const bool changing, const bool geometric, const bool load)
Connect the mesh by adding faces into the nonConformalCyclics.
Mesh stitcher for stationary meshes.
static label addPatch(fvMesh &mesh, const polyPatch &patch)
Add patch. Inserts patch before all processor patches. Returns the.
Definition: fvMeshTools.C:33
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
UPtrList< GeoField > fields(bool strict=false, const HashSet< word > &geometryFields=fvMesh::geometryFields) const
Return the list of fields of type GeoField.
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
A wordList with hashed indices for faster lookup by name.
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
Non-conformal error poly patch. As nonConformalPolyPatch. This patch is a non-coupled non-conformal p...
Non-conformal processor cyclic poly patch. As nonConformalCyclicPolyPatch, but the neighbouring patch...
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
Foam::polyBoundaryMesh.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:270
static bool found(const IOobject &io)
Return whether the given IOobject relates to a mesh on disk.
Definition: polyMesh.C:959
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:215
label nFaces() const
Read and return the specified dictionary from system or from path provided with the -dict option.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const fvPatchList & patches
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:234
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const word & regionName(const solver &region)
Definition: solver.H:218
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:501
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
objects
dictionary dict
word meshPath
Definition: setMeshPath.H:1
const bool overwrite
Definition: setNoOverwrite.H:1
const Foam::wordList regionNames(args.optionFound("allRegions") ? runTime .controlDict().subDict("regionSolvers").toc() :wordList(1, args.optionFound("region") ? args.optionRead< word >("region") :polyMesh::defaultRegion))
Foam::argList args(argc, argv)
Foam::surfaceFields.