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-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  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 "addMeshOption.H"
316  #include "addRegionOption.H"
317  #include "addNoOverwriteOption.H"
318 
319  // Include explicit constant options, have zero from time range
321 
322  #include "setRootCase.H"
323  #include "setMeshPath.H"
324 
325  if (env("FOAM_SIGFPE"))
326  {
328  << "Detected floating point exception trapping (FOAM_SIGFPE)."
329  << " This might give" << nl
330  << " problems when mapping fields. Switch it off in case"
331  << " of problems." << endl;
332  }
333 
334 
335  // Create processor directory if non-existing
336  if (!Pstream::master() && !isDir(args.path()))
337  {
338  Pout<< "Creating case directory " << args.path() << endl;
339  mkDir(args.path());
340  }
341 
342 
343  // Make sure we do not use the master-only reading.
346  // Allow override of time
348  runTime.setTime(times[0], 0);
349 
351  fileName meshSubDir;
352 
353  if (args.optionReadIfPresent("region", regionName))
354  {
355  meshSubDir = regionName/polyMesh::meshSubDir;
356  }
357  else
358  {
359  meshSubDir = polyMesh::meshSubDir;
360  }
361  Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
362 
363  #include "setNoOverwrite.H"
364 
365 
366  // Get time instance directory. Since not all processors have meshes
367  // just use the master one everywhere.
368 
369  fileName masterInstDir;
370  if (Pstream::master())
371  {
372  masterInstDir = runTime.findInstance(meshSubDir, "points");
373  }
374  Pstream::scatter(masterInstDir);
375 
376  // Check who has a mesh
377  const fileName meshAbsolutePath = runTime.path()/masterInstDir/meshSubDir;
378 
379  Info<< "Found points in " << meshAbsolutePath << nl << endl;
380 
381 
382  boolList haveMesh(Pstream::nProcs(), false);
383  haveMesh[Pstream::myProcNo()] = isDir(meshAbsolutePath);
384  Pstream::gatherList(haveMesh);
385  Pstream::scatterList(haveMesh);
386  Info<< "Per processor mesh availability : " << haveMesh << endl;
387  const bool allHaveMesh = (findIndex(haveMesh, false) == -1);
388 
390  (
391  IOobject
392  (
393  regionName,
394  masterInstDir,
395  meshPath,
396  runTime,
398  )
399  );
400 
401  fvMesh& mesh = meshPtr();
402 
403  // Print some statistics
404  Info<< "Before distribution:" << endl;
405  printMeshData(mesh);
406 
407 
408  labelList finalDecomp;
409 
410  // Create decompositionMethod and new decomposition
411  {
412  autoPtr<decompositionMethod> distributor
413  (
415  (
417  )
418  );
419 
420  finalDecomp = distributor().decompose(mesh, mesh.cellCentres());
421  }
422 
423  // Dump decomposition to volScalarField
424  if (!overwrite)
425  {
426  writeDecomposition("decomposition", mesh, finalDecomp);
427  }
428 
429 
430  // Create 0 sized mesh to do all the generation of zero sized
431  // fields on processors that have zero sized meshes. Note that this is
432  // only necessary on master but since polyMesh construction with
433  // Pstream::parRun does parallel comms we have to do it on all
434  // processors
435  autoPtr<fvMeshSubset> subsetterPtr;
436 
437  if (!allHaveMesh)
438  {
439  // Find last non-processor patch.
441 
442  label nonProci = -1;
443 
445  {
446  if (isA<processorPolyPatch>(patches[patchi]))
447  {
448  break;
449  }
450  nonProci++;
451  }
452 
453  if (nonProci == -1)
454  {
456  << "Cannot find non-processor patch on processor "
457  << Pstream::myProcNo() << endl
458  << " Current patches:" << patches.names() << abort(FatalError);
459  }
460 
461  // Subset 0 cells, no parallel comms. This is used to create zero-sized
462  // fields.
463  subsetterPtr.reset(new fvMeshSubset(mesh));
464  subsetterPtr().setLargeCellSubset(labelHashSet(0), nonProci, false);
465  }
466 
467 
468  // Get original objects (before incrementing time!)
469  IOobjectList objects(mesh, runTime.name());
470  // We don't want to map the decomposition (mapping already tested when
471  // mapping the cell centre field)
472  IOobjectList::iterator iter = objects.find("decomposition");
473  if (iter != objects.end())
474  {
475  objects.erase(iter);
476  }
477 
478 
479  // volFields
480 
481  PtrList<volScalarField> volScalarFields;
482  readFields
483  (
484  haveMesh,
485  mesh,
486  subsetterPtr,
487  objects,
488  volScalarFields
489  );
490 
491  PtrList<volVectorField> volVectorFields;
492  readFields
493  (
494  haveMesh,
495  mesh,
496  subsetterPtr,
497  objects,
498  volVectorFields
499  );
500 
501  PtrList<volSphericalTensorField> volSphereTensorFields;
502  readFields
503  (
504  haveMesh,
505  mesh,
506  subsetterPtr,
507  objects,
508  volSphereTensorFields
509  );
510 
511  PtrList<volSymmTensorField> volSymmTensorFields;
512  readFields
513  (
514  haveMesh,
515  mesh,
516  subsetterPtr,
517  objects,
518  volSymmTensorFields
519  );
520 
521  PtrList<volTensorField> volTensorFields;
522  readFields
523  (
524  haveMesh,
525  mesh,
526  subsetterPtr,
527  objects,
528  volTensorFields
529  );
530 
531 
532  // surfaceFields
533 
534  PtrList<surfaceScalarField> surfScalarFields;
535  readFields
536  (
537  haveMesh,
538  mesh,
539  subsetterPtr,
540  objects,
541  surfScalarFields
542  );
543 
544  PtrList<surfaceVectorField> surfVectorFields;
545  readFields
546  (
547  haveMesh,
548  mesh,
549  subsetterPtr,
550  objects,
551  surfVectorFields
552  );
553 
554  PtrList<surfaceSphericalTensorField> surfSphereTensorFields;
555  readFields
556  (
557  haveMesh,
558  mesh,
559  subsetterPtr,
560  objects,
561  surfSphereTensorFields
562  );
563 
564  PtrList<surfaceSymmTensorField> surfSymmTensorFields;
565  readFields
566  (
567  haveMesh,
568  mesh,
569  subsetterPtr,
570  objects,
571  surfSymmTensorFields
572  );
573 
574  PtrList<surfaceTensorField> surfTensorFields;
575  readFields
576  (
577  haveMesh,
578  mesh,
579  subsetterPtr,
580  objects,
581  surfTensorFields
582  );
583 
584 
585  // pointFields
586 
587  pointMesh& pMesh =
588  const_cast<pointMesh&>(pointMesh::New(mesh));
589  PtrList<pointScalarField> pointScalarFields;
590  readFields
591  (
592  haveMesh,
593  pMesh,
594  subsetterPtr,
595  objects,
596  pointScalarFields
597  );
598 
599  PtrList<pointVectorField> pointVectorFields;
600  readFields
601  (
602  haveMesh,
603  pMesh,
604  subsetterPtr,
605  objects,
606  pointVectorFields
607  );
608 
609  PtrList<pointSphericalTensorField> pointSphereTensorFields;
610  readFields
611  (
612  haveMesh,
613  pMesh,
614  subsetterPtr,
615  objects,
616  pointSphereTensorFields
617  );
618 
619  PtrList<pointSymmTensorField> pointSymmTensorFields;
620  readFields
621  (
622  haveMesh,
623  pMesh,
624  subsetterPtr,
625  objects,
626  pointSymmTensorFields
627  );
628 
629  PtrList<pointTensorField> pointTensorFields;
630  readFields
631  (
632  haveMesh,
633  pMesh,
634  subsetterPtr,
635  objects,
636  pointTensorFields
637  );
638 
639  // Debugging: Create additional volField that will be mapped.
640  // Used to test correctness of mapping
641  // volVectorField mapCc("mapCc", 1*mesh.C());
642 
643  // Mesh distribution engine
644  fvMeshDistribute distributor(mesh);
645 
646  // Pout<< "Wanted distribution:"
647  // << distributor.countCells(finalDecomp) << nl << endl;
648 
649  // Do actual sending/receiving of mesh
650  autoPtr<polyDistributionMap> map = distributor.distribute(finalDecomp);
651 
653  // map().distributeFaceData(faceCc);
654 
655 
656  // Print some statistics
657  Info<< "After distribution:" << endl;
658  printMeshData(mesh);
659 
660 
661  if (!overwrite)
662  {
663  runTime++;
664  }
665  else
666  {
667  mesh.setInstance(masterInstDir);
668  }
669  Info<< "Writing redistributed mesh to " << runTime.name() << nl << endl;
670  mesh.write();
671 
672 
673  // Print nice message
674  // ~~~~~~~~~~~~~~~~~~
675 
676  // Determine which processors remain so we can print nice final message.
677  labelList nFaces(Pstream::nProcs());
678  nFaces[Pstream::myProcNo()] = mesh.nFaces();
679  Pstream::gatherList(nFaces);
680  Pstream::scatterList(nFaces);
681 
682  Info<< nl
683  << "You can pick up the redecomposed mesh from the polyMesh directory"
684  << " in " << runTime.name() << "." << nl
685  << "If you redecomposed the mesh to less processors you can delete"
686  << nl
687  << "the processor directories with 0 sized meshes in them." << nl
688  << "Below is a sample set of commands to do this."
689  << " Take care when issuing these" << nl
690  << "commands." << nl << endl;
691 
692  forAll(nFaces, proci)
693  {
694  fileName procDir = "processor" + name(proci);
695 
696  if (nFaces[proci] == 0)
697  {
698  Info<< " rm -r " << procDir.c_str() << nl;
699  }
700  else
701  {
702  fileName timeDir = procDir/runTime.name()/meshSubDir;
703  fileName constDir = procDir/runTime.constant()/meshSubDir;
704 
705  Info<< " rm -r " << constDir.c_str() << nl
706  << " mv " << timeDir.c_str()
707  << ' ' << constDir.c_str() << nl;
708  }
709  }
710  Info<< endl;
711 
712 
713  Info<< "End\n" << endl;
714 
715  return 0;
716 }
717 
718 
719 // ************************************************************************* //
Inter-processor communication reduction functions.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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:367
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:223
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 optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:255
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:284
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:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:426
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
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.
const Time & time() const
Return time.
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:270
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:273
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
label nInternalFaces() const
const vectorField & cellCentres() const
label nCells() 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::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#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:234
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:258
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:218
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)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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:267
autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io)
Load (if it exists) or create zero cell mesh given an IOobject:
objects
word meshPath
Definition: setMeshPath.H:1
const bool overwrite
Definition: setNoOverwrite.H:1
Foam::argList args(argc, argv)