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