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-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  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"
37 #include "domainDecomposition.H"
38 #include "fvFieldReconstructor.H"
40 #include "reconstructLagrangian.H"
41 
42 using namespace Foam;
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 bool haveAllTimes
50 (
51  const HashSet<word>& masterTimeDirSet,
52  const instantList& timeDirs
53 )
54 {
55  // Loop over all times
56  forAll(timeDirs, timei)
57  {
58  if (!masterTimeDirSet.found(timeDirs[timei].name()))
59  {
60  return false;
61  }
62  }
63  return true;
64 }
65 
66 
67 void writeDecomposition(const domainDecomposition& meshes)
68 {
69  // Write as volScalarField::Internal for postprocessing.
71  (
72  IOobject
73  (
74  "cellProc",
75  meshes.completeMesh().time().name(),
76  meshes.completeMesh(),
79  ),
80  meshes.completeMesh(),
81  dimless,
83  );
84 
85  cellProc.write();
86 
87  Info<< "Wrote decomposition as volScalarField::Internal to "
88  << cellProc.name() << " for use in postprocessing."
89  << nl << endl;
90 }
91 
92 
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 
98 int main(int argc, char *argv[])
99 {
101  (
102  "Reconstruct fields of a parallel case"
103  );
104 
105  // Enable -constant ... if someone really wants it
106  // Enable -withZero to prevent accidentally trashing the initial fields
107  timeSelector::addOptions(true, true);
109  #include "addRegionOption.H"
110  #include "addAllRegionsOption.H"
112  (
113  "cellProc",
114  "write cell processor indices as a volScalarField::Internal for "
115  "post-processing."
116  );
118  (
119  "fields",
120  "list",
121  "specify a list of fields to be reconstructed. Eg, '(U T p)' - "
122  "regular expressions not currently supported"
123  );
125  (
126  "noFields",
127  "skip reconstructing fields"
128  );
130  (
131  "lagrangianFields",
132  "list",
133  "specify a list of lagrangian fields to be reconstructed. Eg, '(U d)' -"
134  "regular expressions not currently supported, "
135  "positions always included."
136  );
138  (
139  "noLagrangian",
140  "skip reconstructing lagrangian positions and fields"
141  );
143  (
144  "noSets",
145  "skip reconstructing cellSets, faceSets, pointSets"
146  );
148  (
149  "newTimes",
150  "only reconstruct new times (i.e. that do not exist already)"
151  );
152 
153  #include "setRootCase.H"
154 
155  const bool writeCellProc = args.optionFound("cellProc");
156 
157  HashSet<word> selectedFields;
158  if (args.optionFound("fields"))
159  {
160  args.optionLookup("fields")() >> selectedFields;
161  }
162 
163  const bool noFields = args.optionFound("noFields");
164 
165  if (noFields)
166  {
167  Info<< "Skipping reconstructing fields"
168  << nl << endl;
169  }
170 
171  const bool noLagrangian = args.optionFound("noLagrangian");
172 
173  if (noLagrangian)
174  {
175  Info<< "Skipping reconstructing lagrangian positions and fields"
176  << nl << endl;
177  }
178 
179  const bool noReconstructSets = args.optionFound("noSets");
180 
181  if (noReconstructSets)
182  {
183  Info<< "Skipping reconstructing cellSets, faceSets and pointSets"
184  << nl << endl;
185  }
186 
187  HashSet<word> selectedLagrangianFields;
188  if (args.optionFound("lagrangianFields"))
189  {
190  if (noLagrangian)
191  {
193  << "Cannot specify noLagrangian and lagrangianFields "
194  << "options together."
195  << exit(FatalError);
196  }
197 
198  args.optionLookup("lagrangianFields")() >> selectedLagrangianFields;
199  }
200 
201  // Set time from database
202  Info<< "Create time\n" << endl;
204 
205  // Allow override of time
206  const instantList times = runTimes.selectProc(args);
207 
208  const Time& runTime = runTimes.procTimes()[0];
209 
210  #include "setRegionNames.H"
211 
212  // Determine the processor count
213  const label nProcs = fileHandler().nProcs
214  (
215  args.path(),
217  ? word::null
218  : regionNames[0]
219  );
220 
221  if (!nProcs)
222  {
224  << "No processor* directories found"
225  << exit(FatalError);
226  }
227 
228  // Warn fileHandler of number of processors
229  const_cast<fileOperation&>(fileHandler()).setNProcs(nProcs);
230 
231  // Note that we do not set the runTime time so it is still the
232  // one set through the controlDict. The -time option
233  // only affects the selected set of times from processor0.
234  // - can be illogical
235  // + any point motion handled through mesh.readUpdate
236  if (times.empty())
237  {
238  WarningInFunction << "No times selected" << endl;
239  exit(1);
240  }
241 
242  // Get current times if -newTimes
243  const bool newTimes = args.optionFound("newTimes");
244  instantList masterTimeDirs;
245  if (newTimes)
246  {
247  masterTimeDirs = runTimes.completeTime().times();
248  }
249  HashSet<word> masterTimeDirSet(2*masterTimeDirs.size());
250  forAll(masterTimeDirs, i)
251  {
252  masterTimeDirSet.insert(masterTimeDirs[i].name());
253  }
254  if
255  (
256  newTimes
257  && regionNames.size() == 1
259  && haveAllTimes(masterTimeDirSet, times)
260  )
261  {
262  Info<< "All times already reconstructed.\n\nEnd\n" << endl;
263  return 0;
264  }
265 
266  // Reconstruct all regions
267  forAll(regionNames, regioni)
268  {
269  const word& regionName = regionNames[regioni];
270 
271  const word& regionDir =
273  ? word::null
274  : regionName;
275 
276  // Create meshes
277  Info<< "\n\nReconstructing mesh " << regionName << nl << endl;
278  domainDecomposition meshes(runTimes, regionName);
279  if (meshes.readReconstruct(!noReconstructSets) && writeCellProc)
280  {
281  writeDecomposition(meshes);
282  fileHandler().flush();
283  }
284 
285  // Loop over all times
286  forAll(times, timei)
287  {
288  if (newTimes && masterTimeDirSet.found(times[timei].name()))
289  {
290  Info<< "Skipping time " << times[timei].name()
291  << endl << endl;
292  continue;
293  }
294 
295  // Set the time
296  runTimes.setTime(times[timei], timei);
297 
298  Info<< "Time = " << runTimes.completeTime().userTimeName()
299  << nl << endl;
300 
301  // Update the meshes
302  const fvMesh::readUpdateState state =
303  meshes.readUpdateReconstruct();
304 
305  // Write the mesh out, if necessary
306  if (state != fvMesh::UNCHANGED)
307  {
308  meshes.writeComplete(!noReconstructSets);
309  }
310 
311  // Write the decomposition, if necessary
312  if
313  (
314  writeCellProc
315  && meshes.completeMesh().facesInstance()
316  == runTimes.completeTime().name()
317  )
318  {
319  writeDecomposition(meshes);
320  fileHandler().flush();
321  }
322 
323  // Get list of objects from processor0 database
325  (
326  meshes.procMeshes()[0],
327  runTimes.procTimes()[0].name()
328  );
329 
330  if (!noFields)
331  {
332  // If there are any FV fields, reconstruct them
333  Info<< "Reconstructing FV fields" << nl << endl;
334 
335  fvFieldReconstructor fvReconstructor
336  (
337  meshes.completeMesh(),
338  meshes.procMeshes(),
339  meshes.procFaceAddressing(),
340  meshes.procCellAddressing(),
341  meshes.procFaceAddressingBf()
342  );
343 
344  fvReconstructor.reconstructFvVolumeInternalFields<scalar>
345  (
346  objects,
347  selectedFields
348  );
349  fvReconstructor.reconstructFvVolumeInternalFields<vector>
350  (
351  objects,
352  selectedFields
353  );
354  fvReconstructor.reconstructFvVolumeInternalFields
356  (
357  objects,
358  selectedFields
359  );
360  fvReconstructor.reconstructFvVolumeInternalFields<symmTensor>
361  (
362  objects,
363  selectedFields
364  );
365  fvReconstructor.reconstructFvVolumeInternalFields<tensor>
366  (
367  objects,
368  selectedFields
369  );
370 
371  fvReconstructor.reconstructFvVolumeFields<scalar>
372  (
373  objects,
374  selectedFields
375  );
376  fvReconstructor.reconstructFvVolumeFields<vector>
377  (
378  objects,
379  selectedFields
380  );
381  fvReconstructor.reconstructFvVolumeFields<sphericalTensor>
382  (
383  objects,
384  selectedFields
385  );
386  fvReconstructor.reconstructFvVolumeFields<symmTensor>
387  (
388  objects,
389  selectedFields
390  );
391  fvReconstructor.reconstructFvVolumeFields<tensor>
392  (
393  objects,
394  selectedFields
395  );
396 
397  fvReconstructor.reconstructFvSurfaceFields<scalar>
398  (
399  objects,
400  selectedFields
401  );
402  fvReconstructor.reconstructFvSurfaceFields<vector>
403  (
404  objects,
405  selectedFields
406  );
407  fvReconstructor.reconstructFvSurfaceFields<sphericalTensor>
408  (
409  objects,
410  selectedFields
411  );
412  fvReconstructor.reconstructFvSurfaceFields<symmTensor>
413  (
414  objects,
415  selectedFields
416  );
417  fvReconstructor.reconstructFvSurfaceFields<tensor>
418  (
419  objects,
420  selectedFields
421  );
422 
423  if (fvReconstructor.nReconstructed() == 0)
424  {
425  Info<< "No FV fields" << nl << endl;
426  }
427  }
428 
429  if (!noFields)
430  {
431  Info<< "Reconstructing point fields" << nl << endl;
432 
433  const pointMesh& completePMesh =
434  pointMesh::New(meshes.completeMesh());
435 
436  pointFieldReconstructor pointReconstructor
437  (
438  completePMesh,
439  meshes.procMeshes(),
440  meshes.procPointAddressing()
441  );
442 
443  pointReconstructor.reconstructFields<scalar>
444  (
445  objects,
446  selectedFields
447  );
448  pointReconstructor.reconstructFields<vector>
449  (
450  objects,
451  selectedFields
452  );
453  pointReconstructor.reconstructFields<sphericalTensor>
454  (
455  objects,
456  selectedFields
457  );
458  pointReconstructor.reconstructFields<symmTensor>
459  (
460  objects,
461  selectedFields
462  );
463  pointReconstructor.reconstructFields<tensor>
464  (
465  objects,
466  selectedFields
467  );
468 
469  if (pointReconstructor.nReconstructed() == 0)
470  {
471  Info<< "No point fields" << nl << endl;
472  }
473  }
474 
475 
476  // If there are any clouds, reconstruct them.
477  // The problem is that a cloud of size zero will not get written so
478  // in pass 1 we determine the cloud names and per cloud name the
479  // fields. Note that the fields are stored as IOobjectList from
480  // the first processor that has them. They are in pass2 only used
481  // for name and type (scalar, vector etc).
482 
483  if (!noLagrangian)
484  {
485  HashTable<IOobjectList> cloudObjects;
486 
487  forAll(runTimes.procTimes(), proci)
488  {
489  fileName lagrangianDir
490  (
491  fileHandler().filePath
492  (
493  runTimes.procTimes()[proci].timePath()
494  /regionDir
496  )
497  );
498 
499  fileNameList cloudDirs;
500  if (!lagrangianDir.empty())
501  {
502  cloudDirs = fileHandler().readDir
503  (
504  lagrangianDir,
506  );
507  }
508 
509  forAll(cloudDirs, i)
510  {
511  // Check if we already have cloud objects for this
512  // cloudname
514  cloudObjects.find(cloudDirs[i]);
515 
516  if (iter == cloudObjects.end())
517  {
518  // Do local scan for valid cloud objects
519  IOobjectList sprayObjs
520  (
521  meshes.procMeshes()[proci],
522  runTimes.procTimes()[proci].name(),
523  cloud::prefix/cloudDirs[i]
524  );
525 
526  IOobject* positionsPtr =
527  sprayObjs.lookup(word("positions"));
528 
529  if (positionsPtr)
530  {
531  cloudObjects.insert(cloudDirs[i], sprayObjs);
532  }
533  }
534  }
535  }
536 
537  if (cloudObjects.size())
538  {
539  // Pass2: reconstruct the cloud
540  forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
541  {
542  const word cloudName =
543  string::validate<word>(iter.key());
544 
545  // Objects (on arbitrary processor)
546  const IOobjectList& sprayObjs = iter();
547 
548  Info<< "Reconstructing lagrangian fields for cloud "
549  << cloudName << nl << endl;
550 
552  (
553  meshes.completeMesh(),
554  cloudName,
555  meshes.procMeshes(),
556  meshes.procFaceAddressing(),
557  meshes.procCellAddressing()
558  );
559  reconstructLagrangianFields<label>
560  (
561  cloudName,
562  meshes.completeMesh(),
563  meshes.procMeshes(),
564  sprayObjs,
565  selectedLagrangianFields
566  );
567  reconstructLagrangianFieldFields<label>
568  (
569  cloudName,
570  meshes.completeMesh(),
571  meshes.procMeshes(),
572  sprayObjs,
573  selectedLagrangianFields
574  );
575  reconstructLagrangianFields<scalar>
576  (
577  cloudName,
578  meshes.completeMesh(),
579  meshes.procMeshes(),
580  sprayObjs,
581  selectedLagrangianFields
582  );
583  reconstructLagrangianFieldFields<scalar>
584  (
585  cloudName,
586  meshes.completeMesh(),
587  meshes.procMeshes(),
588  sprayObjs,
589  selectedLagrangianFields
590  );
591  reconstructLagrangianFields<vector>
592  (
593  cloudName,
594  meshes.completeMesh(),
595  meshes.procMeshes(),
596  sprayObjs,
597  selectedLagrangianFields
598  );
599  reconstructLagrangianFieldFields<vector>
600  (
601  cloudName,
602  meshes.completeMesh(),
603  meshes.procMeshes(),
604  sprayObjs,
605  selectedLagrangianFields
606  );
607  reconstructLagrangianFields<sphericalTensor>
608  (
609  cloudName,
610  meshes.completeMesh(),
611  meshes.procMeshes(),
612  sprayObjs,
613  selectedLagrangianFields
614  );
615  reconstructLagrangianFieldFields<sphericalTensor>
616  (
617  cloudName,
618  meshes.completeMesh(),
619  meshes.procMeshes(),
620  sprayObjs,
621  selectedLagrangianFields
622  );
623  reconstructLagrangianFields<symmTensor>
624  (
625  cloudName,
626  meshes.completeMesh(),
627  meshes.procMeshes(),
628  sprayObjs,
629  selectedLagrangianFields
630  );
631  reconstructLagrangianFieldFields<symmTensor>
632  (
633  cloudName,
634  meshes.completeMesh(),
635  meshes.procMeshes(),
636  sprayObjs,
637  selectedLagrangianFields
638  );
639  reconstructLagrangianFields<tensor>
640  (
641  cloudName,
642  meshes.completeMesh(),
643  meshes.procMeshes(),
644  sprayObjs,
645  selectedLagrangianFields
646  );
647  reconstructLagrangianFieldFields<tensor>
648  (
649  cloudName,
650  meshes.completeMesh(),
651  meshes.procMeshes(),
652  sprayObjs,
653  selectedLagrangianFields
654  );
655  }
656  }
657  else
658  {
659  Info<< "No lagrangian fields" << nl << endl;
660  }
661  }
662 
663  // If there is a "uniform" directory in the time region
664  // directory copy from the master processor
665  {
666  fileName uniformDir0
667  (
668  fileHandler().filePath
669  (
670  runTimes.procTimes()[0].timePath()/regionDir/"uniform"
671  )
672  );
673 
674  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
675  {
676  fileHandler().cp
677  (
678  uniformDir0,
679  runTimes.completeTime().timePath()/regionDir
680  );
681  }
682  }
683 
684  // For the first region of a multi-region case additionally
685  // copy the "uniform" directory in the time directory
686  if (regioni == 0 && regionDir != word::null)
687  {
688  fileName uniformDir0
689  (
690  fileHandler().filePath
691  (
692  runTimes.procTimes()[0].timePath()/"uniform"
693  )
694  );
695 
696  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
697  {
698  fileHandler().cp
699  (
700  uniformDir0,
701  runTimes.completeTime().timePath()
702  );
703  }
704  }
705  }
706  }
707 
708  Info<< "\nEnd\n" << endl;
709 
710  return 0;
711 }
712 
713 
714 // ************************************************************************* //
#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
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:111
An STL-conforming const_iterator.
Definition: HashTable.H:484
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:142
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
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 size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:204
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 labelListList & procPointAddressing() const
Access the labels of points for each processor.
const PtrList< fvMesh > & procMeshes() const
Access the processor meshes.
void writeComplete(const bool doSets) const
Write the decomposed meshes and associated data.
bool readReconstruct(const bool doSets)
Read in the processor meshes. Read the complete mesh if it is.
const labelListList & procFaceAddressing() const
Access the labels of faces for each processor (see notes above)
const PtrList< surfaceLabelField::Boundary > & procFaceAddressingBf() const
Access the labels of faces for each processor (see notes above)
const labelList & cellProc() const
Return the distribution as an FV field for writing.
const labelListList & procCellAddressing() const
Access the labels of cells for each processor.
const fvMesh & completeMesh() const
Access the global mesh.
fvMesh::readUpdateState readUpdateReconstruct()
Read-update for reconstruction.
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
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 fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true) const =0
Read a directory and return the entries as a string list.
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
Finite volume reconstructor for volume and surface fields.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:402
Point field reconstructor.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1025
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
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
Foam::word regionName
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
static instantList timeDirs
Definition: globalFoam.H: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
void reconstructLagrangianPositions(const polyMesh &mesh, const word &cloudName, const PtrList< fvMesh > &meshes, const labelListList &faceProcAddressing, const labelListList &cellProcAddressing)
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:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
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
error FatalError
static const char nl
Definition: Ostream.H:260
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"))