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 "volFields.H"
54 #include "fvMeshDistribute.H"
55 #include "polyDistributionMap.H"
56 #include "IOobjectList.H"
57 #include "globalIndex.H"
58 #include "loadOrCreateMesh.H"
60 
61 using namespace Foam;
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 void printMeshData(const polyMesh& mesh)
66 {
67  // Collect all data on master
68 
69  globalIndex globalCells(mesh.nCells());
70  labelListList patchNeiProcNo(Pstream::nProcs());
71  labelListList patchSize(Pstream::nProcs());
72  const labelList& pPatches = mesh.globalData().processorPatches();
73  patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
74  patchSize[Pstream::myProcNo()].setSize(pPatches.size());
75  forAll(pPatches, i)
76  {
77  const processorPolyPatch& ppp = refCast<const processorPolyPatch>
78  (
79  mesh.boundaryMesh()[pPatches[i]]
80  );
81  patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
82  patchSize[Pstream::myProcNo()][i] = ppp.size();
83  }
84  Pstream::gatherList(patchNeiProcNo);
85  Pstream::gatherList(patchSize);
86 
87 
88  // Print stats
89 
90  globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
91 
92  label maxProcCells = 0;
93  label totProcFaces = 0;
94  label maxProcPatches = 0;
95  label totProcPatches = 0;
96  label maxProcFaces = 0;
97 
98  for (label proci = 0; proci < Pstream::nProcs(); proci++)
99  {
100  Info<< endl
101  << "Processor " << proci << nl
102  << " Number of cells = " << globalCells.localSize(proci)
103  << endl;
104 
105  label nProcFaces = 0;
106 
107  const labelList& nei = patchNeiProcNo[proci];
108 
109  forAll(patchNeiProcNo[proci], i)
110  {
111  Info<< " Number of faces shared with processor "
112  << patchNeiProcNo[proci][i] << " = " << patchSize[proci][i]
113  << endl;
114 
115  nProcFaces += patchSize[proci][i];
116  }
117 
118  Info<< " Number of processor patches = " << nei.size() << nl
119  << " Number of processor faces = " << nProcFaces << nl
120  << " Number of boundary faces = "
121  << globalBoundaryFaces.localSize(proci) << endl;
122 
123  maxProcCells = max(maxProcCells, globalCells.localSize(proci));
124  totProcFaces += nProcFaces;
125  totProcPatches += nei.size();
126  maxProcPatches = max(maxProcPatches, nei.size());
127  maxProcFaces = max(maxProcFaces, nProcFaces);
128  }
129 
130  // Stats
131 
132  scalar avgProcCells = scalar(globalCells.size())/Pstream::nProcs();
133  scalar avgProcPatches = scalar(totProcPatches)/Pstream::nProcs();
134  scalar avgProcFaces = scalar(totProcFaces)/Pstream::nProcs();
135 
136  // In case of all faces on one processor. Just to avoid division by 0.
137  if (totProcPatches == 0)
138  {
139  avgProcPatches = 1;
140  }
141  if (totProcFaces == 0)
142  {
143  avgProcFaces = 1;
144  }
145 
146  Info<< nl
147  << "Number of processor faces = " << totProcFaces/2 << nl
148  << "Max number of cells = " << maxProcCells
149  << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
150  << "% above average " << avgProcCells << ")" << nl
151  << "Max number of processor patches = " << maxProcPatches
152  << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
153  << "% above average " << avgProcPatches << ")" << nl
154  << "Max number of faces between processors = " << maxProcFaces
155  << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
156  << "% above average " << avgProcFaces << ")" << nl
157  << endl;
158 }
159 
160 
161 // Debugging: write volScalarField with decomposition for post processing.
162 void writeDecomposition
163 (
164  const word& name,
165  const fvMesh& mesh,
166  const labelList& decomp
167 )
168 {
169  Info<< "Writing wanted cell distribution to volScalarField " << name
170  << " for postprocessing purposes." << nl << endl;
171 
172  volScalarField procCells
173  (
174  IOobject
175  (
176  name,
177  mesh.time().name(),
178  mesh,
181  false // do not register
182  ),
183  mesh,
185  extrapolatedCalculatedFvPatchScalarField::typeName
186  );
187 
188  forAll(procCells, cI)
189  {
190  procCells[cI] = decomp[cI];
191  }
192  procCells.correctBoundaryConditions();
193 
194  procCells.write();
195 }
196 
197 
198 template<class GeoField>
199 void readFields
200 (
201  const boolList& haveMesh,
202  const typename GeoField::Mesh& mesh,
203  const autoPtr<fvMeshSubset>& subsetterPtr,
204  IOobjectList& allObjects,
206 )
207 {
208  // Get my objects of type
209  IOobjectList objects(allObjects.lookupClass(GeoField::typeName));
210 
211  // Check that we all have all objects
212  wordList objectNames = objects.toc();
213 
214  // Get master names
215  wordList masterNames(objectNames);
216  Pstream::scatter(masterNames);
217 
218  if (haveMesh[Pstream::myProcNo()] && objectNames != masterNames)
219  {
221  << "differing fields of type " << GeoField::typeName
222  << " on processors." << endl
223  << "Master has:" << masterNames << endl
224  << Pstream::myProcNo() << " has:" << objectNames
225  << abort(FatalError);
226  }
227 
228  fields.setSize(masterNames.size());
229 
230  // Have master send all fields to processors that don't have a mesh
231  if (Pstream::master())
232  {
233  forAll(masterNames, i)
234  {
235  const word& name = masterNames[i];
236  IOobject& io = *objects[name];
238 
239  // Load field
240  fields.set(i, new GeoField(io, mesh));
241 
242  // Create zero sized field and send
243  if (subsetterPtr.valid())
244  {
245  tmp<GeoField> tsubfld = subsetterPtr().interpolate(fields[i]);
246 
247  // Send to all processors that don't have a mesh
248  for (label proci = 1; proci < Pstream::nProcs(); proci++)
249  {
250  if (!haveMesh[proci])
251  {
253  toProc<< tsubfld();
254  }
255  }
256  }
257  }
258  }
259  else if (!haveMesh[Pstream::myProcNo()])
260  {
261  // Don't have mesh (nor fields). Receive empty field from master.
262 
263  forAll(masterNames, i)
264  {
265  const word& name = masterNames[i];
266 
267  // Receive field
268  IPstream fromMaster
269  (
272  );
273  dictionary fieldDict(fromMaster);
274 
275  fields.set
276  (
277  i,
278  new GeoField
279  (
280  IOobject
281  (
282  name,
283  mesh.thisDb().time().name(),
284  mesh.thisDb(),
287  ),
288  mesh,
289  fieldDict
290  )
291  );
292 
294  // fields[i].write();
295  }
296  }
297  else
298  {
299  // Have mesh so just try to load
300  forAll(masterNames, i)
301  {
302  const word& name = masterNames[i];
303  IOobject& io = *objects[name];
305 
306  // Load field
307  fields.set(i, new GeoField(io, mesh));
308  }
309  }
310 }
311 
312 
313 int main(int argc, char *argv[])
314 {
315  #include "addRegionOption.H"
316  #include "addOverwriteOption.H"
317 
318  // Include explicit constant options, have zero from time range
320 
321  #include "setRootCase.H"
322 
323  if (env("FOAM_SIGFPE"))
324  {
326  << "Detected floating point exception trapping (FOAM_SIGFPE)."
327  << " This might give" << nl
328  << " problems when mapping fields. Switch it off in case"
329  << " of problems." << endl;
330  }
331 
332 
333  // Create processor directory if non-existing
334  if (!Pstream::master() && !isDir(args.path()))
335  {
336  Pout<< "Creating case directory " << args.path() << endl;
337  mkDir(args.path());
338  }
339 
340 
341  // Make sure we do not use the master-only reading.
344  // Allow override of time
346  runTime.setTime(times[0], 0);
347 
349  fileName meshSubDir;
350 
351  if (args.optionReadIfPresent("region", regionName))
352  {
353  meshSubDir = regionName/polyMesh::meshSubDir;
354  }
355  else
356  {
357  meshSubDir = polyMesh::meshSubDir;
358  }
359  Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
360 
361  const bool overwrite = args.optionFound("overwrite");
362 
363 
364  // Get time instance directory. Since not all processors have meshes
365  // just use the master one everywhere.
366 
367  fileName masterInstDir;
368  if (Pstream::master())
369  {
370  masterInstDir = runTime.findInstance(meshSubDir, "points");
371  }
372  Pstream::scatter(masterInstDir);
373 
374  // Check who has a mesh
375  const fileName meshPath = runTime.path()/masterInstDir/meshSubDir;
376 
377  Info<< "Found points in " << meshPath << nl << endl;
378 
379 
380  boolList haveMesh(Pstream::nProcs(), false);
381  haveMesh[Pstream::myProcNo()] = isDir(meshPath);
382  Pstream::gatherList(haveMesh);
383  Pstream::scatterList(haveMesh);
384  Info<< "Per processor mesh availability : " << haveMesh << endl;
385  const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
386 
388  (
389  IOobject
390  (
391  regionName,
392  masterInstDir,
393  runTime,
395  )
396  );
397 
398  fvMesh& mesh = meshPtr();
399 
400  // Print some statistics
401  Info<< "Before distribution:" << endl;
402  printMeshData(mesh);
403 
404 
405  labelList finalDecomp;
406 
407  // Create decompositionMethod and new decomposition
408  {
409  autoPtr<decompositionMethod> distributor
410  (
412  (
414  )
415  );
416 
417  finalDecomp = distributor().decompose(mesh, mesh.cellCentres());
418  }
419 
420  // Dump decomposition to volScalarField
421  if (!overwrite)
422  {
423  writeDecomposition("decomposition", mesh, finalDecomp);
424  }
425 
426 
427  // Create 0 sized mesh to do all the generation of zero sized
428  // fields on processors that have zero sized meshes. Note that this is
429  // only necessary on master but since polyMesh construction with
430  // Pstream::parRun does parallel comms we have to do it on all
431  // processors
432  autoPtr<fvMeshSubset> subsetterPtr;
433 
434  if (!allHaveMesh)
435  {
436  // Find last non-processor patch.
437  const polyBoundaryMesh& patches = mesh.boundaryMesh();
438 
439  label nonProci = -1;
440 
442  {
443  if (isA<processorPolyPatch>(patches[patchi]))
444  {
445  break;
446  }
447  nonProci++;
448  }
449 
450  if (nonProci == -1)
451  {
453  << "Cannot find non-processor patch on processor "
454  << Pstream::myProcNo() << endl
455  << " Current patches:" << patches.names() << abort(FatalError);
456  }
457 
458  // Subset 0 cells, no parallel comms. This is used to create zero-sized
459  // fields.
460  subsetterPtr.reset(new fvMeshSubset(mesh));
461  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProci, false);
462  }
463 
464 
465  // Get original objects (before incrementing time!)
466  IOobjectList objects(mesh, runTime.name());
467  // We don't want to map the decomposition (mapping already tested when
468  // mapping the cell centre field)
469  IOobjectList::iterator iter = objects.find("decomposition");
470  if (iter != objects.end())
471  {
472  objects.erase(iter);
473  }
474 
475 
476  // volFields
477 
478  PtrList<volScalarField> volScalarFields;
479  readFields
480  (
481  haveMesh,
482  mesh,
483  subsetterPtr,
484  objects,
485  volScalarFields
486  );
487 
488  PtrList<volVectorField> volVectorFields;
489  readFields
490  (
491  haveMesh,
492  mesh,
493  subsetterPtr,
494  objects,
495  volVectorFields
496  );
497 
498  PtrList<volSphericalTensorField> volSphereTensorFields;
499  readFields
500  (
501  haveMesh,
502  mesh,
503  subsetterPtr,
504  objects,
505  volSphereTensorFields
506  );
507 
508  PtrList<volSymmTensorField> volSymmTensorFields;
509  readFields
510  (
511  haveMesh,
512  mesh,
513  subsetterPtr,
514  objects,
515  volSymmTensorFields
516  );
517 
518  PtrList<volTensorField> volTensorFields;
519  readFields
520  (
521  haveMesh,
522  mesh,
523  subsetterPtr,
524  objects,
525  volTensorFields
526  );
527 
528 
529  // surfaceFields
530 
531  PtrList<surfaceScalarField> surfScalarFields;
532  readFields
533  (
534  haveMesh,
535  mesh,
536  subsetterPtr,
537  objects,
538  surfScalarFields
539  );
540 
541  PtrList<surfaceVectorField> surfVectorFields;
542  readFields
543  (
544  haveMesh,
545  mesh,
546  subsetterPtr,
547  objects,
548  surfVectorFields
549  );
550 
551  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
552  readFields
553  (
554  haveMesh,
555  mesh,
556  subsetterPtr,
557  objects,
558  surfSphereTensorFields
559  );
560 
561  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
562  readFields
563  (
564  haveMesh,
565  mesh,
566  subsetterPtr,
567  objects,
568  surfSymmTensorFields
569  );
570 
571  PtrList<surfaceTensorField> surfTensorFields;
572  readFields
573  (
574  haveMesh,
575  mesh,
576  subsetterPtr,
577  objects,
578  surfTensorFields
579  );
580 
581 
582  // pointFields
583 
584  pointMesh& pMesh =
585  const_cast<pointMesh&>(pointMesh::New(mesh));
586  PtrList<pointScalarField> pointScalarFields;
587  readFields
588  (
589  haveMesh,
590  pMesh,
591  subsetterPtr,
592  objects,
593  pointScalarFields
594  );
595 
596  PtrList<pointVectorField> pointVectorFields;
597  readFields
598  (
599  haveMesh,
600  pMesh,
601  subsetterPtr,
602  objects,
603  pointVectorFields
604  );
605 
606  PtrList<pointSphericalTensorField> pointSphereTensorFields;
607  readFields
608  (
609  haveMesh,
610  pMesh,
611  subsetterPtr,
612  objects,
613  pointSphereTensorFields
614  );
615 
616  PtrList<pointSymmTensorField> pointSymmTensorFields;
617  readFields
618  (
619  haveMesh,
620  pMesh,
621  subsetterPtr,
622  objects,
623  pointSymmTensorFields
624  );
625 
626  PtrList<pointTensorField> pointTensorFields;
627  readFields
628  (
629  haveMesh,
630  pMesh,
631  subsetterPtr,
632  objects,
633  pointTensorFields
634  );
635 
636  // Debugging: Create additional volField that will be mapped.
637  // Used to test correctness of mapping
638  // volVectorField mapCc("mapCc", 1*mesh.C());
639 
640  // Mesh distribution engine
641  fvMeshDistribute distributor(mesh);
642 
643  // Pout<< "Wanted distribution:"
644  // << distributor.countCells(finalDecomp) << nl << endl;
645 
646  // Do actual sending/receiving of mesh
647  autoPtr<polyDistributionMap> map = distributor.distribute(finalDecomp);
648 
650  // map().distributeFaceData(faceCc);
651 
652 
653  // Print some statistics
654  Info<< "After distribution:" << endl;
655  printMeshData(mesh);
656 
657 
658  if (!overwrite)
659  {
660  runTime++;
661  }
662  else
663  {
664  mesh.setInstance(masterInstDir);
665  }
666  Info<< "Writing redistributed mesh to " << runTime.name() << nl << endl;
667  mesh.write();
668 
669 
670  // Print nice message
671  // ~~~~~~~~~~~~~~~~~~
672 
673  // Determine which processors remain so we can print nice final message.
674  labelList nFaces(Pstream::nProcs());
675  nFaces[Pstream::myProcNo()] = mesh.nFaces();
676  Pstream::gatherList(nFaces);
677  Pstream::scatterList(nFaces);
678 
679  Info<< nl
680  << "You can pick up the redecomposed mesh from the polyMesh directory"
681  << " in " << runTime.name() << "." << nl
682  << "If you redecomposed the mesh to less processors you can delete"
683  << nl
684  << "the processor directories with 0 sized meshes in them." << nl
685  << "Below is a sample set of commands to do this."
686  << " Take care when issuing these" << nl
687  << "commands." << nl << endl;
688 
689  forAll(nFaces, proci)
690  {
691  fileName procDir = "processor" + name(proci);
692 
693  if (nFaces[proci] == 0)
694  {
695  Info<< " rm -r " << procDir.c_str() << nl;
696  }
697  else
698  {
699  fileName timeDir = procDir/runTime.name()/meshSubDir;
700  fileName constDir = procDir/runTime.constant()/meshSubDir;
701 
702  Info<< " rm -r " << constDir.c_str() << nl
703  << " mv " << timeDir.c_str()
704  << ' ' << constDir.c_str() << nl;
705  }
706  }
707  Info<< endl;
708 
709 
710  Info<< "End\n" << endl;
711 
712  return 0;
713 }
714 
715 
716 // ************************************************************************* //
Inter-processor communication reduction functions.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
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
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
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 IOdictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
static autoPtr< decompositionMethod > NewDistributor(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
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:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1737
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:53
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1515
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:270
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
#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
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:228
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:257
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
const word & regionName(const solver &region)
Definition: solver.H:209
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:213
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:266
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
objects
Foam::argList args(argc, argv)