decomposePar.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  decomposePar
26 
27 Description
28  Automatically decomposes a mesh and fields of a case for parallel
29  execution of OpenFOAM.
30 
31 Usage
32  \b decomposePar [OPTION]
33 
34  Options:
35  - \par -cellProc
36  Write cell processor indices as a volScalarField::Internal for
37  post-processing.
38 
39  - \par -region <regionName> \n
40  Decompose named region. Does not check for existence of processor*.
41 
42  - \par -allRegions \n
43  Decompose all regions in regionSolvers. Does not check for
44  existence of processor*.
45 
46  - \par -copyZero \n
47  Copy \a 0 directory to processor* rather than decompose the fields.
48 
49  - \par -copyUniform \n
50  Copy any \a uniform directories too.
51 
52  - \par -constant
53  Decompose mesh and fields in the constant directory.
54 
55  - \par -time xxx:yyy \n
56  Override controlDict settings and decompose selected times.
57 
58  - \par -fields \n
59  Use existing geometry decomposition and convert fields only.
60 
61  - \par -noSets \n
62  Skip decomposing cellSets, faceSets, pointSets.
63 
64  - \par -force \n
65  Remove any existing \a processor subdirectories before decomposing the
66  geometry.
67 
68 \*---------------------------------------------------------------------------*/
69 
70 #include "argList.H"
71 #include "timeSelector.H"
72 #include "IOobjectList.H"
73 #include "processorRunTimes.H"
75 #include "decompositionMethod.H"
76 #include "fvFieldDecomposer.H"
77 #include "pointFieldDecomposer.H"
79 
80 using namespace Foam;
81 
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84 namespace Foam
85 {
86 
87 bool haveUniform
88 (
89  const processorRunTimes& runTimes,
90  const word& regionDir = word::null
91 )
92 {
93  return
95  (
96  runTimes.completeTime().timePath()/regionDir/"uniform"
97  );
98 }
99 
100 
101 void decomposeUniform
102 (
103  const bool copyUniform,
104  const processorRunTimes& runTimes,
105  const word& regionDir = word::null
106 )
107 {
108  const fileName uniformDir(regionDir/"uniform");
109 
110  forAll(runTimes.procTimes(), proci)
111  {
112  const fileName procTimePath =
113  fileHandler().filePath(runTimes.procTimes()[proci].timePath());
114 
115  if (!fileHandler().isDir(procTimePath))
116  {
117  fileHandler().mkDir(procTimePath);
118  }
119 
120  if (copyUniform)
121  {
122  if (!fileHandler().exists(procTimePath/uniformDir))
123  {
124  fileHandler().cp
125  (
126  runTimes.completeTime().timePath()/uniformDir,
127  procTimePath/uniformDir
128  );
129  }
130  }
131  else
132  {
133  // Link with relative paths
134  string parentPath = string("..")/"..";
135 
136  if (regionDir != word::null)
137  {
138  parentPath = parentPath/"..";
139  }
140 
141  fileName currentDir(cwd());
142 
143  chDir(procTimePath);
144 
145  if (!fileHandler().exists(uniformDir))
146  {
147  fileHandler().ln
148  (
149  parentPath/runTimes.completeTime().name()/uniformDir,
150  uniformDir
151  );
152  }
153 
154  chDir(currentDir);
155  }
156  }
157 }
158 
159 
160 void writeDecomposition(const domainDecomposition& meshes)
161 {
162  // Write as volScalarField::Internal for postprocessing.
163  volScalarField::Internal cellProc
164  (
165  IOobject
166  (
167  "cellProc",
168  meshes.completeMesh().time().name(),
169  meshes.completeMesh(),
172  ),
173  meshes.completeMesh(),
174  dimless,
175  scalarField(scalarList(meshes.cellProc()))
176  );
177 
178  cellProc.write();
179 
180  Info<< "Wrote decomposition as volScalarField::Internal to "
181  << cellProc.name() << " for use in postprocessing"
182  << endl;
183 }
184 
185 
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 namespace Foam
192 {
193 
194 class delayedNewLine
195 {
196  mutable bool first_;
197 
198 public:
199 
200  delayedNewLine()
201  :
202  first_(true)
203  {}
204 
205  friend Ostream& operator<<(Ostream& os, const delayedNewLine& dnl)
206  {
207  if (!dnl.first_) os << nl;
208  dnl.first_ = false;
209  return os;
210  }
211 };
212 
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 int main(int argc, char *argv[])
219 {
221  (
222  "decompose a mesh and fields of a case for parallel execution"
223  );
224 
226  #include "addRegionOption.H"
227  #include "addAllRegionsOption.H"
229  (
230  "cellProc",
231  "write cell processor indices as a volScalarField::Internal for "
232  "post-processing"
233  );
235  (
236  "copyZero",
237  "Copy \a 0 directory to processor* rather than decompose the fields"
238  );
240  (
241  "copyUniform",
242  "copy any uniform/ directories too"
243  );
245  (
246  "fields",
247  "use existing geometry decomposition and convert fields only"
248  );
250  (
251  "noFields",
252  "opposite of -fields; only decompose geometry"
253  );
255  (
256  "noSets",
257  "skip decomposing cellSets, faceSets, pointSets"
258  );
260  (
261  "force",
262  "remove existing processor*/ subdirs before decomposing the geometry"
263  );
264 
265  // Include explicit constant option, execute from zero by default
266  timeSelector::addOptions(true, false);
267 
268  #include "setRootCase.H"
269 
270  const bool region = args.optionFound("region");
271  const bool writeCellProc = args.optionFound("cellProc");
272  const bool copyZero = args.optionFound("copyZero");
273  const bool copyUniform = args.optionFound("copyUniform");
274  const bool decomposeFieldsOnly = args.optionFound("fields");
275  const bool decomposeGeomOnly = args.optionFound("noFields");
276  const bool decomposeSets = !args.optionFound("noSets");
277  const bool forceOverwrite = args.optionFound("force");
278 
279  if (decomposeGeomOnly)
280  {
281  Info<< "Skipping decomposing fields" << nl << endl;
282 
283  if (decomposeFieldsOnly || copyZero)
284  {
286  << "Cannot combine geometry-only decomposition (-noFields)"
287  << " with field decomposition (-fields or -copyZero)"
288  << exit(FatalError);
289  }
290  }
291 
292  // Set time from database
293  Info<< "Create time" << nl << endl;
295  const Time& runTime = runTimes.completeTime();
296 
297  // Allow override of time
298  const instantList times = runTimes.selectComplete(args);
299 
300  #include "setRegionNames.H"
301 
302  // Remove existing processor directories if requested
303  if (forceOverwrite)
304  {
305  if (region)
306  {
308  << "Cannot force the decomposition of a single region"
309  << exit(FatalError);
310  }
311 
312  const label nProcs0 =
313  fileHandler().nProcs(runTimes.completeTime().path());
314 
315  Info<< "Removing " << nProcs0
316  << " existing processor directories" << nl << endl;
317 
318  // Remove existing processor directories
319  const fileNameList dirs
320  (
322  (
323  runTimes.completeTime().path(),
325  )
326  );
327  forAllReverse(dirs, diri)
328  {
329  const fileName& d = dirs[diri];
330 
331  // Starts with 'processors'
332  if (d.find("processors") == 0)
333  {
334  if (fileHandler().exists(d))
335  {
336  fileHandler().rmDir(d);
337  }
338  }
339 
340  // Starts with 'processor'
341  if (d.find("processor") == 0)
342  {
343  // Check that integer after processor
344  fileName num(d.substr(9));
345  label proci = -1;
346  if (Foam::read(num.c_str(), proci))
347  {
348  if (fileHandler().exists(d))
349  {
350  fileHandler().rmDir(d);
351  }
352  }
353  }
354  }
355 
356  // Flush file handler to clear any detected processor directories
357  fileHandler().flush();
358  }
359 
360  // Check the specified number of processes is consistent with any existing
361  // processor directories
362  {
363  const label nProcs0 =
364  fileHandler().nProcs(runTimes.completeTime().path());
365 
366  if (nProcs0 && nProcs0 != runTimes.nProcs())
367  {
369  << "Case is already decomposed with " << nProcs0
370  << " domains, use the -force option or manually" << nl
371  << "remove processor directories before decomposing. e.g.,"
372  << nl
373  << " rm -rf " << runTimes.completeTime().path().c_str()
374  << "/processor*"
375  << nl
376  << exit(FatalError);
377  }
378  }
379 
380  // Get the decomposition dictionary
381  const dictionary decomposeParDict =
383 
384  // Check existing decomposition
385  forAll(regionNames, regioni)
386  {
387  const word& regionName = regionNames[regioni];
388  const word regionDir =
390 
391  // Determine the existing processor count directly
392  const label nProcs =
393  fileHandler().nProcs(runTimes.completeTime().path(), regionDir);
394 
395  // Get requested numberOfSubdomains
396  const label nDomains =
397  decomposeParDict.lookup<label>("numberOfSubdomains");
398 
399  // Give file handler a chance to determine the output directory
400  const_cast<fileOperation&>(fileHandler()).setNProcs(nDomains);
401 
402  // Sanity check number of processors in a previously decomposed case
403  if (decomposeFieldsOnly && nProcs != nDomains)
404  {
406  << "Specified -fields, but the case was decomposed with "
407  << nProcs << " domains" << nl << "instead of " << nDomains
408  << " domains as specified in decomposeParDict" << nl
409  << exit(FatalError);
410  }
411  }
412 
413  // Create meshes
414  multiDomainDecomposition regionMeshes(runTimes, regionNames);
415  if
416  (
417  !(decomposeFieldsOnly && copyZero)
418  && regionMeshes.readDecompose(decomposeSets)
419  )
420  {
421  Info<< endl;
422 
423  if (writeCellProc)
424  {
425  forAll(regionNames, regioni)
426  {
427  writeDecomposition(regionMeshes[regioni]());
428  Info<< endl;
429  fileHandler().flush();
430  }
431  }
432  }
433 
434  // Get flag to determine whether or not to distribute uniform data
435  const bool distributed =
436  decomposeParDict.lookupOrDefault<bool>("distributed", false);
437 
438  // Loop over all times
439  forAll(times, timei)
440  {
441  // Set the time
442  runTimes.setTime(times[timei], timei);
443 
444  Info<< "Time = " << runTimes.completeTime().userTimeName()
445  << nl << endl;
446 
447  // Update the meshes, if necessary
448  const fvMesh::readUpdateState stat =
449  !(decomposeFieldsOnly && copyZero)
450  ? regionMeshes.readUpdateDecompose()
452  if (stat >= fvMesh::TOPO_CHANGE) Info<< endl;
453 
454  // Write the mesh out (if anything has changed), if necessary
455  if (!decomposeFieldsOnly)
456  {
457  regionMeshes.writeProcs(decomposeSets);
458  }
459 
460  // Write the decomposition, if necessary
461  forAll(regionNames, regioni)
462  {
463  if (writeCellProc && stat >= fvMesh::TOPO_CHANGE)
464  {
465  writeDecomposition(regionMeshes[regioni]());
466  Info<< endl;
467  fileHandler().flush();
468  }
469  }
470 
471  // If only decomposing geometry then there is no more to do
472  if (decomposeGeomOnly)
473  {
474  continue;
475  }
476 
477  // If copying from zero then just copy everything from the <time>
478  // directory to the processor*/<time> directories without altering them
479  if (copyZero)
480  {
481  const fileName completeTimePath =
482  runTimes.completeTime().timePath();
483 
484  fileName prevProcTimePath;
485  for (label proci = 0; proci < runTimes.nProcs(); proci++)
486  {
487  const Time& procRunTime = runTimes.procTimes()[proci];
488 
489  if (fileHandler().isDir(completeTimePath))
490  {
491  const fileName procTimePath
492  (
493  fileHandler().objectPath
494  (
495  IOobject
496  (
497  "",
498  procRunTime.name(),
499  procRunTime
500  ),
501  word::null
502  )
503  );
504 
505  if (procTimePath != prevProcTimePath)
506  {
507  Info<< "Processor " << proci
508  << ": copying " << completeTimePath << nl
509  << " to " << procTimePath << endl;
510  fileHandler().cp(completeTimePath, procTimePath);
511  prevProcTimePath = procTimePath;
512  }
513  }
514  }
515 
516  Info<< endl;
517 
518  continue;
519  }
520 
521  // Otherwise we are doing full region-by-region decomposition of all
522  // the available fields
523  forAll(regionNames, regioni)
524  {
525  const word& regionName = regionNames[regioni];
526  const word regionDir =
528 
529  const delayedNewLine dnl;
530 
531  // Prefixed scope
532  {
533  const RegionRef<domainDecomposition> meshes =
534  regionMeshes[regioni];
535 
536  // Search for objects at this time
538  (
539  meshes().completeMesh(),
540  runTimes.completeTime().name()
541  );
542 
543  {
544  Info<< dnl << "Decomposing FV fields" << endl;
545 
547  {
548  fvFieldDecomposer fvDecomposer
549  (
550  meshes().completeMesh(),
551  meshes().procMeshes(),
552  meshes().procFaceAddressing(),
553  meshes().procCellAddressing(),
554  meshes().procFaceAddressingBf()
555  );
556 
557  #define DO_FV_VOL_INTERNAL_FIELDS_TYPE(Type, nullArg) \
558  fvDecomposer.decomposeVolInternalFields<Type> \
559  (objects);
560  FOR_ALL_FIELD_TYPES(DO_FV_VOL_INTERNAL_FIELDS_TYPE)
561  #undef DO_FV_VOL_INTERNAL_FIELDS_TYPE
562 
563  #define DO_FV_VOL_FIELDS_TYPE(Type, nullArg) \
564  fvDecomposer.decomposeVolFields<Type> \
565  (objects);
566  FOR_ALL_FIELD_TYPES(DO_FV_VOL_FIELDS_TYPE)
567  #undef DO_FV_VOL_FIELDS_TYPE
568 
569  #define DO_FV_SURFACE_FIELDS_TYPE(Type, nullArg) \
570  fvDecomposer.decomposeFvSurfaceFields<Type> \
571  (objects);
572  FOR_ALL_FIELD_TYPES(DO_FV_SURFACE_FIELDS_TYPE)
573  #undef DO_FV_SURFACE_FIELDS_TYPE
574  }
575  else
576  {
577  Info<< dnl << " (no FV fields)" << endl;
578  }
579  }
580 
581  {
582  Info<< dnl << "Decomposing point fields" << endl;
583 
585  {
586  pointFieldDecomposer pointDecomposer
587  (
588  pointMesh::New(meshes().completeMesh()),
589  meshes().procMeshes(),
590  meshes().procPointAddressing()
591  );
592 
593  #define DO_POINT_FIELDS_TYPE(Type, nullArg) \
594  pointDecomposer.decomposeFields<Type> \
595  (objects);
596  FOR_ALL_FIELD_TYPES(DO_POINT_FIELDS_TYPE)
597  #undef DO_POINT_FIELDS_TYPE
598  }
599  else
600  {
601  Info<< dnl << " (no point fields)" << endl;
602  }
603  }
604 
605  {
606  // Find cloud directories
607  fileNameList cloudDirs
608  (
610  (
611  runTimes.completeTime().timePath()
612  /regionDir
613  /cloud::prefix,
615  )
616  );
617 
618  // Add objects in any found cloud directories
619  HashTable<IOobjectList> cloudsObjects;
620  forAll(cloudDirs, i)
621  {
622  // Do local scan for valid cloud objects
623  IOobjectList cloudObjs
624  (
625  meshes().completeMesh(),
626  runTimes.completeTime().name(),
627  cloud::prefix/cloudDirs[i],
630  false
631  );
632 
633  // If "positions" is present, then add to the table
634  if (cloudObjs.lookup(word("positions")))
635  {
636  cloudsObjects.insert(cloudDirs[i], cloudObjs);
637  }
638  }
639 
640  // Decompose the objects found above
641  if (cloudsObjects.size())
642  {
644  (
646  cloudsObjects,
647  iter
648  )
649  {
650  const word cloudName =
651  string::validate<word>(iter.key());
652 
653  const IOobjectList& cloudObjects = iter();
654 
655  Info<< dnl << "Decomposing lagrangian fields for "
656  << "cloud " << cloudName << endl;
657 
658  if
659  (
661  (
662  cloudObjects
663  )
664  )
665  {
667  lagrangianDecomposer
668  (
669  meshes().completeMesh(),
670  meshes().procMeshes(),
671  meshes().procFaceAddressing(),
672  meshes().procCellAddressing(),
673  cloudName
674  );
675 
676  #define DO_CLOUD_FIELDS_TYPE(Type, nullArg) \
677  lagrangianDecomposer.decomposeFields<Type> \
678  (cloudObjects);
679  DO_CLOUD_FIELDS_TYPE(label, );
680  FOR_ALL_FIELD_TYPES(DO_CLOUD_FIELDS_TYPE)
681  #undef DO_CLOUD_FIELDS_TYPE
682  }
683  else
684  {
685  Info<< dnl << " (no lagrangian fields)"
686  << endl;
687  }
688  }
689  }
690  }
691  }
692 
693  Info<< dnl;
694  }
695 
696  // Distribute the uniform directory
697  if (haveUniform(runTimes))
698  {
699  Info<< "Distributing uniform files" << endl;
700 
701  decomposeUniform
702  (
703  copyUniform || distributed,
704  runTimes
705  );
706 
707  Info<< endl;
708  }
709 
710  if (regionNames == wordList(1, polyMesh::defaultRegion)) continue;
711 
712  // Distribute the region uniform directories
713  forAll(regionNames, regioni)
714  {
715  const word& regionName = regionNames[regioni];
716  const word regionDir =
718 
719  if (haveUniform(runTimes, regionDir))
720  {
721  // Prefixed scope
722  {
723  const RegionRef<domainDecomposition> meshes =
724  regionMeshes[regioni];
725 
726  Info<< "Distributing uniform files" << endl;
727 
728  decomposeUniform
729  (
730  copyUniform || distributed,
731  runTimes,
732  regionDir
733  );
734  }
735 
736  Info<< endl;
737  }
738  }
739  }
740 
741  Info<< "End" << nl << endl;
742 
743  return 0;
744 }
745 
746 
747 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
#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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An STL-conforming hash table.
Definition: HashTable.H:127
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
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
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
fileName path() const
Explicitly inherit path from TimePaths to disambiguate from.
Definition: TimePaths.H:138
word userTimeName() const
Return current user time name with units.
Definition: Time.C:856
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:208
fileName timePath() const
Return current time path.
Definition: Time.H:281
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
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 void noParallel()
Remove the parallel options.
Definition: argList.C:175
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
static IOdictionary 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:162
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 word & name() const
Return const reference to name.
Automatic domain decomposition class for finite-volume meshes.
const labelList & cellProc() const
Return the distribution as an FV field for writing.
const fvMesh & completeMesh() const
Access the global mesh.
A class for handling file names.
Definition: fileName.H:82
virtual bool rmDir(const fileName &) const =0
Remove a directory and its contents.
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
virtual bool ln(const fileName &src, const fileName &dst) const =0
Create a softlink. dst should not exist. Returns true if.
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
virtual bool isDir(const fileName &, const bool followLink=true) const =0
Does the name exist as a directory in the file system?
virtual bool mkDir(const fileName &, mode_t=0777) const =0
Make directory.
virtual fileName filePath(const bool globalFile, const IOobject &, const word &typeName) const =0
Search for an object. globalFile : also check undecomposed case.
Finite Volume volume and surface field decomposer.
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:108
@ TOPO_CHANGE
Definition: fvMesh.H:112
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
Point field decomposer.
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
instantList selectComplete(const argList &args)
Select the time.
label nProcs() const
Return the number of processors.
const PtrList< Time > & procTimes() const
Access the processor run times.
void setTime(const instant &inst, const label newIndex)
Set the time.
const Time & completeTime() const
Access the complete run time.
A class for handling character strings derived from std::string.
Definition: string.H:79
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
fileName cwd()
Return current working directory path name.
Definition: POSIX.C:241
List< word > wordList
A List of words.
Definition: fileName.H:54
bool read(const char *, int32_t &)
Definition: int32IO.C:85
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 dimensionSet dimless
const word & regionName(const solver &region)
Definition: solver.H:209
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
messageStream Info
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
static const char nl
Definition: Ostream.H:266
fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:662
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
Definition: POSIX.C:284
objects
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)
const word cloudName(propsDict.lookup("cloudName"))