reconstructPar.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  reconstructPar
26 
27 Description
28  Reconstructs fields of a case that is decomposed for parallel
29  execution of OpenFOAM.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "timeSelector.H"
35 #include "IOobjectList.H"
36 #include "processorRunTimes.H"
38 #include "fvFieldReconstructor.H"
42 
43 using namespace Foam;
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 bool haveUniform
51 (
52  const processorRunTimes& runTimes,
53  const word& regionDir = word::null
54 )
55 {
56  return
58  (
59  fileHandler().filePath
60  (
61  runTimes.procTimes()[0].timePath()/regionDir/"uniform"
62  )
63  );
64 }
65 
66 
67 void reconstructUniform
68 (
69  const processorRunTimes& runTimes,
70  const word& regionDir = word::null
71 )
72 {
73  fileHandler().cp
74  (
75  fileHandler().filePath
76  (
77  runTimes.procTimes()[0].timePath()/regionDir/"uniform"
78  ),
79  runTimes.completeTime().timePath()/regionDir
80  );
81 }
82 
83 
84 void writeDecomposition(const domainDecomposition& meshes)
85 {
86  // Write as volScalarField::Internal for postprocessing.
88  (
89  IOobject
90  (
91  "cellProc",
92  meshes.completeMesh().time().name(),
93  meshes.completeMesh(),
96  ),
97  meshes.completeMesh(),
98  dimless,
100  );
101 
102  cellProc.write();
103 
104  Info<< "Wrote decomposition as volScalarField::Internal to "
105  << cellProc.name() << " for use in postprocessing"
106  << endl;
107 }
108 
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113 
114 namespace Foam
115 {
116 
117 class delayedNewLine
118 {
119  mutable bool first_;
120 
121 public:
122 
123  delayedNewLine()
124  :
125  first_(true)
126  {}
127 
128  friend Ostream& operator<<(Ostream& os, const delayedNewLine& dnl)
129  {
130  if (!dnl.first_) os << nl;
131  dnl.first_ = false;
132  return os;
133  }
134 };
135 
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 int main(int argc, char *argv[])
142 {
144  (
145  "Reconstruct fields of a parallel case"
146  );
147 
149  #include "addMeshOption.H"
150  #include "addRegionOption.H"
151  #include "addAllRegionsOption.H"
153  (
154  "cellProc",
155  "write cell processor indices as a volScalarField::Internal for "
156  "post-processing"
157  );
159  (
160  "fields",
161  "list",
162  "specify a list of fields to be reconstructed. Eg, '(U T p)' - "
163  "regular expressions not currently supported"
164  );
166  (
167  "noFields",
168  "skip reconstructing fields"
169  );
171  (
172  "lagrangianFields",
173  "list",
174  "specify a list of lagrangian fields to be reconstructed. Eg, '(U d)' -"
175  "regular expressions not currently supported, "
176  "positions always included"
177  );
179  (
180  "noLagrangian",
181  "skip reconstructing lagrangian positions and fields"
182  );
184  (
185  "noSets",
186  "skip reconstructing cellSets, faceSets, pointSets"
187  );
189  (
190  "newTimes",
191  "only reconstruct new times (i.e. that do not exist already)"
192  );
194  (
195  "rm",
196  "remove processor time directories after reconstruction"
197  );
198 
199  // Include explicit constant options, and explicit zero option (to prevent
200  // the user accidentally trashing the initial fields)
201  timeSelector::addOptions(true, true);
202 
203  #include "setRootCase.H"
204  #include "setMeshPath.H"
205 
206  const bool writeCellProc = args.optionFound("cellProc");
207 
208  HashSet<word> selectedFields;
209  if (args.optionFound("fields"))
210  {
211  args.optionLookup("fields")() >> selectedFields;
212  }
213 
214  const bool noFields = args.optionFound("noFields");
215 
216  if (noFields)
217  {
218  Info<< "Skipping reconstructing fields" << nl << endl;
219  }
220 
221  const bool noLagrangian = args.optionFound("noLagrangian");
222 
223  if (noLagrangian)
224  {
225  Info<< "Skipping reconstructing lagrangian positions and fields"
226  << nl << endl;
227  }
228 
229  const bool noReconstructSets = args.optionFound("noSets");
230 
231  if (noReconstructSets)
232  {
233  Info<< "Skipping reconstructing cellSets, faceSets and pointSets"
234  << nl << endl;
235  }
236 
237  HashSet<word> selectedLagrangianFields;
238  if (args.optionFound("lagrangianFields"))
239  {
240  if (noLagrangian)
241  {
243  << "Cannot specify noLagrangian and lagrangianFields "
244  << "options together"
245  << exit(FatalError);
246  }
247 
248  args.optionLookup("lagrangianFields")() >> selectedLagrangianFields;
249  }
250 
251  // Set time from database
252  Info<< "Create time" << nl << endl;
254 
255  // Get the times to reconstruct
256  instantList times = runTimes.selectProc(args);
257 
258  const Time& runTime = runTimes.procTimes()[0];
259 
260  #include "setRegionNames.H"
261 
262  // Determine the processor count
263  const label nProcs = fileHandler().nProcs
264  (
265  args.path(),
267  ? word::null
268  : regionNames[0]
269  );
270 
271  if (!nProcs)
272  {
274  << "No processor* directories found"
275  << exit(FatalError);
276  }
277 
278  // Warn fileHandler of number of processors
279  const_cast<fileOperation&>(fileHandler()).setNProcs(nProcs);
280 
281  // Quit if no times
282  if (times.empty())
283  {
284  WarningInFunction << "No times selected" << nl << endl;
285  exit(1);
286  }
287 
288  // If only reconstructing new times then filter out existing times
289  if (args.optionFound("newTimes"))
290  {
291  // Get all existing times
292  const instantList existingTimes = runTimes.completeTime().times();
293 
294  // Put into a set
295  HashSet<word> existingTimesSet;
296  existingTimesSet.resize(2*existingTimes.size());
297  forAll(existingTimes, i)
298  {
299  existingTimesSet.insert(existingTimes[i].name());
300  }
301 
302  // Remove times from the existing time set by shuffling up
303  label timei = 0;
304  forAll(times, timej)
305  {
306  if (!existingTimesSet.found(times[timej].name()))
307  {
308  times[timei ++] = times[timej];
309  }
310  }
311  times.resize(timei);
312  }
313 
314  // Quit if no times
315  if (times.empty())
316  {
317  Info<< "All times already reconstructed" << nl << nl
318  << "End" << nl << endl;
319  return 0;
320  }
321 
322  // Create meshes
323  multiDomainDecomposition regionMeshes(runTimes, meshPath, regionNames);
324  if (regionMeshes.readReconstruct(!noReconstructSets))
325  {
326  Info<< endl;
327 
328  if (writeCellProc)
329  {
330  forAll(regionNames, regioni)
331  {
332  writeDecomposition(regionMeshes[regioni]());
333  Info<< endl;
334  fileHandler().flush();
335  }
336  }
337  }
338 
339  // Loop over all times
340  forAll(times, timei)
341  {
342  // Set the time
343  runTimes.setTime(times[timei], timei);
344 
345  Info<< "Time = " << runTimes.completeTime().userTimeName()
346  << nl << endl;
347 
348  // Update the meshes
349  const fvMesh::readUpdateState stat =
350  regionMeshes.readUpdateReconstruct();
351  if (stat >= fvMesh::TOPO_CHANGE) Info<< endl;
352 
353  // Write the mesh out (if anything has changed)
354  regionMeshes.writeComplete(!noReconstructSets);
355 
356  // Write the decomposition, if necessary
357  forAll(regionNames, regioni)
358  {
359  if (writeCellProc && stat >= fvMesh::TOPO_CHANGE)
360  {
361  writeDecomposition(regionMeshes[regioni]());
362  Info<< endl;
363  fileHandler().flush();
364  }
365  }
366 
367  // Do a region-by-region reconstruction of all the available fields
368  forAll(regionNames, regioni)
369  {
370  const word& regionName = regionNames[regioni];
371  const word regionDir =
373 
374  const delayedNewLine dnl;
375 
376  // Prefixed scope
377  {
378  const RegionRef<domainDecomposition> meshes =
379  regionMeshes[regioni];
380 
381  // Search for objects at this time
383  (
384  meshes().procMeshes()[0],
385  runTimes.procTimes()[0].name()
386  );
387 
388  if (!noFields)
389  {
390  Info<< dnl << "Reconstructing FV fields" << endl;
391 
392  if
393  (
394  fvFieldReconstructor::reconstructs
395  (
396  objects,
397  selectedFields
398  )
399  )
400  {
401  fvFieldReconstructor fvReconstructor
402  (
403  meshes().completeMesh(),
404  meshes().procMeshes(),
405  meshes().procFaceAddressing(),
406  meshes().procCellAddressing(),
407  meshes().procFaceAddressingBf()
408  );
409 
410  #define DO_FV_VOL_INTERNAL_FIELDS_TYPE(Type, nullArg) \
411  fvReconstructor.reconstructVolInternalFields<Type> \
412  (objects, selectedFields);
413  FOR_ALL_FIELD_TYPES(DO_FV_VOL_INTERNAL_FIELDS_TYPE)
414  #undef DO_FV_VOL_INTERNAL_FIELDS_TYPE
415 
416  #define DO_FV_VOL_FIELDS_TYPE(Type, nullArg) \
417  fvReconstructor.reconstructVolFields<Type> \
418  (objects, selectedFields);
419  FOR_ALL_FIELD_TYPES(DO_FV_VOL_FIELDS_TYPE)
420  #undef DO_FV_VOL_FIELDS_TYPE
421 
422  #define DO_FV_SURFACE_FIELDS_TYPE(Type, nullArg) \
423  fvReconstructor.reconstructFvSurfaceFields<Type> \
424  (objects, selectedFields);
425  FOR_ALL_FIELD_TYPES(DO_FV_SURFACE_FIELDS_TYPE)
426  #undef DO_FV_SURFACE_FIELDS_TYPE
427  }
428  else
429  {
430  Info<< dnl << " (no FV fields)" << endl;
431  }
432  }
433 
434  if (!noFields)
435  {
436  Info<< dnl << "Reconstructing point fields" << endl;
437 
438  if
439  (
440  pointFieldReconstructor::reconstructs
441  (
442  objects,
443  selectedFields
444  )
445  )
446  {
447  pointFieldReconstructor pointReconstructor
448  (
449  pointMesh::New(meshes().completeMesh()),
450  meshes().procMeshes(),
451  meshes().procPointAddressing()
452  );
453 
454  #define DO_POINT_FIELDS_TYPE(Type, nullArg) \
455  pointReconstructor.reconstructFields<Type> \
456  (objects, selectedFields);
458  #undef DO_POINT_FIELDS_TYPE
459  }
460  else
461  {
462  Info<< dnl << " (no point fields)" << endl;
463  }
464  }
465 
466  if (!noLagrangian)
467  {
468  // Search for clouds that exist on any processor and add
469  // them into this table of cloud objects
470  HashTable<IOobjectList> cloudsObjects;
471  forAll(runTimes.procTimes(), proci)
472  {
473  // Find cloud directories
474  fileNameList cloudDirs
475  (
477  (
478  fileHandler().filePath
479  (
480  runTimes.procTimes()[proci].timePath()
481  /regionDir
483  ),
485  )
486  );
487 
488  // Add objects in any found cloud directories
489  forAll(cloudDirs, i)
490  {
491  // Pass if we already have an objects for this name
493  cloudsObjects.find(cloudDirs[i]);
494  if (iter != cloudsObjects.end()) continue;
495 
496  // Do local scan for valid cloud objects
497  IOobjectList cloudObjs
498  (
499  meshes().procMeshes()[proci],
500  runTimes.procTimes()[proci].name(),
501  lagrangian::cloud::prefix/cloudDirs[i],
504  false
505  );
506 
507  // If "positions" is present, then add to the table
508  if (cloudObjs.lookup(word("positions")))
509  {
510  cloudsObjects.insert(cloudDirs[i], cloudObjs);
511  }
512  }
513  }
514 
515  // Reconstruct the objects found above
516  if (cloudsObjects.size())
517  {
519  (
521  cloudsObjects,
522  iter
523  )
524  {
525  const word cloudName =
526  string::validate<word>(iter.key());
527 
528  const IOobjectList& cloudObjects = iter();
529 
530  Info<< dnl << "Reconstructing lagrangian fields "
531  << "for cloud " << cloudName << endl;
532 
533  if
534  (
535  lagrangianFieldReconstructor::reconstructs
536  (
537  cloudObjects,
538  selectedLagrangianFields
539  )
540  )
541  {
543  lagrangianReconstructor
544  (
545  meshes().completeMesh(),
546  meshes().procMeshes(),
547  meshes().procFaceAddressing(),
548  meshes().procCellAddressing(),
549  cloudName
550  );
551 
552  lagrangianReconstructor.reconstructPositions();
553 
554  #define DO_CLOUD_FIELDS_TYPE(Type, nullArg) \
555  lagrangianReconstructor \
556  .reconstructFields<Type> \
557  (cloudObjects, selectedLagrangianFields);
558  DO_CLOUD_FIELDS_TYPE(label, )
559  FOR_ALL_FIELD_TYPES(DO_CLOUD_FIELDS_TYPE)
560  #undef DO_CLOUD_FIELDS_TYPE
561  }
562  else
563  {
564  Info<< dnl << " (no lagrangian fields)"
565  << endl;
566  }
567  }
568  }
569  }
570 
571  if (!noLagrangian)
572  {
573  // Search for Lagrangian meshes that exist on any processor
574  // and add them into this table of objects
575  HashTable<IOobjectList> LagrangianObjects;
576  forAll(runTimes.procTimes(), proci)
577  {
578  // Find Lagrangian directories
579  fileNameList LagrangianDirs
580  (
582  (
583  fileHandler().filePath
584  (
585  runTimes.procTimes()[proci].timePath()
586  /regionDir
588  ),
590  )
591  );
592 
593  // Add objects in any found Lagrangian directories
594  forAll(LagrangianDirs, i)
595  {
596  // Pass if we already have an objects for this name
597  if
598  (
599  LagrangianObjects.find(LagrangianDirs[i])
600  != LagrangianObjects.end()
601  ) continue;
602 
603  // Do local scan for valid Lagrangian objects
605  (
606  meshes().procMeshes()[proci],
607  runTimes.procTimes()[proci].name(),
608  LagrangianMesh::prefix/LagrangianDirs[i],
611  false
612  );
613 
614  // If coordinates or fields are present then add
615  // this set of objects to the table
616  if
617  (
619  || LagrangianFieldReconstructor::reconstructs
620  (
621  objects,
622  selectedLagrangianFields
623  )
624  )
625  {
626  LagrangianObjects.insert
627  (
628  LagrangianDirs[i],
629  objects
630  );
631  }
632  }
633  }
634 
635  // Reconstruct the objects found above
636  if (LagrangianObjects.size())
637  {
639  (
641  LagrangianObjects,
642  iter
643  )
644  {
645  const word LagrangianName =
646  string::validate<word>(iter.key());
647 
648  Info<< dnl << "Reconstructing Lagrangian fields "
649  << "for " << LagrangianName << endl;
650 
652  LagrangianReconstructor
653  (
654  meshes().completeMesh(),
655  meshes().procMeshes(),
656  meshes().procFaceAddressing(),
657  meshes().procCellAddressing(),
658  LagrangianName
659  );
660 
661  LagrangianReconstructor.reconstructPositions();
662 
663  if
664  (
665  LagrangianFieldReconstructor::reconstructs
666  (
667  iter(),
668  selectedLagrangianFields
669  )
670  )
671  {
672  #define DO_LAGRANGIAN_FIELDS_TYPE( \
673  Type, GeoField) \
674  LagrangianReconstructor \
675  .reconstructFields<GeoField<Type>> \
676  (iter(), selectedLagrangianFields);
678  (
679  label,
681  )
683  (
686  );
688  (
689  label,
691  )
693  (
696  );
697  #undef DO_LAGRANGIAN_FIELDS_TYPE
698 
699  // --> Note we don't have to explicitly
700  // reconstruct the dynamic variants of these
701  // fields as they are IO compatible with the
702  // non-dynamic fields
703  }
704  else
705  {
706  Info<< dnl << " (no Lagrangian fields)"
707  << endl;
708  }
709  }
710  }
711  }
712  }
713 
714  Info<< dnl;
715  }
716 
717  // Collect the uniform directory
718  if (haveUniform(runTimes))
719  {
720  Info<< "Collecting uniform files" << endl;
721 
722  reconstructUniform(runTimes);
723 
724  Info<< endl;
725  }
726 
728  {
729  // Collect the region uniform directories
730  forAll(regionNames, regioni)
731  {
732  const word& regionName = regionNames[regioni];
733  const word regionDir =
735  ? word::null
736  : regionName;
737 
738  if (haveUniform(runTimes, regionDir))
739  {
740  // Prefixed scope
741  {
742  const RegionRef<domainDecomposition> meshes =
743  regionMeshes[regioni];
744 
745  Info<< "Collecting uniform files" << endl;
746 
747  reconstructUniform(runTimes, regionDir);
748  }
749 
750  Info<< endl;
751  }
752  }
753  }
754 
755  if (args.optionFound("rm") && times[timei].name() != Time::constantName)
756  {
757  const bool allRegions = args.optionFound("allRegions");
758 
759  forAll(regionNames, regioni)
760  {
761  const word& regionName = regionNames[regioni];
762  const word regionDir =
763  allRegions || regionName == polyMesh::defaultRegion
764  ? word::null
765  : regionName;
766 
767  const RegionRef<domainDecomposition> meshes =
768  regionMeshes[regioni];
769 
770  Info<< "Removing processors time directory" << endl;
771 
772  for (label proci=0; proci<nProcs; proci++)
773  {
774  const fileName procTimePath = fileHandler().filePath
775  (
776  runTimes.procTimes()[proci].timePath()/regionDir
777  );
778 
779  if (isDir(procTimePath))
780  {
781  rmDir(procTimePath);
782  }
783  }
784  }
785 
786  Info<< endl;
787  }
788  }
789 
790  Info<< "End" << nl << endl;
791 
792  return 0;
793 }
794 
795 
796 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
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...
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
An STL-conforming const_iterator.
Definition: HashTable.H:498
const Key & key() const
Return the Key corresponding to the iterator.
Definition: HashTableI.H:254
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
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:167
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:506
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static const word coordinatesName
Name of the coordinates field.
static const word prefix
Instance prefix.
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static const word constantName
Definition: TimePaths.H:62
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
instantList times() const
Search the case for valid time directories.
Definition: Time.C:640
word userTimeName() const
Return current user time name with units.
Definition: Time.C:845
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:208
fileName timePath() const
Return current time path.
Definition: Time.H:275
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
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
fileName path() const
Return the path to the caseName.
Definition: argListI.H:66
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
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 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 fileName filePath(const bool globalFile, const IOobject &) const =0
Search for an object. globalFile : also check undecomposed case.
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?
Finite volume reconstructor for volume and surface fields.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:105
@ TOPO_CHANGE
Definition: fvMesh.H:109
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:69
Point field reconstructor.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:270
instantList selectProc(const argList &args)
Select the time.
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.
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
#define DO_LAGRANGIAN_FIELDS_TYPE(Type, nullArg)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
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:258
const dimensionSet dimless
const word & regionName(const solver &region)
Definition: solver.H:218
messageStream Info
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
static const char nl
Definition: Ostream.H:267
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
#define DO_POINT_FIELDS_TYPE(Type, nullArg)
objects
word meshPath
Definition: setMeshPath.H:1
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)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
const word cloudName(propsDict.lookup("cloudName"))