redistributePar.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-2023 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  redistributePar
26 
27 Description
28  Redistributes existing decomposed mesh and fields according to the current
29  settings in the decomposeParDict file.
30 
31  Must be run on maximum number of source and destination processors.
32  Balances mesh and writes new mesh to new time directory.
33 
34  Can also work like decomposePar:
35  \verbatim
36  # Create empty processor directories (have to exist for argList)
37  mkdir processor0
38  ..
39  mkdir processorN
40 
41  # Copy undecomposed polyMesh
42  cp -R constant processor0
43 
44  # Distribute
45  mpirun -np ddd redistributePar -parallel
46  \endverbatim
47 \*---------------------------------------------------------------------------*/
48 
49 #include "argList.H"
50 #include "timeSelector.H"
51 #include "decompositionMethod.H"
52 #include "PstreamReduceOps.H"
53 #include "fvMeshDistribute.H"
54 #include "polyDistributionMap.H"
55 #include "IOobjectList.H"
56 #include "globalIndex.H"
57 #include "loadOrCreateMesh.H"
59 
60 using namespace Foam;
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 void printMeshData(const polyMesh& mesh)
65 {
66  // Collect all data on master
67 
68  globalIndex globalCells(mesh.nCells());
69  labelListList patchNeiProcNo(Pstream::nProcs());
70  labelListList patchSize(Pstream::nProcs());
71  const labelList& pPatches = mesh.globalData().processorPatches();
72  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
73  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
74  forAll(pPatches, i)
75  {
76  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
77  (
78  mesh.boundaryMesh()[pPatches[i]]
79  );
80  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
81  patchSize[Pstream::myProcNo()][i] = ppp.size();
82  }
83  Pstream::gatherList(patchNeiProcNo);
84  Pstream::gatherList(patchSize);
85 
86 
87  // Print stats
88 
89  globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
90 
91  label maxProcCells = 0;
92  label totProcFaces = 0;
93  label maxProcPatches = 0;
94  label totProcPatches = 0;
95  label maxProcFaces = 0;
96 
97  for (label proci = 0; proci < Pstream::nProcs(); proci++)
98  {
99  Info<< endl
100  << "Processor " << proci << nl
101  << " Number of cells = " << globalCells.localSize(proci)
102  << endl;
103 
104  label nProcFaces = 0;
105 
106  const labelList& nei = patchNeiProcNo[proci];
107 
108  forAll(patchNeiProcNo[proci], i)
109  {
110  Info<< " Number of faces shared with processor "
111  << patchNeiProcNo[proci][i] << " = " << patchSize[proci][i]
112  << endl;
113 
114  nProcFaces += patchSize[proci][i];
115  }
116 
117  Info<< " Number of processor patches = " << nei.size() << nl
118  << " Number of processor faces = " << nProcFaces << nl
119  << " Number of boundary faces = "
120  << globalBoundaryFaces.localSize(proci) << endl;
121 
122  maxProcCells = max(maxProcCells, globalCells.localSize(proci));
123  totProcFaces += nProcFaces;
124  totProcPatches += nei.size();
125  maxProcPatches = max(maxProcPatches, nei.size());
126  maxProcFaces = max(maxProcFaces, nProcFaces);
127  }
128 
129  // Stats
130 
131  scalar avgProcCells = scalar(globalCells.size())/Pstream::nProcs();
132  scalar avgProcPatches = scalar(totProcPatches)/Pstream::nProcs();
133  scalar avgProcFaces = scalar(totProcFaces)/Pstream::nProcs();
134 
135  // In case of all faces on one processor. Just to avoid division by 0.
136  if (totProcPatches == 0)
137  {
138  avgProcPatches = 1;
139  }
140  if (totProcFaces == 0)
141  {
142  avgProcFaces = 1;
143  }
144 
145  Info<< nl
146  << "Number of processor faces = " << totProcFaces/2 << nl
147  << "Max number of cells = " << maxProcCells
148  << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
149  << "% above average " << avgProcCells << ")" << nl
150  << "Max number of processor patches = " << maxProcPatches
151  << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
152  << "% above average " << avgProcPatches << ")" << nl
153  << "Max number of faces between processors = " << maxProcFaces
154  << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
155  << "% above average " << avgProcFaces << ")" << nl
156  << endl;
157 }
158 
159 
160 // Debugging: write volScalarField with decomposition for post processing.
161 void writeDecomposition
162 (
163  const word& name,
164  const fvMesh& mesh,
165  const labelList& decomp
166 )
167 {
168  Info<< "Writing wanted cell distribution to volScalarField " << name
169  << " for postprocessing purposes." << nl << endl;
170 
171  volScalarField procCells
172  (
173  IOobject
174  (
175  name,
176  mesh.time().name(),
177  mesh,
180  false // do not register
181  ),
182  mesh,
184  extrapolatedCalculatedFvPatchScalarField::typeName
185  );
186 
187  forAll(procCells, cI)
188  {
189  procCells[cI] = decomp[cI];
190  }
191  procCells.correctBoundaryConditions();
192 
193  procCells.write();
194 }
195 
196 
197 template<class GeoField>
198 void readFields
199 (
200  const boolList& haveMesh,
201  const typename GeoField::Mesh& mesh,
202  const autoPtr<fvMeshSubset>& subsetterPtr,
203  IOobjectList& allObjects,
205 )
206 {
207  // Get my objects of type
208  IOobjectList objects(allObjects.lookupClass(GeoField::typeName));
209 
210  // Check that we all have all objects
211  wordList objectNames = objects.toc();
212 
213  // Get master names
214  wordList masterNames(objectNames);
215  Pstream::scatter(masterNames);
216 
217  if (haveMesh[Pstream::myProcNo()] && objectNames != masterNames)
218  {
220  << "differing fields of type " << GeoField::typeName
221  << " on processors." << endl
222  << "Master has:" << masterNames << endl
223  << Pstream::myProcNo() << " has:" << objectNames
224  << abort(FatalError);
225  }
226 
227  fields.setSize(masterNames.size());
228 
229  // Have master send all fields to processors that don't have a mesh
230  if (Pstream::master())
231  {
232  forAll(masterNames, i)
233  {
234  const word& name = masterNames[i];
235  IOobject& io = *objects[name];
237 
238  // Load field
239  fields.set(i, new GeoField(io, mesh));
240 
241  // Create zero sized field and send
242  if (subsetterPtr.valid())
243  {
244  tmp<GeoField> tsubfld = subsetterPtr().interpolate(fields[i]);
245 
246  // Send to all processors that don't have a mesh
247  for (label proci = 1; proci < Pstream::nProcs(); proci++)
248  {
249  if (!haveMesh[proci])
250  {
252  toProc<< tsubfld();
253  }
254  }
255  }
256  }
257  }
258  else if (!haveMesh[Pstream::myProcNo()])
259  {
260  // Don't have mesh (nor fields). Receive empty field from master.
261 
262  forAll(masterNames, i)
263  {
264  const word& name = masterNames[i];
265 
266  // Receive field
267  IPstream fromMaster
268  (
271  );
272  dictionary fieldDict(fromMaster);
273 
274  fields.set
275  (
276  i,
277  new GeoField
278  (
279  IOobject
280  (
281  name,
282  mesh.thisDb().time().name(),
283  mesh.thisDb(),
286  ),
287  mesh,
288  fieldDict
289  )
290  );
291 
293  // fields[i].write();
294  }
295  }
296  else
297  {
298  // Have mesh so just try to load
299  forAll(masterNames, i)
300  {
301  const word& name = masterNames[i];
302  IOobject& io = *objects[name];
304 
305  // Load field
306  fields.set(i, new GeoField(io, mesh));
307  }
308  }
309 }
310 
311 
312 int main(int argc, char *argv[])
313 {
314  #include "addRegionOption.H"
315  #include "addOverwriteOption.H"
316 
317  // Include explicit constant options, have zero from time range
319 
320  #include "setRootCase.H"
321 
322  if (env("FOAM_SIGFPE"))
323  {
325  << "Detected floating point exception trapping (FOAM_SIGFPE)."
326  << " This might give" << nl
327  << " problems when mapping fields. Switch it off in case"
328  << " of problems." << endl;
329  }
330 
331 
332  // Create processor directory if non-existing
333  if (!Pstream::master() && !isDir(args.path()))
334  {
335  Pout<< "Creating case directory " << args.path() << endl;
336  mkDir(args.path());
337  }
338 
339 
340  // Make sure we do not use the master-only reading.
343  // Allow override of time
345  runTime.setTime(times[0], 0);
346 
348  fileName meshSubDir;
349 
350  if (args.optionReadIfPresent("region", regionName))
351  {
352  meshSubDir = regionName/polyMesh::meshSubDir;
353  }
354  else
355  {
356  meshSubDir = polyMesh::meshSubDir;
357  }
358  Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
359 
360  const bool overwrite = args.optionFound("overwrite");
361 
362 
363  // Get time instance directory. Since not all processors have meshes
364  // just use the master one everywhere.
365 
366  fileName masterInstDir;
367  if (Pstream::master())
368  {
369  masterInstDir = runTime.findInstance(meshSubDir, "points");
370  }
371  Pstream::scatter(masterInstDir);
372 
373  // Check who has a mesh
374  const fileName meshPath = runTime.path()/masterInstDir/meshSubDir;
375 
376  Info<< "Found points in " << meshPath << nl << endl;
377 
378 
379  boolList haveMesh(Pstream::nProcs(), false);
380  haveMesh[Pstream::myProcNo()] = isDir(meshPath);
381  Pstream::gatherList(haveMesh);
382  Pstream::scatterList(haveMesh);
383  Info<< "Per processor mesh availability : " << haveMesh << endl;
384  const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
385 
387  (
388  IOobject
389  (
390  regionName,
391  masterInstDir,
392  runTime,
394  )
395  );
396 
397  fvMesh& mesh = meshPtr();
398 
399  // Print some statistics
400  Info<< "Before distribution:" << endl;
401  printMeshData(mesh);
402 
403 
404  labelList finalDecomp;
405 
406  // Create decompositionMethod and new decomposition
407  {
408  autoPtr<decompositionMethod> distributor
409  (
411  (
413  )
414  );
415 
416  finalDecomp = distributor().decompose(mesh, mesh.cellCentres());
417  }
418 
419  // Dump decomposition to volScalarField
420  if (!overwrite)
421  {
422  writeDecomposition("decomposition", mesh, finalDecomp);
423  }
424 
425 
426  // Create 0 sized mesh to do all the generation of zero sized
427  // fields on processors that have zero sized meshes. Note that this is
428  // only necessary on master but since polyMesh construction with
429  // Pstream::parRun does parallel comms we have to do it on all
430  // processors
431  autoPtr<fvMeshSubset> subsetterPtr;
432 
433  if (!allHaveMesh)
434  {
435  // Find last non-processor patch.
436  const polyBoundaryMesh& patches = mesh.boundaryMesh();
437 
438  label nonProci = -1;
439 
441  {
442  if (isA<processorPolyPatch>(patches[patchi]))
443  {
444  break;
445  }
446  nonProci++;
447  }
448 
449  if (nonProci == -1)
450  {
452  << "Cannot find non-processor patch on processor "
453  << Pstream::myProcNo() << endl
454  << " Current patches:" << patches.names() << abort(FatalError);
455  }
456 
457  // Subset 0 cells, no parallel comms. This is used to create zero-sized
458  // fields.
459  subsetterPtr.reset(new fvMeshSubset(mesh));
460  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProci, false);
461  }
462 
463 
464  // Get original objects (before incrementing time!)
465  IOobjectList objects(mesh, runTime.name());
466  // We don't want to map the decomposition (mapping already tested when
467  // mapping the cell centre field)
468  IOobjectList::iterator iter = objects.find("decomposition");
469  if (iter != objects.end())
470  {
471  objects.erase(iter);
472  }
473 
474 
475  // volFields
476 
477  PtrList<volScalarField> volScalarFields;
478  readFields
479  (
480  haveMesh,
481  mesh,
482  subsetterPtr,
483  objects,
484  volScalarFields
485  );
486 
487  PtrList<volVectorField> volVectorFields;
488  readFields
489  (
490  haveMesh,
491  mesh,
492  subsetterPtr,
493  objects,
494  volVectorFields
495  );
496 
497  PtrList<volSphericalTensorField> volSphereTensorFields;
498  readFields
499  (
500  haveMesh,
501  mesh,
502  subsetterPtr,
503  objects,
504  volSphereTensorFields
505  );
506 
507  PtrList<volSymmTensorField> volSymmTensorFields;
508  readFields
509  (
510  haveMesh,
511  mesh,
512  subsetterPtr,
513  objects,
514  volSymmTensorFields
515  );
516 
517  PtrList<volTensorField> volTensorFields;
518  readFields
519  (
520  haveMesh,
521  mesh,
522  subsetterPtr,
523  objects,
524  volTensorFields
525  );
526 
527 
528  // surfaceFields
529 
530  PtrList<surfaceScalarField> surfScalarFields;
531  readFields
532  (
533  haveMesh,
534  mesh,
535  subsetterPtr,
536  objects,
537  surfScalarFields
538  );
539 
540  PtrList<surfaceVectorField> surfVectorFields;
541  readFields
542  (
543  haveMesh,
544  mesh,
545  subsetterPtr,
546  objects,
547  surfVectorFields
548  );
549 
550  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
551  readFields
552  (
553  haveMesh,
554  mesh,
555  subsetterPtr,
556  objects,
557  surfSphereTensorFields
558  );
559 
560  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
561  readFields
562  (
563  haveMesh,
564  mesh,
565  subsetterPtr,
566  objects,
567  surfSymmTensorFields
568  );
569 
570  PtrList<surfaceTensorField> surfTensorFields;
571  readFields
572  (
573  haveMesh,
574  mesh,
575  subsetterPtr,
576  objects,
577  surfTensorFields
578  );
579 
580 
581  // pointFields
582 
583  pointMesh& pMesh =
584  const_cast<pointMesh&>(pointMesh::New(mesh));
585  PtrList<pointScalarField> pointScalarFields;
586  readFields
587  (
588  haveMesh,
589  pMesh,
590  subsetterPtr,
591  objects,
592  pointScalarFields
593  );
594 
595  PtrList<pointVectorField> pointVectorFields;
596  readFields
597  (
598  haveMesh,
599  pMesh,
600  subsetterPtr,
601  objects,
602  pointVectorFields
603  );
604 
605  PtrList<pointSphericalTensorField> pointSphereTensorFields;
606  readFields
607  (
608  haveMesh,
609  pMesh,
610  subsetterPtr,
611  objects,
612  pointSphereTensorFields
613  );
614 
615  PtrList<pointSymmTensorField> pointSymmTensorFields;
616  readFields
617  (
618  haveMesh,
619  pMesh,
620  subsetterPtr,
621  objects,
622  pointSymmTensorFields
623  );
624 
625  PtrList<pointTensorField> pointTensorFields;
626  readFields
627  (
628  haveMesh,
629  pMesh,
630  subsetterPtr,
631  objects,
632  pointTensorFields
633  );
634 
635  // Debugging: Create additional volField that will be mapped.
636  // Used to test correctness of mapping
637  // volVectorField mapCc("mapCc", 1*mesh.C());
638 
639  // Mesh distribution engine
640  fvMeshDistribute distributor(mesh);
641 
642  // Pout<< "Wanted distribution:"
643  // << distributor.countCells(finalDecomp) << nl << endl;
644 
645  // Do actual sending/receiving of mesh
646  autoPtr<polyDistributionMap> map = distributor.distribute(finalDecomp);
647 
649  // map().distributeFaceData(faceCc);
650 
651 
652  // Print some statistics
653  Info<< "After distribution:" << endl;
654  printMeshData(mesh);
655 
656 
657  if (!overwrite)
658  {
659  runTime++;
660  }
661  else
662  {
663  mesh.setInstance(masterInstDir);
664  }
665  Info<< "Writing redistributed mesh to " << runTime.name() << nl << endl;
666  mesh.write();
667 
668 
669  // Print nice message
670  // ~~~~~~~~~~~~~~~~~~
671 
672  // Determine which processors remain so we can print nice final message.
673  labelList nFaces(Pstream::nProcs());
674  nFaces[Pstream::myProcNo()] = mesh.nFaces();
675  Pstream::gatherList(nFaces);
676  Pstream::scatterList(nFaces);
677 
678  Info<< nl
679  << "You can pick up the redecomposed mesh from the polyMesh directory"
680  << " in " << runTime.name() << "." << nl
681  << "If you redecomposed the mesh to less processors you can delete"
682  << nl
683  << "the processor directories with 0 sized meshes in them." << nl
684  << "Below is a sample set of commands to do this."
685  << " Take care when issuing these" << nl
686  << "commands." << nl << endl;
687 
688  forAll(nFaces, proci)
689  {
690  fileName procDir = "processor" + name(proci);
691 
692  if (nFaces[proci] == 0)
693  {
694  Info<< " rm -r " << procDir.c_str() << nl;
695  }
696  else
697  {
698  fileName timeDir = procDir/runTime.name()/meshSubDir;
699  fileName constDir = procDir/runTime.constant()/meshSubDir;
700 
701  Info<< " rm -r " << constDir.c_str() << nl
702  << " mv " << timeDir.c_str()
703  << ' ' << constDir.c_str() << nl;
704  }
705  }
706  Info<< endl;
707 
708 
709  Info<< "End\n" << endl;
710 
711  return 0;
712 }
713 
714 
715 // ************************************************************************* //
Inter-processor communication reduction functions.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
HashTable< IOobject *, word, string::hash >::iterator iterator
Definition: HashPtrTable.H:82
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:190
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
writeOption writeOpt() const
Definition: IOobject.H:370
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:226
Input inter-processor communications stream.
Definition: IPstream.H:54
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Output inter-processor communications stream.
Definition: OPstream.H:54
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
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
static int masterNo()
Process index of the master.
Definition: UPstream.H:417
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
fileName path() const
Return the path to the caseName.
Definition: argListI.H:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
static autoPtr< decompositionMethod > NewDistributor(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
static dictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const word & name() const
Return const reference to name.
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:265
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:402
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1736
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
const labelList & processorPatches() const
Return list of processor patch labels.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1562
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:271
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
label nCells() const
const vectorField & cellCentres() const
label nInternalFaces() const
label nFaces() const
Neighbour processor patch.
int neighbProcNo() const
Return neighbour processor number.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList selectIfPresent(Time &runTime, const argList &args)
If any time option provided return the set of times (as select0)
Definition: timeSelector.C:283
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::word regionName
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static fvMesh * meshPtr
Definition: globalFoam.H:52
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:230
Load or create (0 size) a mesh. Used in distributing meshes to a larger number of processors.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
bool env(const word &)
Return true if environment variable of given name is defined.
Definition: POSIX.C:91
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
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
messageStream Info
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
static const char nl
Definition: Ostream.H:260
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
objects
Foam::argList args(argc, argv)