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-2024 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 -overwrite \n
35  Replace the old mesh with the new one, rather than 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"
59 #include "fvMeshTools.H"
60 #include "hashedWordList.H"
61 #include "IOobjectList.H"
62 #include "MultiRegionList.H"
67 #include "polyMesh.H"
68 #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> ncPatchFieldDicts;
94 
96 
97 
98  // Constructors
99 
100  // Construct null (for List)
101  nonConformalCouple()
102  {}
103 
104  // Construct from arguments
105  nonConformalCouple
106  (
107  const word& primaryRegionName,
108  const argList& args
109  )
110  :
111  regionNames(primaryRegionName, primaryRegionName),
112  origPatchNames(args[1], args[2]),
113  ncPatchType(nonConformalCyclicPolyPatch::typeName),
114  ncPatchNames
115  (
116  ncPatchType + "_on_" + args[1],
117  ncPatchType + "_on_" + args[2]
118  ),
119  ncPatchFieldDicts(),
120  transform(true)
121  {}
122 
123  //- Construct from dictionary
124  nonConformalCouple
125  (
126  const word& primaryRegionName,
127  const dictionary& dict
128  )
129  {
130  const bool haveRegion = dict.found("region");
131  const bool haveRegions = dict.found("regions");
132 
133  if (haveRegion && haveRegions)
134  {
136  << "Regions should be specified either by a \"region\" "
137  << "entry with a single region name (for in-region "
138  << "cyclics), or by a \"regions\" entry with a pair of "
139  << "names (for inter-region mapped walls)"
140  << exit(FatalIOError);
141  }
142 
143  const bool havePatches = dict.found("patches");
144  const bool haveOwnerNeighbour =
145  dict.found("owner") || dict.found("neighbour");
146 
147  if (havePatches == haveOwnerNeighbour)
148  {
150  << "Patches should be specified with either a single "
151  << "\"patches\" entry with a pair of patch names, "
152  << "or with two sub-dictionaries named \"owner\" and "
153  << "\"neighbour\"" << exit(FatalIOError);
154  }
155 
156  if (havePatches)
157  {
158  const word thisRegionName =
159  haveRegion ? dict.lookup<word>("region") : word::null;
160 
161  regionNames =
162  haveRegion ? Pair<word>(thisRegionName, thisRegionName)
163  : haveRegions ? dict.lookup<Pair<word>>("regions")
164  : Pair<word>(primaryRegionName, primaryRegionName);
165  origPatchNames =
166  dict.lookup<Pair<word>>("patches");
167  ncPatchType =
169  (
170  "type",
171  regionNames.first() == regionNames.second()
172  ? nonConformalCyclicPolyPatch::typeName
173  : nonConformalMappedWallPolyPatch::typeName
174  );
175  ncPatchNames =
176  Pair<word>
177  (
178  dict.dictName() + "_on_" + origPatchNames[0],
179  dict.dictName() + "_on_" + origPatchNames[1]
180  );
181  ncPatchFieldDicts = Pair<dictionary>();
182  transform = cyclicTransform(dict, true);
183  }
184  else
185  {
186  const dictionary& ownerDict = dict.subDict("owner");
187  const dictionary& neighbourDict = dict.subDict("neighbour");
188 
189  regionNames =
190  Pair<word>
191  (
192  ownerDict.lookupOrDefault<word>
193  (
194  "region",
195  primaryRegionName
196  ),
197  neighbourDict.lookupOrDefault<word>
198  (
199  "region",
200  primaryRegionName
201  )
202  );
203  origPatchNames =
204  Pair<word>
205  (
206  ownerDict.lookup<word>("patch"),
207  neighbourDict.lookup<word>("patch")
208  );
209  ncPatchType =
211  (
212  "type",
213  regionNames.first() == regionNames.second()
214  ? nonConformalCyclicPolyPatch::typeName
215  : nonConformalMappedWallPolyPatch::typeName
216  );
217  ncPatchNames =
218  Pair<word>
219  (
220  ownerDict.lookupOrDefault<word>
221  (
222  "name",
223  dict.dictName() + "_on_" + origPatchNames[0]
224  ),
225  neighbourDict.lookupOrDefault<word>
226  (
227  "name",
228  dict.dictName() + "_on_" + origPatchNames[1]
229  )
230  );
231  ncPatchFieldDicts =
233  (
234  ownerDict.subOrEmptyDict("patchFields"),
235  neighbourDict.subOrEmptyDict("patchFields")
236  );
237  transform = cyclicTransform(dict, true);
238  }
239  }
240 };
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 template<class Type>
246 void evaluateNonConformalProcessorCyclics(const fvMesh& mesh)
247 {
249 
250  forAll(fields, i)
251  {
252  const label nReq = Pstream::nRequests();
253 
254  forAll(mesh.boundary(), patchi)
255  {
256  typename VolField<Type>::Patch& pf =
257  fields[i].boundaryFieldRef()[patchi];
258 
259  if (isA<nonConformalProcessorCyclicPolyPatch>(pf.patch().patch()))
260  {
261  pf.initEvaluate(Pstream::defaultCommsType);
262  }
263  }
264 
265  if
266  (
269  )
270  {
271  Pstream::waitRequests(nReq);
272  }
273 
274  forAll(mesh.boundary(), patchi)
275  {
276  typename VolField<Type>::Patch& pf =
277  fields[i].boundaryFieldRef()[patchi];
278 
279  if (isA<nonConformalProcessorCyclicPolyPatch>(pf.patch().patch()))
280  {
281  pf.evaluate(Pstream::defaultCommsType);
282  }
283  }
284  }
285 }
286 
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 int main(int argc, char *argv[])
291 {
292  #include "addOverwriteOption.H"
293  #include "addRegionOption.H"
294  #include "addDictOption.H"
295 
296  const bool haveArgs = argList::hasArgs(argc, argv);
297  if (haveArgs)
298  {
299  argList::validArgs.append("patch1");
300  argList::validArgs.append("patch2");
302  (
303  "fields",
304  "add non-conformal boundary conditions to the fields"
305  );
306  }
307 
308  #include "setRootCase.H"
310 
311  const Foam::word primaryRegionName =
313 
314  // Flag to determine whether or not patches are added to fields
315  bool fields;
316 
317  // Read the couples to be created from arguments or a system dictionary
318  List<nonConformalCouple> couples;
319  if (haveArgs)
320  {
321  fields = args.optionFound("fields");
322 
323  couples.append(nonConformalCouple(primaryRegionName, args));
324  }
325  else
326  {
327  const dictionary dict
328  (
329  systemDict
330  (
331  "createNonConformalCouplesDict",
332  args,
333  runTime,
334  primaryRegionName
335  )
336  );
337 
338  fields = dict.lookupOrDefault<bool>("fields", false);
339 
340  const dictionary& couplesDict =
341  dict.optionalSubDict("nonConformalCouples");
342 
343  forAllConstIter(dictionary, couplesDict, iter)
344  {
345  if (!iter().isDict()) continue;
346 
347  const dictionary& subDict = iter().dict();
348 
349  const bool havePatches = subDict.found("patches");
350  const bool haveOwnerNeighbour =
351  subDict.found("owner") || subDict.found("neighbour");
352 
353  if (havePatches == haveOwnerNeighbour)
354  {
355  FatalIOErrorInFunction(subDict)
356  << "Patches should be specified with either a single "
357  << "\"patches\" entry with a pair of patch names, "
358  << "or with two sub-dictionaries named \"owner\" and "
359  << "\"neighbour\"." << exit(FatalIOError);
360  }
361 
362  couples.append(nonConformalCouple(primaryRegionName, subDict));
363  }
364  }
365 
366  // Create all the meshes needed
368  PtrList<fvMesh> regionMeshes;
369  forAll(couples, i)
370  {
371  forAll(couples[i].regionNames, sidei)
372  {
373  const word& regionName = couples[i].regionNames[sidei];
374 
375  if (regionNames.found(regionName)) continue;
376 
378  regionMeshes.append
379  (
380  new fvMesh
381  (
382  IOobject
383  (
384  regionName,
385  runTime.name(),
386  runTime,
388  ),
389  false
390  )
391  );
392  }
393  }
394 
395  const bool overwrite = args.optionFound("overwrite");
396 
397  const word oldInstance = regionMeshes[0].pointsInstance();
398 
399  // Read the fields
400  if (fields) Info<< "Reading geometric fields" << nl << endl;
401 
402  #define DeclareRegionGeoTypeFields(Type, Geo) \
403  PtrList<PtrList<Geo##Field<Type>>> \
404  CAT4(region, Geo, CAPITALIZE(Type), Fields)(regionMeshes.size()); \
405  forAll(regionMeshes, regioni) \
406  { \
407  CAT4(region, Geo, CAPITALIZE(Type), Fields).set \
408  ( \
409  regioni, \
410  new PtrList<Geo##Field<Type>> \
411  ); \
412  }
413  FOR_ALL_FIELD_TYPES(DeclareRegionGeoTypeFields, Vol);
414  FOR_ALL_FIELD_TYPES(DeclareRegionGeoTypeFields, Surface);
415  FOR_ALL_FIELD_TYPES(DeclareRegionGeoTypeFields, Point);
416  #undef DeclareRegionGeoTypeFields
417 
418  if (fields)
419  {
420  MultiRegionUList<fvMesh> multiRegionMeshes(regionMeshes, false);
421 
422  forAll(regionMeshes, regioni)
423  {
424  RegionRef<fvMesh> mesh = multiRegionMeshes[regioni];
425 
426  IOobjectList objects(mesh, runTime.name());
427 
428  #define ReadRegionGeoTypeFields(Type, Geo, mesh) \
429  ReadFields \
430  ( \
431  mesh, \
432  objects, \
433  CAT4(region, Geo, CAPITALIZE(Type), Fields)[regioni] \
434  );
435  FOR_ALL_FIELD_TYPES(ReadRegionGeoTypeFields, Vol, mesh);
436  FOR_ALL_FIELD_TYPES(ReadRegionGeoTypeFields, Surface, mesh);
437  const pointMesh& pMesh = pointMesh::New(mesh);
438  FOR_ALL_FIELD_TYPES(ReadRegionGeoTypeFields, Point, pMesh);
439  #undef ReadRegionGeoTypeFields
440  }
441 
442  Info<< endl;
443  }
444 
445  if (!overwrite)
446  {
447  runTime++;
448  }
449 
450  // Make sure the meshes are not connected before couples are added
451  forAll(regionMeshes, regioni)
452  {
453  regionMeshes[regioni].conform();
454  }
455 
456  // Find the first processor patch and face
457  labelList regionFirstProcPatchis(regionMeshes.size());
458  labelList regionFirstProcFaceis(regionMeshes.size());;
459  forAll(regionMeshes, regioni)
460  {
461  const fvMesh& mesh = regionMeshes[regioni];
462  const polyBoundaryMesh& patches = mesh.boundaryMesh();
463  label& firstProcPatchi = regionFirstProcPatchis[regioni];
464  label& firstProcFacei = regionFirstProcFaceis[regioni];
465 
466  firstProcPatchi = patches.size();
467  firstProcFacei = mesh.nFaces();
468 
470  {
471  const polyPatch& pp = patches[patchi];
472 
473  const bool isProcPp = isA<processorPolyPatch>(pp);
474 
475  if (isProcPp && firstProcPatchi == patches.size())
476  {
477  firstProcPatchi = patchi;
478  firstProcFacei = pp.start();
479  }
480 
481  if (!isProcPp && firstProcPatchi != patches.size())
482  {
484  << "Processor patches of region " << regionNames[regioni]
485  << " do not follow boundary patches"
486  << exit(FatalError);
487  }
488  }
489  }
490 
491  // Start building lists of patches and patch-fields to add
492  List<List<polyPatch*>> newPatches(regionMeshes.size());
493  List<boolList> newPatchIsCouple(regionMeshes.size());
494  List<List<dictionary>> newPatchFieldDicts(regionMeshes.size());
495 
496  // Clone the non-processor patches
497  forAll(regionMeshes, regioni)
498  {
499  const fvMesh& mesh = regionMeshes[regioni];
500  const polyBoundaryMesh& patches = mesh.boundaryMesh();
501  const label firstProcPatchi = regionFirstProcPatchis[regioni];
502 
503  for (label patchi = 0; patchi < firstProcPatchi; ++ patchi)
504  {
505  const polyPatch& pp = patches[patchi];
506 
507  newPatches[regioni].append
508  (
509  pp.clone
510  (
511  patches,
512  patchi,
513  pp.size(),
514  pp.start()
515  ).ptr()
516  );
517  newPatchIsCouple[regioni].append(false);
518  newPatchFieldDicts[regioni].append(dictionary());
519  }
520  }
521 
522  // Add the non-processor coupled patches
523  forAll(couples, couplei)
524  {
525  const nonConformalCouple& couple = couples[couplei];
526 
527  Info<< indent << "Adding " << couple.ncPatchType
528  << " interfaces between patches: " << incrIndent << nl
529  << indent << couple.origPatchNames << decrIndent << nl
530  << "In regions: " << incrIndent << nl
531  << indent << couple.regionNames << decrIndent << nl
532  << indent << "Named:" << incrIndent << nl
533  << indent << couple.ncPatchNames << decrIndent << nl
534  << indent << "With transform: " << incrIndent << nl;
535  couple.transform.write(Info);
536  Info<< decrIndent << nl;
537 
538  auto appendPatch = [&](const bool owner)
539  {
540  const label regioni =
541  regionNames[couple.regionNames[!owner]];
542 
543  dictionary patchDict
544  (
545  "type", couple.ncPatchType,
546  "nFaces", 0,
547  "startFace", regionFirstProcFaceis[regioni],
548  "originalPatch", couple.origPatchNames[!owner],
549  "neighbourRegion", couple.regionNames[owner],
550  "neighbourPatch", couple.ncPatchNames[owner]
551  );
552 
553  {
554  OStringStream oss;
555  (owner ? couple.transform : inv(couple.transform)).write(oss);
556  patchDict.merge(IStringStream(oss.str())());
557  }
558 
559  newPatches[regioni].append
560  (
562  (
563  couple.ncPatchNames[!owner],
564  patchDict,
565  newPatches[regioni].size(),
566  regionMeshes[regioni].boundaryMesh()
567  ).ptr()
568  );
569  newPatchIsCouple[regioni].append(true);
570  newPatchFieldDicts[regioni].append
571  (
572  couple.ncPatchFieldDicts[!owner]
573  );
574  };
575  appendPatch(true);
576  appendPatch(false);
577  }
578 
579  // Add the error patches. Note there is only one for each original patch,
580  // regardless of how many couplings are attached to that patch.
581  {
582  // Create a table of unique original patches
583  HashTable<label, Pair<word>, Hash<Pair<word>>> regionOrigPatchToIndex;
584  DynamicList<Pair<word>> regionOrigPatches;
585  forAll(couples, couplei)
586  {
587  const nonConformalCouple& couple = couples[couplei];
588  forAll(couple.regionNames, i)
589  {
590  const Pair<word> regionOrigPatch
591  (
592  couple.regionNames[i],
593  couple.origPatchNames[i]
594  );
595 
596  if (!regionOrigPatchToIndex.found(regionOrigPatch))
597  {
598  regionOrigPatchToIndex.insert
599  (
600  regionOrigPatch,
601  regionOrigPatches.size()
602  );
603  regionOrigPatches.append(regionOrigPatch);
604  }
605  }
606  }
607 
608  // Add an error patch for each unique original patch
609  forAll(regionOrigPatches, i)
610  {
611  const word& regionName = regionOrigPatches[i].first();
612  const label regioni = regionNames[regionName];
613  const word& origPatchName = regionOrigPatches[i].second();
614 
615  newPatches[regioni].append
616  (
618  (
619  nonConformalErrorPolyPatch::typeName
620  + "_on_"
621  + origPatchName,
622  0,
623  regionFirstProcFaceis[regioni],
624  newPatches[regioni].size(),
625  regionMeshes[regioni].boundaryMesh(),
626  nonConformalErrorPolyPatch::typeName,
627  origPatchName
628  )
629  );
630  newPatchIsCouple[regioni].append(false);
631  newPatchFieldDicts[regioni].append(dictionary());
632  }
633  }
634 
635  // Clone the processor patches
636  forAll(regionMeshes, regioni)
637  {
638  const fvMesh& mesh = regionMeshes[regioni];
639  const polyBoundaryMesh& patches = mesh.boundaryMesh();
640  const label firstProcPatchi = regionFirstProcPatchis[regioni];
641 
642  for (label patchi = firstProcPatchi; patchi < patches.size(); ++ patchi)
643  {
644  const polyPatch& pp = patches[patchi];
645 
646  newPatches[regioni].append
647  (
648  pp.clone
649  (
650  patches,
651  newPatches[regioni].size(),
652  pp.size(),
653  pp.start()
654  ).ptr()
655  );
656  newPatchIsCouple[regioni].append(false);
657  newPatchFieldDicts[regioni].append(dictionary());
658  }
659  }
660 
661  // Add the processor cyclic patches
662  if (Pstream::parRun())
663  {
664  forAll(couples, couplei)
665  {
666  const nonConformalCouple& couple = couples[couplei];
667 
668  if (couple.ncPatchType != nonConformalCyclicPolyPatch::typeName)
669  continue;
670 
671  const label regioni = regionNames[couples[couplei].regionNames[0]];
672  const fvMesh& mesh = regionMeshes[regioni];
673  const polyBoundaryMesh& patches = mesh.boundaryMesh();
674 
675  const polyPatch& patch1 = patches[couple.origPatchNames[0]];
676  const polyPatch& patch2 = patches[couple.origPatchNames[1]];
677 
678  boolList procHasPatch1(Pstream::nProcs(), false);
679  procHasPatch1[Pstream::myProcNo()] = !patch1.empty();
680  Pstream::gatherList(procHasPatch1);
681  Pstream::scatterList(procHasPatch1);
682 
683  boolList procHasPatch2(Pstream::nProcs(), false);
684  procHasPatch2[Pstream::myProcNo()] = !patch2.empty();
685  Pstream::gatherList(procHasPatch2);
686  Pstream::scatterList(procHasPatch2);
687 
688  // Multiple cyclic interfaces must be ordered in a specific way for
689  // processor communication to function correctly.
690  //
691  // A communication that is sent from the cyclic owner is received
692  // on the cyclic neighbour and vice versa. Therefore, in a coupled
693  // pair of processors if one sends the owner first the other must
694  // receive the neighbour first.
695  //
696  // We ensure the above by ordering the patches so that for the
697  // lower indexed processor the owner interface comes first, and for
698  // the higher indexed processor the neighbour comes first.
699 
700  auto appendProcPatches = [&](const bool owner, const bool first)
701  {
702  const boolList& procHasPatchA =
703  owner ? procHasPatch1 : procHasPatch2;
704  const boolList& procHasPatchB =
705  owner ? procHasPatch2 : procHasPatch1;
706 
707  if (procHasPatchA[Pstream::myProcNo()])
708  {
709  forAll(procHasPatchB, proci)
710  {
711  if
712  (
713  (
714  (first && proci > Pstream::myProcNo())
715  || (!first && proci < Pstream::myProcNo())
716  )
717  && procHasPatchB[proci]
718  )
719  {
720  newPatches[regioni].append
721  (
723  (
724  0,
725  mesh.nFaces(),
726  newPatches.size(),
727  patches,
729  proci,
730  couple.ncPatchNames[!owner],
731  couple.origPatchNames[!owner]
732  )
733  );
734  newPatchIsCouple[regioni].append(true);
735  newPatchFieldDicts[regioni].append
736  (
737  couple.ncPatchFieldDicts[!owner]
738  );
739  }
740  }
741  }
742  };
743 
744  appendProcPatches(true, true);
745  appendProcPatches(false, true);
746  appendProcPatches(false, false);
747  appendProcPatches(true, false);
748  }
749  }
750 
751  // Re-patch the meshes. Create constraint or calculated patch fields at
752  // this stage; don't apply the patch field dictionaries. The patch fields
753  // will be specified properly later after all patches have been added and
754  // the meshes have been stitched. That way the geometry is fully available
755  // for construction of the patch fields.
756  forAll(regionMeshes, regioni)
757  {
758  forAll(newPatches[regioni], newPatchi)
759  {
761  (
762  regionMeshes[regioni],
763  *newPatches[regioni][newPatchi],
764  dictionary(),
766  false
767  );
768  }
769  }
770 
771  // Connect the meshes
772  {
773  MultiRegionUList<fvMesh> multiRegionMeshes(regionMeshes, false);
774  forAll(regionMeshes, regioni)
775  {
776  RegionRef<fvMesh> mesh = multiRegionMeshes[regioni];
777  fvMeshStitchers::stationary(mesh).connect(false, false, false);
778  }
779  }
780 
781  // Set the fields on the new patches. This allows constraint types (e.g.,
782  // nonConformalCyclic) to be set by default, but requires a dictionary
783  // setting for non-constrained types (e.g., nonConformalMappedWall). It
784  // also allows constraint types to be overridden (e.g., with jumpCyclic) if
785  // there is a dictionary present.
786  forAll(regionMeshes, regioni)
787  {
788  forAll(newPatches[regioni], newPatchi)
789  {
790  if (newPatchIsCouple[regioni][newPatchi])
791  {
792  fvMeshTools::setPatchFields
793  (
794  regionMeshes[regioni],
795  newPatchi,
796  newPatchFieldDicts[regioni][newPatchi]
797  );
798  }
799  }
800  }
801 
802  // Communicate values across non-conformal processor cyclics so that they
803  // contain valid values that can be written to disk
804  if (Pstream::parRun())
805  {
806  forAll(regionMeshes, regioni)
807  {
808  const fvMesh& mesh = regionMeshes[regioni];
809 
810  #define EVALUATE_NON_CONFORMAL_PROCESSOR_CYCLICS(Type, nullArg) \
811  evaluateNonConformalProcessorCyclics<Type>(mesh);
812  FOR_ALL_FIELD_TYPES(EVALUATE_NON_CONFORMAL_PROCESSOR_CYCLICS)
813  #undef EVALUATE_NON_CONFORMAL_PROCESSOR_CYCLICS
814  }
815  }
816 
817  // Set the precision of the points data to 10
819 
820  // Set the instance so that the mesh writes
821  const word newInstance = overwrite ? oldInstance : runTime.name();
822  forAll(regionMeshes, regioni)
823  {
824  regionMeshes[regioni].setInstance(newInstance);
825  regionMeshes[regioni].setPolyFacesBfInstance(newInstance);
826  }
827 
828  // Write resulting mesh
829  Info<< nl << "Writing mesh to " << runTime.name() << nl << endl;
830  forAll(regionMeshes, regioni)
831  {
832  regionMeshes[regioni].write();
833  }
834 
835  Info<< "End\n" << endl;
836 
837  return 0;
838 }
839 
840 
841 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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.
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:458
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
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:243
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
Cyclic plane transformation.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:894
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:921
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
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, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:34
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
UPtrList< GeoField > fields(const bool strict=false) const
Return the list of fields of type GeoField.
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:267
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
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
#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:228
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:241
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:257
const word & regionName(const solver &region)
Definition: solver.H:209
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
objects
dictionary dict
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.