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-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  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"
80 
81 using namespace Foam;
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 namespace Foam
86 {
87 
88 bool haveUniform
89 (
90  const processorRunTimes& runTimes,
91  const word& regionDir = word::null
92 )
93 {
94  return
96  (
97  runTimes.completeTime().timePath()/regionDir/"uniform"
98  );
99 }
100 
101 
102 void decomposeUniform
103 (
104  const bool copyUniform,
105  const processorRunTimes& runTimes,
106  const word& regionDir = word::null
107 )
108 {
109  const fileName uniformDir(regionDir/"uniform");
110 
111  forAll(runTimes.procTimes(), proci)
112  {
113  const fileName procTimePath =
114  fileHandler().filePath(runTimes.procTimes()[proci].timePath());
115 
116  if (!fileHandler().isDir(procTimePath))
117  {
118  fileHandler().mkDir(procTimePath);
119  }
120 
121  if (copyUniform)
122  {
123  if (!fileHandler().exists(procTimePath/uniformDir))
124  {
125  fileHandler().cp
126  (
127  runTimes.completeTime().timePath()/uniformDir,
128  procTimePath/uniformDir
129  );
130  }
131  }
132  else
133  {
134  // Link with relative paths
135  string parentPath = string("..")/"..";
136 
137  if (regionDir != word::null)
138  {
139  parentPath = parentPath/"..";
140  }
141 
142  fileName currentDir(cwd());
143 
144  chDir(procTimePath);
145 
146  if (!fileHandler().exists(uniformDir))
147  {
148  fileHandler().ln
149  (
150  parentPath/runTimes.completeTime().name()/uniformDir,
151  uniformDir
152  );
153  }
154 
155  chDir(currentDir);
156  }
157  }
158 }
159 
160 
161 void writeDecomposition(const domainDecomposition& meshes)
162 {
163  // Write as volScalarField::Internal for postprocessing.
164  volScalarField::Internal cellProc
165  (
166  IOobject
167  (
168  "cellProc",
169  meshes.completeMesh().time().name(),
170  meshes.completeMesh(),
173  ),
174  meshes.completeMesh(),
175  dimless,
176  scalarField(scalarList(meshes.cellProc()))
177  );
178 
179  cellProc.write();
180 
181  Info<< "Wrote decomposition as volScalarField::Internal to "
182  << cellProc.name() << " for use in postprocessing"
183  << endl;
184 }
185 
186 
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 namespace Foam
193 {
194 
195 class delayedNewLine
196 {
197  mutable bool first_;
198 
199 public:
200 
201  delayedNewLine()
202  :
203  first_(true)
204  {}
205 
206  friend Ostream& operator<<(Ostream& os, const delayedNewLine& dnl)
207  {
208  if (!dnl.first_) os << nl;
209  dnl.first_ = false;
210  return os;
211  }
212 };
213 
214 }
215 
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 int main(int argc, char *argv[])
220 {
222  (
223  "decompose a mesh and fields of a case for parallel execution"
224  );
225 
227  #include "addMeshOption.H"
228  #include "addRegionOption.H"
229  #include "addAllRegionsOption.H"
231  (
232  "cellProc",
233  "write cell processor indices as a volScalarField::Internal for "
234  "post-processing"
235  );
237  (
238  "copyZero",
239  "Copy \a 0 directory to processor* rather than decompose the fields"
240  );
242  (
243  "copyUniform",
244  "copy any uniform/ directories too"
245  );
247  (
248  "fields",
249  "use existing geometry decomposition and convert fields only"
250  );
252  (
253  "noFields",
254  "opposite of -fields; only decompose geometry"
255  );
257  (
258  "noSets",
259  "skip decomposing cellSets, faceSets, pointSets"
260  );
262  (
263  "force",
264  "remove existing processor*/ subdirs before decomposing the geometry"
265  );
266 
267  // Include explicit constant option, execute from zero by default
268  timeSelector::addOptions(true, false);
269 
270  #include "setRootCase.H"
271  #include "setMeshPath.H"
272 
273  const bool region = args.optionFound("region");
274  const bool writeCellProc = args.optionFound("cellProc");
275  const bool copyZero = args.optionFound("copyZero");
276  const bool copyUniform = args.optionFound("copyUniform");
277  const bool decomposeFieldsOnly = args.optionFound("fields");
278  const bool decomposeGeomOnly = args.optionFound("noFields");
279  const bool decomposeSets = !args.optionFound("noSets");
280  const bool forceOverwrite = args.optionFound("force");
281 
282  if (decomposeGeomOnly)
283  {
284  Info<< "Skipping decomposing fields" << nl << endl;
285 
286  if (decomposeFieldsOnly || copyZero)
287  {
289  << "Cannot combine geometry-only decomposition (-noFields)"
290  << " with field decomposition (-fields or -copyZero)"
291  << exit(FatalError);
292  }
293  }
294 
295  // Set time from database
296  Info<< "Create time" << nl << endl;
298  const Time& runTime = runTimes.completeTime();
299 
300  // Allow override of time
301  const instantList times = runTimes.selectComplete(args);
302 
303  #include "setRegionNames.H"
304 
305  // Remove existing processor directories if requested
306  if (forceOverwrite)
307  {
308  if (region)
309  {
311  << "Cannot force the decomposition of a single region"
312  << exit(FatalError);
313  }
314 
315  const label nProcs0 =
316  fileHandler().nProcs(runTimes.completeTime().path());
317 
318  Info<< "Removing " << nProcs0
319  << " existing processor directories" << nl << endl;
320 
321  // Remove existing processor directories
322  const fileNameList dirs
323  (
325  (
326  runTimes.completeTime().path(),
328  )
329  );
330  forAllReverse(dirs, diri)
331  {
332  const fileName& d = dirs[diri];
333 
334  // Starts with 'processors'
335  if (d.find("processors") == 0)
336  {
337  if (fileHandler().exists(d))
338  {
339  fileHandler().rmDir(d);
340  }
341  }
342 
343  // Starts with 'processor'
344  if (d.find("processor") == 0)
345  {
346  // Check that integer after processor
347  fileName num(d.substr(9));
348  label proci = -1;
349  if (Foam::read(num.c_str(), proci))
350  {
351  if (fileHandler().exists(d))
352  {
353  fileHandler().rmDir(d);
354  }
355  }
356  }
357  }
358 
359  // Flush file handler to clear any detected processor directories
360  fileHandler().flush();
361  }
362 
363  // Check the specified number of processes is consistent with any existing
364  // processor directories
365  {
366  const label nProcs0 =
367  fileHandler().nProcs(runTimes.completeTime().path());
368 
369  if (nProcs0 && nProcs0 != runTimes.nProcs())
370  {
372  << "Case is already decomposed with " << nProcs0
373  << " domains, use the -force option or manually" << nl
374  << "remove processor directories before decomposing. e.g.,"
375  << nl
376  << " rm -rf " << runTimes.completeTime().path().c_str()
377  << "/processor*"
378  << nl
379  << exit(FatalError);
380  }
381  }
382 
383  // Get the decomposition dictionary
384  const dictionary decomposeParDict =
386 
387  // Check existing decomposition
388  forAll(regionNames, regioni)
389  {
390  const word& regionName = regionNames[regioni];
391  const word regionDir =
393 
394  // Determine the existing processor count directly
395  const label nProcs =
396  fileHandler().nProcs(runTimes.completeTime().path(), regionDir);
397 
398  // Get requested numberOfSubdomains
399  const label nDomains =
400  decomposeParDict.lookup<label>("numberOfSubdomains");
401 
402  // Give file handler a chance to determine the output directory
403  const_cast<fileOperation&>(fileHandler()).setNProcs(nDomains);
404 
405  // Sanity check number of processors in a previously decomposed case
406  if (decomposeFieldsOnly && nProcs != nDomains)
407  {
409  << "Specified -fields, but the case was decomposed with "
410  << nProcs << " domains" << nl << "instead of " << nDomains
411  << " domains as specified in decomposeParDict" << nl
412  << exit(FatalError);
413  }
414  }
415 
416  // Create meshes
417  multiDomainDecomposition regionMeshes(runTimes, meshPath, regionNames);
418  if
419  (
420  !(decomposeFieldsOnly && copyZero)
421  && regionMeshes.readDecompose(decomposeSets)
422  )
423  {
424  Info<< endl;
425 
426  if (writeCellProc)
427  {
428  forAll(regionNames, regioni)
429  {
430  writeDecomposition(regionMeshes[regioni]());
431  Info<< endl;
432  fileHandler().flush();
433  }
434  }
435  }
436 
437  // Get flag to determine whether or not to distribute uniform data
438  const bool distributed =
439  decomposeParDict.lookupOrDefault<bool>("distributed", false);
440 
441  // Loop over all times
442  forAll(times, timei)
443  {
444  // Set the time
445  runTimes.setTime(times[timei], timei);
446 
447  Info<< "Time = " << runTimes.completeTime().userTimeName()
448  << nl << endl;
449 
450  // Update the meshes, if necessary
451  const fvMesh::readUpdateState stat =
452  !(decomposeFieldsOnly && copyZero)
453  ? regionMeshes.readUpdateDecompose()
455  if (stat >= fvMesh::TOPO_CHANGE) Info<< endl;
456 
457  // Write the mesh out (if anything has changed), if necessary
458  if (!decomposeFieldsOnly)
459  {
460  regionMeshes.writeProcs(decomposeSets);
461  }
462 
463  // Write the decomposition, if necessary
464  forAll(regionNames, regioni)
465  {
466  if (writeCellProc && stat >= fvMesh::TOPO_CHANGE)
467  {
468  writeDecomposition(regionMeshes[regioni]());
469  Info<< endl;
470  fileHandler().flush();
471  }
472  }
473 
474  // If only decomposing geometry then there is no more to do
475  if (decomposeGeomOnly)
476  {
477  continue;
478  }
479 
480  // If copying from zero then just copy everything from the <time>
481  // directory to the processor*/<time> directories without altering them
482  if (copyZero)
483  {
484  const fileName completeTimePath =
485  runTimes.completeTime().timePath();
486 
487  fileName prevProcTimePath;
488  for (label proci = 0; proci < runTimes.nProcs(); proci++)
489  {
490  const Time& procRunTime = runTimes.procTimes()[proci];
491 
492  if (fileHandler().isDir(completeTimePath))
493  {
494  const fileName procTimePath
495  (
496  fileHandler().objectPath
497  (
498  IOobject(word::null, fileName::null, procRunTime)
499  )
500  );
501 
502  if (procTimePath != prevProcTimePath)
503  {
504  Info<< "Processor " << proci
505  << ": copying " << completeTimePath << nl
506  << " to " << procTimePath << endl;
507  fileHandler().cp(completeTimePath, procTimePath);
508  prevProcTimePath = procTimePath;
509  }
510  }
511  }
512 
513  Info<< endl;
514 
515  continue;
516  }
517 
518  // Otherwise we are doing full region-by-region decomposition of all
519  // the available fields
520  forAll(regionNames, regioni)
521  {
522  const word& regionName = regionNames[regioni];
523  const word regionDir =
525 
526  const delayedNewLine dnl;
527 
528  // Prefixed scope
529  {
530  const RegionRef<domainDecomposition> meshes =
531  regionMeshes[regioni];
532 
533  // Search for objects at this time
535  (
536  meshes().completeMesh(),
537  runTimes.completeTime().name()
538  );
539 
540  {
541  Info<< dnl << "Decomposing FV fields" << endl;
542 
544  {
545  fvFieldDecomposer fvDecomposer
546  (
547  meshes().completeMesh(),
548  meshes().procMeshes(),
549  meshes().procFaceAddressing(),
550  meshes().procCellAddressing(),
551  meshes().procFaceAddressingBf()
552  );
553 
554  #define DO_FV_VOL_INTERNAL_FIELDS_TYPE(Type, nullArg) \
555  fvDecomposer.decomposeVolInternalFields<Type> \
556  (objects);
557  FOR_ALL_FIELD_TYPES(DO_FV_VOL_INTERNAL_FIELDS_TYPE)
558  #undef DO_FV_VOL_INTERNAL_FIELDS_TYPE
559 
560  #define DO_FV_VOL_FIELDS_TYPE(Type, nullArg) \
561  fvDecomposer.decomposeVolFields<Type> \
562  (objects);
563  FOR_ALL_FIELD_TYPES(DO_FV_VOL_FIELDS_TYPE)
564  #undef DO_FV_VOL_FIELDS_TYPE
565 
566  #define DO_FV_SURFACE_FIELDS_TYPE(Type, nullArg) \
567  fvDecomposer.decomposeFvSurfaceFields<Type> \
568  (objects);
569  FOR_ALL_FIELD_TYPES(DO_FV_SURFACE_FIELDS_TYPE)
570  #undef DO_FV_SURFACE_FIELDS_TYPE
571  }
572  else
573  {
574  Info<< dnl << " (no FV fields)" << endl;
575  }
576  }
577 
578  {
579  Info<< dnl << "Decomposing point fields" << endl;
580 
582  {
583  pointFieldDecomposer pointDecomposer
584  (
585  pointMesh::New(meshes().completeMesh()),
586  meshes().procMeshes(),
587  meshes().procPointAddressing()
588  );
589 
590  #define DO_POINT_FIELDS_TYPE(Type, nullArg) \
591  pointDecomposer.decomposeFields<Type> \
592  (objects);
594  #undef DO_POINT_FIELDS_TYPE
595  }
596  else
597  {
598  Info<< dnl << " (no point fields)" << endl;
599  }
600  }
601 
602  {
603  // Find cloud directories
604  fileNameList cloudDirs
605  (
607  (
608  runTimes.completeTime().timePath()
609  /regionDir
612  )
613  );
614 
615  // Add objects in any found cloud directories
616  HashTable<IOobjectList> cloudsObjects;
617  forAll(cloudDirs, i)
618  {
619  // Do local scan for valid cloud objects
620  IOobjectList cloudObjs
621  (
622  meshes().completeMesh(),
623  runTimes.completeTime().name(),
624  lagrangian::cloud::prefix/cloudDirs[i],
627  false
628  );
629 
630  // If "positions" is present, then add to the table
631  if (cloudObjs.lookup(word("positions")))
632  {
633  cloudsObjects.insert(cloudDirs[i], cloudObjs);
634  }
635  }
636 
637  // Decompose the objects found above
638  if (cloudsObjects.size())
639  {
641  (
643  cloudsObjects,
644  iter
645  )
646  {
647  const word cloudName =
648  string::validate<word>(iter.key());
649 
650  const IOobjectList& cloudObjects = iter();
651 
652  Info<< dnl << "Decomposing lagrangian fields for "
653  << "cloud " << cloudName << endl;
654 
655  if
656  (
658  (
659  cloudObjects
660  )
661  )
662  {
664  lagrangianDecomposer
665  (
666  meshes().completeMesh(),
667  meshes().procMeshes(),
668  meshes().procFaceAddressing(),
669  meshes().procCellAddressing(),
670  cloudName
671  );
672 
673  lagrangianDecomposer.decomposePositions();
674 
675  #define DO_CLOUD_FIELDS_TYPE(Type, nullArg) \
676  lagrangianDecomposer.decomposeFields<Type> \
677  (cloudObjects);
678  DO_CLOUD_FIELDS_TYPE(label, )
679  FOR_ALL_FIELD_TYPES(DO_CLOUD_FIELDS_TYPE)
680  #undef DO_CLOUD_FIELDS_TYPE
681  }
682  else
683  {
684  Info<< dnl << " (no lagrangian fields)"
685  << endl;
686  }
687  }
688  }
689  }
690 
691  {
692  // Find Lagrangian directories
693  fileNameList LagrangianDirs
694  (
696  (
697  runTimes.completeTime().timePath()
698  /regionDir
701  )
702  );
703 
704  // Add objects in any found Lagrangian directories
705  HashTable<IOobjectList> LagrangianObjects;
706  forAll(LagrangianDirs, i)
707  {
708  // Do local scan for valid Lagrangian objects
710  (
711  meshes().completeMesh(),
712  runTimes.completeTime().name(),
713  LagrangianMesh::prefix/LagrangianDirs[i],
716  false
717  );
718 
719  // If coordinates or fields are present then add
720  // this set of objects to the table
721  if
722  (
724  || LagrangianFieldDecomposer::decomposes
725  (
726  objects
727  )
728  )
729  {
730  LagrangianObjects.insert
731  (
732  LagrangianDirs[i],
733  objects
734  );
735  }
736  }
737 
738  // Decompose the objects found above
739  if (LagrangianObjects.size())
740  {
742  (
744  LagrangianObjects,
745  iter
746  )
747  {
748  const word LagrangianName =
749  string::validate<word>(iter.key());
750 
751  Info<< dnl << "Decomposing Lagrangian fields "
752  << "for " << LagrangianName << endl;
753 
755  LagrangianDecomposer
756  (
757  meshes().completeMesh(),
758  meshes().procMeshes(),
759  meshes().procFaceAddressing(),
760  meshes().procCellAddressing(),
761  LagrangianName
762  );
763 
764  LagrangianDecomposer.decomposePositions();
765 
766  if
767  (
768  LagrangianFieldDecomposer::decomposes
769  (
770  iter()
771  )
772  )
773  {
774  #define DO_LAGRANGIAN_FIELDS_TYPE( \
775  Type, GeoField) \
776  LagrangianDecomposer \
777  .decomposeFields<GeoField<Type>> \
778  (iter());
780  (
781  label,
783  )
785  (
788  )
790  (
791  label,
793  )
795  (
798  )
799  #undef DO_LAGRANGIAN_FIELDS_TYPE
800 
801  // --> Note we don't have to explicitly
802  // decompose the dynamic variants of these
803  // fields as they are IO compatible with the
804  // non-dynamic fields
805  }
806  else
807  {
808  Info<< dnl << " (no lagrangian fields)"
809  << endl;
810  }
811  }
812  }
813  }
814  }
815 
816  Info<< dnl;
817  }
818 
819  // Distribute the uniform directory
820  if (haveUniform(runTimes))
821  {
822  Info<< "Distributing uniform files" << endl;
823 
824  decomposeUniform
825  (
826  copyUniform || distributed,
827  runTimes
828  );
829 
830  Info<< endl;
831  }
832 
833  if (regionNames == wordList(1, polyMesh::defaultRegion)) continue;
834 
835  // Distribute the region uniform directories
836  forAll(regionNames, regioni)
837  {
838  const word& regionName = regionNames[regioni];
839  const word regionDir =
841 
842  if (haveUniform(runTimes, regionDir))
843  {
844  // Prefixed scope
845  {
846  const RegionRef<domainDecomposition> meshes =
847  regionMeshes[regioni];
848 
849  Info<< "Distributing uniform files" << endl;
850 
851  decomposeUniform
852  (
853  copyUniform || distributed,
854  runTimes,
855  regionDir
856  );
857  }
858 
859  Info<< endl;
860  }
861  }
862  }
863 
864  Info<< "End" << nl << endl;
865 
866  return 0;
867 }
868 
869 
870 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:445
#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.
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
static const word coordinatesName
Name of the coordinates field.
static const word prefix
Instance prefix.
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: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
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 IOdictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
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
static const fileName null
An empty fileName.
Definition: fileName.H:97
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 fileName filePath(const bool globalFile, const IOobject &) const =0
Search for an object. globalFile : also check undecomposed case.
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.
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:420
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:105
@ TOPO_CHANGE
Definition: fvMesh.H:109
static bool decomposes(const IOobjectList &objects)
Return whether anything in the object list gets decomposed.
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:69
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:270
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
#define DO_LAGRANGIAN_FIELDS_TYPE(Type, nullArg)
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:258
const dimensionSet dimless
const word & regionName(const solver &region)
Definition: solver.H:218
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: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
bool chDir(const fileName &dir)
Change the current directory to the one given and return true,.
Definition: POSIX.C:284
#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)
const word cloudName(propsDict.lookup("cloudName"))