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-2024 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"
41 
42 using namespace Foam;
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 bool haveUniform
50 (
51  const processorRunTimes& runTimes,
52  const word& regionDir = word::null
53 )
54 {
55  return
57  (
58  fileHandler().filePath
59  (
60  runTimes.procTimes()[0].timePath()/regionDir/"uniform"
61  )
62  );
63 }
64 
65 
66 void reconstructUniform
67 (
68  const processorRunTimes& runTimes,
69  const word& regionDir = word::null
70 )
71 {
72  fileHandler().cp
73  (
74  fileHandler().filePath
75  (
76  runTimes.procTimes()[0].timePath()/regionDir/"uniform"
77  ),
78  runTimes.completeTime().timePath()/regionDir
79  );
80 }
81 
82 
83 void writeDecomposition(const domainDecomposition& meshes)
84 {
85  // Write as volScalarField::Internal for postprocessing.
87  (
88  IOobject
89  (
90  "cellProc",
91  meshes.completeMesh().time().name(),
92  meshes.completeMesh(),
95  ),
96  meshes.completeMesh(),
97  dimless,
99  );
100 
101  cellProc.write();
102 
103  Info<< "Wrote decomposition as volScalarField::Internal to "
104  << cellProc.name() << " for use in postprocessing"
105  << endl;
106 }
107 
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 
113 namespace Foam
114 {
115 
116 class delayedNewLine
117 {
118  mutable bool first_;
119 
120 public:
121 
122  delayedNewLine()
123  :
124  first_(true)
125  {}
126 
127  friend Ostream& operator<<(Ostream& os, const delayedNewLine& dnl)
128  {
129  if (!dnl.first_) os << nl;
130  dnl.first_ = false;
131  return os;
132  }
133 };
134 
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 int main(int argc, char *argv[])
141 {
143  (
144  "Reconstruct fields of a parallel case"
145  );
146 
148  #include "addRegionOption.H"
149  #include "addAllRegionsOption.H"
151  (
152  "cellProc",
153  "write cell processor indices as a volScalarField::Internal for "
154  "post-processing"
155  );
157  (
158  "fields",
159  "list",
160  "specify a list of fields to be reconstructed. Eg, '(U T p)' - "
161  "regular expressions not currently supported"
162  );
164  (
165  "noFields",
166  "skip reconstructing fields"
167  );
169  (
170  "lagrangianFields",
171  "list",
172  "specify a list of lagrangian fields to be reconstructed. Eg, '(U d)' -"
173  "regular expressions not currently supported, "
174  "positions always included"
175  );
177  (
178  "noLagrangian",
179  "skip reconstructing lagrangian positions and fields"
180  );
182  (
183  "noSets",
184  "skip reconstructing cellSets, faceSets, pointSets"
185  );
187  (
188  "newTimes",
189  "only reconstruct new times (i.e. that do not exist already)"
190  );
192  (
193  "rm",
194  "remove processor time directories after reconstruction"
195  );
196 
197  // Include explicit constant options, and explicit zero option (to prevent
198  // the user accidentally trashing the initial fields)
199  timeSelector::addOptions(true, true);
200 
201  #include "setRootCase.H"
202 
203  const bool writeCellProc = args.optionFound("cellProc");
204 
205  HashSet<word> selectedFields;
206  if (args.optionFound("fields"))
207  {
208  args.optionLookup("fields")() >> selectedFields;
209  }
210 
211  const bool noFields = args.optionFound("noFields");
212 
213  if (noFields)
214  {
215  Info<< "Skipping reconstructing fields" << nl << endl;
216  }
217 
218  const bool noLagrangian = args.optionFound("noLagrangian");
219 
220  if (noLagrangian)
221  {
222  Info<< "Skipping reconstructing lagrangian positions and fields"
223  << nl << endl;
224  }
225 
226  const bool noReconstructSets = args.optionFound("noSets");
227 
228  if (noReconstructSets)
229  {
230  Info<< "Skipping reconstructing cellSets, faceSets and pointSets"
231  << nl << endl;
232  }
233 
234  HashSet<word> selectedLagrangianFields;
235  if (args.optionFound("lagrangianFields"))
236  {
237  if (noLagrangian)
238  {
240  << "Cannot specify noLagrangian and lagrangianFields "
241  << "options together"
242  << exit(FatalError);
243  }
244 
245  args.optionLookup("lagrangianFields")() >> selectedLagrangianFields;
246  }
247 
248  // Set time from database
249  Info<< "Create time" << nl << endl;
251 
252  // Get the times to reconstruct
253  instantList times = runTimes.selectProc(args);
254 
255  const Time& runTime = runTimes.procTimes()[0];
256 
257  #include "setRegionNames.H"
258 
259  // Determine the processor count
260  const label nProcs = fileHandler().nProcs
261  (
262  args.path(),
264  ? word::null
265  : regionNames[0]
266  );
267 
268  if (!nProcs)
269  {
271  << "No processor* directories found"
272  << exit(FatalError);
273  }
274 
275  // Warn fileHandler of number of processors
276  const_cast<fileOperation&>(fileHandler()).setNProcs(nProcs);
277 
278  // Quit if no times
279  if (times.empty())
280  {
281  WarningInFunction << "No times selected" << nl << endl;
282  exit(1);
283  }
284 
285  // If only reconstructing new times then filter out existing times
286  if (args.optionFound("newTimes"))
287  {
288  // Get all existing times
289  const instantList existingTimes = runTimes.completeTime().times();
290 
291  // Put into a set
292  HashSet<word> existingTimesSet;
293  existingTimesSet.resize(2*existingTimes.size());
294  forAll(existingTimes, i)
295  {
296  existingTimesSet.insert(existingTimes[i].name());
297  }
298 
299  // Remove times from the existing time set by shuffling up
300  label timei = 0;
301  forAll(times, timej)
302  {
303  if (!existingTimesSet.found(times[timej].name()))
304  {
305  times[timei ++] = times[timej];
306  }
307  }
308  times.resize(timei);
309  }
310 
311  // Quit if no times
312  if (times.empty())
313  {
314  Info<< "All times already reconstructed" << nl << nl
315  << "End" << nl << endl;
316  return 0;
317  }
318 
319  // Create meshes
320  multiDomainDecomposition regionMeshes(runTimes, regionNames);
321  if (regionMeshes.readReconstruct(!noReconstructSets))
322  {
323  Info<< endl;
324 
325  if (writeCellProc)
326  {
327  forAll(regionNames, regioni)
328  {
329  writeDecomposition(regionMeshes[regioni]());
330  Info<< endl;
331  fileHandler().flush();
332  }
333  }
334  }
335 
336  // Loop over all times
337  forAll(times, timei)
338  {
339  // Set the time
340  runTimes.setTime(times[timei], timei);
341 
342  Info<< "Time = " << runTimes.completeTime().userTimeName()
343  << nl << endl;
344 
345  // Update the meshes
346  const fvMesh::readUpdateState stat =
347  regionMeshes.readUpdateReconstruct();
348  if (stat >= fvMesh::TOPO_CHANGE) Info<< endl;
349 
350  // Write the mesh out (if anything has changed)
351  regionMeshes.writeComplete(!noReconstructSets);
352 
353  // Write the decomposition, if necessary
354  forAll(regionNames, regioni)
355  {
356  if (writeCellProc && stat >= fvMesh::TOPO_CHANGE)
357  {
358  writeDecomposition(regionMeshes[regioni]());
359  Info<< endl;
360  fileHandler().flush();
361  }
362  }
363 
364  // Do a region-by-region reconstruction of all the available fields
365  forAll(regionNames, regioni)
366  {
367  const word& regionName = regionNames[regioni];
368  const word regionDir =
370 
371  const delayedNewLine dnl;
372 
373  // Prefixed scope
374  {
375  const RegionRef<domainDecomposition> meshes =
376  regionMeshes[regioni];
377 
378  // Search for objects at this time
380  (
381  meshes().procMeshes()[0],
382  runTimes.procTimes()[0].name()
383  );
384 
385  if (!noFields)
386  {
387  Info<< dnl << "Reconstructing FV fields" << endl;
388 
389  if
390  (
391  fvFieldReconstructor::reconstructs
392  (
393  objects,
394  selectedFields
395  )
396  )
397  {
398  fvFieldReconstructor fvReconstructor
399  (
400  meshes().completeMesh(),
401  meshes().procMeshes(),
402  meshes().procFaceAddressing(),
403  meshes().procCellAddressing(),
404  meshes().procFaceAddressingBf()
405  );
406 
407  #define DO_FV_VOL_INTERNAL_FIELDS_TYPE(Type, nullArg) \
408  fvReconstructor.reconstructVolInternalFields<Type> \
409  (objects, selectedFields);
410  FOR_ALL_FIELD_TYPES(DO_FV_VOL_INTERNAL_FIELDS_TYPE)
411  #undef DO_FV_VOL_INTERNAL_FIELDS_TYPE
412 
413  #define DO_FV_VOL_FIELDS_TYPE(Type, nullArg) \
414  fvReconstructor.reconstructVolFields<Type> \
415  (objects, selectedFields);
416  FOR_ALL_FIELD_TYPES(DO_FV_VOL_FIELDS_TYPE)
417  #undef DO_FV_VOL_FIELDS_TYPE
418 
419  #define DO_FV_SURFACE_FIELDS_TYPE(Type, nullArg) \
420  fvReconstructor.reconstructFvSurfaceFields<Type> \
421  (objects, selectedFields);
422  FOR_ALL_FIELD_TYPES(DO_FV_SURFACE_FIELDS_TYPE)
423  #undef DO_FV_SURFACE_FIELDS_TYPE
424  }
425  else
426  {
427  Info<< dnl << " (no FV fields)" << endl;
428  }
429  }
430 
431  if (!noFields)
432  {
433  Info<< dnl << "Reconstructing point fields" << endl;
434 
435  if
436  (
437  pointFieldReconstructor::reconstructs
438  (
439  objects,
440  selectedFields
441  )
442  )
443  {
444  pointFieldReconstructor pointReconstructor
445  (
446  pointMesh::New(meshes().completeMesh()),
447  meshes().procMeshes(),
448  meshes().procPointAddressing()
449  );
450 
451  #define DO_POINT_FIELDS_TYPE(Type, nullArg) \
452  pointReconstructor.reconstructFields<Type> \
453  (objects, selectedFields);
454  FOR_ALL_FIELD_TYPES(DO_POINT_FIELDS_TYPE)
455  #undef DO_POINT_FIELDS_TYPE
456  }
457  else
458  {
459  Info<< dnl << " (no point fields)" << endl;
460  }
461  }
462 
463  if (!noLagrangian)
464  {
465  // Search for clouds that exist on any processor and add
466  // them into this table of cloud objects
467  HashTable<IOobjectList> cloudsObjects;
468  forAll(runTimes.procTimes(), proci)
469  {
470  // Find cloud directories
471  fileNameList cloudDirs
472  (
474  (
475  fileHandler().filePath
476  (
477  runTimes.procTimes()[proci].timePath()
478  /regionDir
480  ),
482  )
483  );
484 
485  // Add objects in any found cloud directories
486  forAll(cloudDirs, i)
487  {
488  // Pass if we already have an objects for this name
490  cloudsObjects.find(cloudDirs[i]);
491  if (iter != cloudsObjects.end()) continue;
492 
493  // Do local scan for valid cloud objects
494  IOobjectList cloudObjs
495  (
496  meshes().procMeshes()[proci],
497  runTimes.procTimes()[proci].name(),
498  cloud::prefix/cloudDirs[i],
501  false
502  );
503 
504  // If "positions" is present, then add to the table
505  if (cloudObjs.lookup(word("positions")))
506  {
507  cloudsObjects.insert(cloudDirs[i], cloudObjs);
508  }
509  }
510  }
511 
512  // Reconstruct the objects found above
513  if (cloudsObjects.size())
514  {
516  (
518  cloudsObjects,
519  iter
520  )
521  {
522  const word cloudName =
523  string::validate<word>(iter.key());
524 
525  const IOobjectList& cloudObjects = iter();
526 
527  Info<< dnl << "Reconstructing lagrangian fields "
528  << "for cloud " << cloudName << endl;
529 
530  if
531  (
532  lagrangianFieldReconstructor::reconstructs
533  (
534  cloudObjects,
535  selectedLagrangianFields
536  )
537  )
538  {
540  lagrangianReconstructor
541  (
542  meshes().completeMesh(),
543  meshes().procMeshes(),
544  meshes().procFaceAddressing(),
545  meshes().procCellAddressing(),
546  cloudName
547  );
548 
549  #define DO_CLOUD_FIELDS_TYPE(Type, nullArg) \
550  lagrangianReconstructor \
551  .reconstructFields<Type> \
552  (cloudObjects, selectedLagrangianFields);
553  DO_CLOUD_FIELDS_TYPE(label, );
554  FOR_ALL_FIELD_TYPES(DO_CLOUD_FIELDS_TYPE)
555  #undef DO_CLOUD_FIELDS_TYPE
556  }
557  else
558  {
559  Info<< dnl << " (no lagrangian fields)"
560  << endl;
561  }
562  }
563  }
564  }
565  }
566 
567  Info<< dnl;
568  }
569 
570  // Collect the uniform directory
571  if (haveUniform(runTimes))
572  {
573  Info<< "Collecting uniform files" << endl;
574 
575  reconstructUniform(runTimes);
576 
577  Info<< endl;
578  }
579 
581  {
582  // Collect the region uniform directories
583  forAll(regionNames, regioni)
584  {
585  const word& regionName = regionNames[regioni];
586  const word regionDir =
588  ? word::null
589  : regionName;
590 
591  if (haveUniform(runTimes, regionDir))
592  {
593  // Prefixed scope
594  {
595  const RegionRef<domainDecomposition> meshes =
596  regionMeshes[regioni];
597 
598  Info<< "Collecting uniform files" << endl;
599 
600  reconstructUniform(runTimes, regionDir);
601  }
602 
603  Info<< endl;
604  }
605  }
606  }
607 
608  if (args.optionFound("rm") && times[timei].name() != Time::constantName)
609  {
610  const bool allRegions = args.optionFound("allRegions");
611 
612  forAll(regionNames, regioni)
613  {
614  const word& regionName = regionNames[regioni];
615  const word regionDir =
616  allRegions || regionName == polyMesh::defaultRegion
617  ? word::null
618  : regionName;
619 
620  const RegionRef<domainDecomposition> meshes =
621  regionMeshes[regioni];
622 
623  Info<< "Removing processors time directory" << endl;
624 
625  for (label proci=0; proci<nProcs; proci++)
626  {
627  const fileName procTimePath = fileHandler().filePath
628  (
629  runTimes.procTimes()[proci].timePath()/regionDir
630  );
631 
632  if (isDir(procTimePath))
633  {
634  rmDir(procTimePath);
635  }
636  }
637  }
638 
639  Info<< endl;
640  }
641  }
642 
643  Info<< "End" << nl << endl;
644 
645  return 0;
646 }
647 
648 
649 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#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...
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
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:651
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
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
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
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 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 fileName filePath(const bool globalFile, const IOobject &, const word &typeName) const =0
Search for an object. globalFile : also check undecomposed case.
Finite volume reconstructor for volume and surface fields.
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
Point field reconstructor.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
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 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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const word & regionName(const solver &region)
Definition: solver.H:209
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
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
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)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
const word cloudName(propsDict.lookup("cloudName"))