reconstructPar.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  reconstructPar
29 
30 Group
31  grpParallelUtilities
32 
33 Description
34  Reconstructs fields of a case that is decomposed for parallel
35  execution of OpenFOAM.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "argList.H"
40 #include "timeSelector.H"
41 
42 #include "fvCFD.H"
43 #include "IOobjectList.H"
44 #include "processorMeshes.H"
45 #include "regionProperties.H"
46 #include "fvFieldReconstructor.H"
49 
50 #include "faCFD.H"
51 #include "faMesh.H"
52 #include "processorFaMeshes.H"
53 #include "faFieldReconstructor.H"
54 
55 #include "cellSet.H"
56 #include "faceSet.H"
57 #include "pointSet.H"
58 
59 #include "hexRef8Data.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 bool haveAllTimes
64 (
65  const wordHashSet& masterTimeDirSet,
66  const instantList& timeDirs
67 )
68 {
69  // Loop over all times
70  for (const instant& t : timeDirs)
71  {
72  if (!masterTimeDirSet.found(t.name()))
73  {
74  return false;
75  }
76  }
77  return true;
78 }
79 
80 
81 int main(int argc, char *argv[])
82 {
83  argList::addNote
84  (
85  "Reconstruct fields of a parallel case"
86  );
87 
88  // Enable -constant ... if someone really wants it
89  // Enable -withZero to prevent accidentally trashing the initial fields
90  timeSelector::addOptions(true, true); // constant(true), zero(true)
91  argList::noParallel();
92 
93  #include "addAllRegionOptions.H"
94 
95  argList::addVerboseOption();
96  argList::addOption
97  (
98  "fields",
99  "wordRes",
100  "Specify single or multiple fields to reconstruct (all by default)."
101  " Eg, 'T' or '(p T U \"alpha.*\")'"
102  );
103  argList::addBoolOption
104  (
105  "no-fields", // noFields
106  "Skip reconstructing fields"
107  );
108  argList::addOptionCompat("no-fields", {"noFields", 2106});
109  argList::addOption
110  (
111  "lagrangianFields",
112  "wordRes",
113  "Specify single or multiple lagrangian fields to reconstruct"
114  " (all by default)."
115  " Eg, '(U d)'"
116  " - Positions are always included."
117  );
118  argList::addBoolOption
119  (
120  "no-lagrangian", // noLagrangian
121  "Skip reconstructing lagrangian positions and fields"
122  );
123  argList::addOptionCompat("no-lagrangian", {"noLagrangian", 2106});
124 
125  argList::addBoolOption
126  (
127  "no-sets",
128  "Skip reconstructing cellSets, faceSets, pointSets"
129  );
130  argList::addOptionCompat("no-sets", {"noSets", 2106});
131 
132  argList::addBoolOption
133  (
134  "newTimes",
135  "Only reconstruct new times (i.e. that do not exist already)"
136  );
137 
138  #include "setRootCase.H"
139  #include "createTime.H"
140 
141 
142  const bool doFields = !args.found("no-fields");
143  wordRes selectedFields;
144 
145  if (doFields)
146  {
147  args.readListIfPresent<wordRe>("fields", selectedFields);
148  }
149  else
150  {
151  Info<< "Skipping reconstructing fields";
152  if (args.found("fields"))
153  {
154  Info<< ". Ignore -fields option";
155  }
156  Info<< nl << endl;
157  }
158 
159 
160  const bool doFiniteArea = !args.found("no-finite-area");
161  if (!doFiniteArea)
162  {
163  Info<< "Skipping reconstructing finiteArea mesh/fields"
164  << nl << endl;
165  }
166 
167 
168  const bool doLagrangian = !args.found("no-lagrangian");
169  wordRes selectedLagrangianFields;
170 
171  if (doLagrangian)
172  {
173  args.readListIfPresent<wordRe>
174  (
175  "lagrangianFields", selectedLagrangianFields
176  );
177  }
178  else
179  {
180  Info<< "Skipping reconstructing lagrangian positions/fields";
181  if (args.found("lagrangianFields"))
182  {
183  Info<< ". Ignore -lagrangianFields option";
184  }
185  Info<< nl << endl;
186  }
187 
188 
189  const bool doReconstructSets = !args.found("no-sets");
190 
191  if (!doReconstructSets)
192  {
193  Info<< "Skipping reconstructing cellSets, faceSets and pointSets"
194  << nl << endl;
195  }
196 
197  const bool newTimes = args.found("newTimes");
198 
199  // Get region names
200  #include "getAllRegionOptions.H"
201 
202  // Determine the processor count
203  label nProcs{0};
204 
205  if (regionNames.empty())
206  {
208  << "No regions specified or detected."
209  << exit(FatalError);
210  }
211  else if (regionNames[0] == polyMesh::defaultRegion)
212  {
213  nProcs = fileHandler().nProcs(args.path());
214  }
215  else
216  {
217  nProcs = fileHandler().nProcs(args.path(), regionNames[0]);
218 
219  if (regionNames.size() == 1)
220  {
221  Info<< "Using region: " << regionNames[0] << nl << endl;
222  }
223  }
224 
225  if (!nProcs)
226  {
228  << "No processor* directories found"
229  << exit(FatalError);
230  }
231 
232  // Warn fileHandler of number of processors
233  const_cast<fileOperation&>(fileHandler()).nProcs(nProcs);
234 
235  // Create the processor databases
236  PtrList<Time> databases(nProcs);
237 
238  forAll(databases, proci)
239  {
240  databases.set
241  (
242  proci,
243  new Time
244  (
245  Time::controlDictName,
246  args.rootPath(),
247  args.caseName()/("processor" + Foam::name(proci)),
249  args.allowLibs()
250  )
251  );
252  }
253 
254  // Use the times list from the master processor
255  // and select a subset based on the command-line options
257  (
258  databases[0].times(),
259  args
260  );
261 
262  // Note that we do not set the runTime time so it is still the
263  // one set through the controlDict. The -time option
264  // only affects the selected set of times from processor0.
265  // - can be illogical
266  // + any point motion handled through mesh.readUpdate
267 
268 
269  if (timeDirs.empty())
270  {
271  WarningInFunction << "No times selected";
272  exit(1);
273  }
274 
275 
276  // Get current times if -newTimes
277  instantList masterTimeDirs;
278  if (newTimes)
279  {
280  masterTimeDirs = runTime.times();
281  }
282  wordHashSet masterTimeDirSet(2*masterTimeDirs.size());
283  for (const instant& t : masterTimeDirs)
284  {
285  masterTimeDirSet.insert(t.name());
286  }
287 
288 
289  // Set all times on processor meshes equal to reconstructed mesh
290  forAll(databases, proci)
291  {
292  databases[proci].setTime(runTime);
293  }
294 
295 
296  forAll(regionNames, regioni)
297  {
298  const word& regionName = regionNames[regioni];
300 
301  Info<< "\n\nReconstructing fields" << nl
302  << "region=" << regionName << nl << endl;
303 
304  if
305  (
306  newTimes
307  && regionNames.size() == 1
308  && regionDir.empty()
309  && haveAllTimes(masterTimeDirSet, timeDirs)
310  )
311  {
312  Info<< "Skipping region " << regionName
313  << " since already have all times"
314  << endl << endl;
315  continue;
316  }
317 
318 
319  fvMesh mesh
320  (
321  IOobject
322  (
323  regionName,
324  runTime.timeName(),
325  runTime,
327  )
328  );
329 
330 
331  // Read all meshes and addressing to reconstructed mesh
332  processorMeshes procMeshes(databases, regionName);
333 
334  // Loop over all times
335  forAll(timeDirs, timei)
336  {
337  if (newTimes && masterTimeDirSet.found(timeDirs[timei].name()))
338  {
339  Info<< "Skipping time " << timeDirs[timei].name()
340  << endl << endl;
341  continue;
342  }
343 
344 
345  // Set time for global database
346  runTime.setTime(timeDirs[timei], timei);
347 
348  Info<< "Time = " << runTime.timeName() << endl << endl;
349 
350  // Set time for all databases
351  forAll(databases, proci)
352  {
353  databases[proci].setTime(timeDirs[timei], timei);
354  }
355 
356  // Check if any new meshes need to be read.
357  polyMesh::readUpdateState meshStat = mesh.readUpdate();
358 
359  polyMesh::readUpdateState procStat = procMeshes.readUpdate();
360 
361  if (procStat == polyMesh::POINTS_MOVED)
362  {
363  // Reconstruct the points for moving mesh cases and write
364  // them out
365  procMeshes.reconstructPoints(mesh);
366  }
367  else if (meshStat != procStat)
368  {
370  << "readUpdate for the reconstructed mesh:"
371  << meshStat << nl
372  << "readUpdate for the processor meshes :"
373  << procStat << nl
374  << "These should be equal or your addressing"
375  << " might be incorrect."
376  << " Please check your time directories for any "
377  << "mesh directories." << endl;
378  }
379 
380 
381  // Get list of objects from processor0 database
382  IOobjectList objects
383  (
384  procMeshes.meshes()[0],
385  databases[0].timeName(),
386  IOobjectOption::NO_REGISTER
387  );
388 
389  IOobjectList faObjects;
390 
391  if (doFiniteArea && doFields)
392  {
393  // List of area mesh objects (assuming single region)
394  // - scan on processor0
395  faObjects = IOobjectList
396  (
397  procMeshes.meshes()[0],
398  databases[0].timeName(),
399  faMesh::dbDir(word::null), // local relative to mesh
400  IOobjectOption::NO_REGISTER
401  );
402  }
403 
404  if (doFields)
405  {
406  // If there are any FV fields, reconstruct them
407  Info<< "Reconstructing FV fields" << nl << endl;
408 
409  fvFieldReconstructor reconstructor
410  (
411  mesh,
412  procMeshes.meshes(),
413  procMeshes.faceProcAddressing(),
414  procMeshes.cellProcAddressing(),
415  procMeshes.boundaryProcAddressing()
416  );
417 
418  reconstructor.reconstructAllFields(objects, selectedFields);
419 
420  if (reconstructor.nReconstructed() == 0)
421  {
422  Info<< "No FV fields" << nl << endl;
423  }
424  }
425 
426  if (doFields)
427  {
428  Info<< "Reconstructing point fields" << nl << endl;
429 
430  const pointMesh& pMesh = pointMesh::New(mesh);
431  PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
432 
433  forAll(pMeshes, proci)
434  {
435  pMeshes.set
436  (
437  proci,
438  new pointMesh(procMeshes.meshes()[proci])
439  );
440  }
441 
442  pointFieldReconstructor reconstructor
443  (
444  pMesh,
445  pMeshes,
446  procMeshes.pointProcAddressing(),
447  procMeshes.boundaryProcAddressing()
448  );
449 
450  reconstructor.reconstructAllFields(objects, selectedFields);
451 
452  if (reconstructor.nReconstructed() == 0)
453  {
454  Info<< "No point fields" << nl << endl;
455  }
456  }
457 
458 
459  // If there are any clouds, reconstruct them.
460  // The problem is that a cloud of size zero will not get written so
461  // in pass 1 we determine the cloud names and per cloud name the
462  // fields. Note that the fields are stored as IOobjectList from
463  // the first processor that has them. They are in pass2 only used
464  // for name and type (scalar, vector etc).
465 
466  if (doLagrangian)
467  {
468  HashTable<IOobjectList> allCloudObjects;
469 
470  forAll(databases, proci)
471  {
472  fileName lagrangianDir
473  (
474  fileHandler().filePath
475  (
476  databases[proci].timePath()
477  / regionDir
478  / cloud::prefix
479  )
480  );
481 
482  fileNameList cloudDirs;
483  if (!lagrangianDir.empty())
484  {
485  cloudDirs = fileHandler().readDir
486  (
487  lagrangianDir,
488  fileName::DIRECTORY
489  );
490  }
491 
492  for (const fileName& cloudDir : cloudDirs)
493  {
494  // Check if we already have cloud objects for this
495  // cloudname
496  if (!allCloudObjects.found(cloudDir))
497  {
498  // Do local scan for valid cloud objects
499  IOobjectList localObjs
500  (
501  procMeshes.meshes()[proci],
502  databases[proci].timeName(),
503  cloud::prefix/cloudDir
504  );
505 
506  if
507  (
508  localObjs.found("coordinates")
509  || localObjs.found("positions")
510  )
511  {
512  allCloudObjects.insert(cloudDir, localObjs);
513  }
514  }
515  }
516  }
517 
518 
519  if (allCloudObjects.size())
520  {
521  lagrangianReconstructor reconstructor
522  (
523  mesh,
524  procMeshes.meshes(),
525  procMeshes.faceProcAddressing(),
526  procMeshes.cellProcAddressing()
527  );
528 
529  // Pass2: reconstruct the cloud
530  forAllConstIters(allCloudObjects, iter)
531  {
532  const word cloudName = word::validate(iter.key());
533 
534  // Objects (on arbitrary processor)
535  const IOobjectList& cloudObjs = iter.val();
536 
537  Info<< "Reconstructing lagrangian fields for cloud "
538  << cloudName << nl << endl;
539 
540  reconstructor.reconstructPositions(cloudName);
541 
542  reconstructor.reconstructAllFields
543  (
544  cloudName,
545  cloudObjs,
546  selectedLagrangianFields
547  );
548  }
549  }
550  else
551  {
552  Info<< "No lagrangian fields" << nl << endl;
553  }
554  }
555 
556 
557  // If there are any FA fields, reconstruct them
558 
559  if (!doFiniteArea)
560  {
561  }
562  else if
563  (
564  faObjects.count<areaScalarField>()
565  || faObjects.count<areaVectorField>()
566  || faObjects.count<areaSphericalTensorField>()
567  || faObjects.count<areaSymmTensorField>()
568  || faObjects.count<areaTensorField>()
569  || faObjects.count<edgeScalarField>()
570  )
571  {
572  Info << "Reconstructing FA fields" << nl << endl;
573 
574  faMesh aMesh(mesh);
575 
576  processorFaMeshes procFaMeshes(procMeshes.meshes());
577 
578  faFieldReconstructor reconstructor
579  (
580  aMesh,
581  procFaMeshes.meshes(),
582  procFaMeshes.edgeProcAddressing(),
583  procFaMeshes.faceProcAddressing(),
584  procFaMeshes.boundaryProcAddressing()
585  );
586 
587  reconstructor.reconstructAllFields(faObjects);
588  }
589  else
590  {
591  Info << "No FA fields" << nl << endl;
592  }
593 
594  if (doReconstructSets)
595  {
596  // Scan to find all sets
597  HashTable<label> cSetNames;
598  HashTable<label> fSetNames;
599  HashTable<label> pSetNames;
600 
601  forAll(procMeshes.meshes(), proci)
602  {
603  const fvMesh& procMesh = procMeshes.meshes()[proci];
604 
605  // Note: look at sets in current time only or between
606  // mesh and current time?. For now current time. This will
607  // miss out on sets in intermediate times that have not
608  // been reconstructed.
609  IOobjectList objects
610  (
611  procMesh,
612  databases[0].timeName(), //procMesh.facesInstance()
613  polyMesh::meshSubDir/"sets"
614  );
615 
616  for (const IOobject& io : objects.csorted<cellSet>())
617  {
618  cSetNames.insert(io.name(), cSetNames.size());
619  }
620 
621  for (const IOobject& io : objects.csorted<faceSet>())
622  {
623  fSetNames.insert(io.name(), fSetNames.size());
624  }
625 
626  for (const IOobject& io : objects.csorted<pointSet>())
627  {
628  pSetNames.insert(io.name(), pSetNames.size());
629  }
630  }
631 
632  if (cSetNames.size() || fSetNames.size() || pSetNames.size())
633  {
634  // Construct all sets
635  PtrList<cellSet> cellSets(cSetNames.size());
636  PtrList<faceSet> faceSets(fSetNames.size());
637  PtrList<pointSet> pointSets(pSetNames.size());
638 
639  Info<< "Reconstructing sets:" << endl;
640  if (cSetNames.size())
641  {
642  Info<< " cellSets "
643  << cSetNames.sortedToc() << endl;
644  }
645  if (fSetNames.size())
646  {
647  Info<< " faceSets "
648  << fSetNames.sortedToc() << endl;
649  }
650  if (pSetNames.size())
651  {
652  Info<< " pointSets "
653  << pSetNames.sortedToc() << endl;
654  }
655 
656  // Load sets
657  forAll(procMeshes.meshes(), proci)
658  {
659  const fvMesh& procMesh = procMeshes.meshes()[proci];
660 
661  IOobjectList objects
662  (
663  procMesh,
664  databases[0].timeName(),
665  polyMesh::meshSubDir/"sets"
666  );
667 
668  // cellSets
669  const labelList& cellMap =
670  procMeshes.cellProcAddressing()[proci];
671 
672  for (const IOobject& io : objects.csorted<cellSet>())
673  {
674  // Load cellSet
675  const cellSet procSet(io);
676  const label seti = cSetNames[io.name()];
677  if (!cellSets.set(seti))
678  {
679  cellSets.set
680  (
681  seti,
682  new cellSet
683  (
684  mesh,
685  io.name(),
686  procSet.size()
687  )
688  );
689  }
690  cellSet& cSet = cellSets[seti];
691  cSet.instance() = runTime.timeName();
692 
693  for (const label celli : procSet)
694  {
695  cSet.insert(cellMap[celli]);
696  }
697  }
698 
699  // faceSets
700  const labelList& faceMap =
701  procMeshes.faceProcAddressing()[proci];
702 
703  for (const IOobject& io : objects.csorted<faceSet>())
704  {
705  // Load faceSet
706  const faceSet procSet(io);
707  const label seti = fSetNames[io.name()];
708  if (!faceSets.set(seti))
709  {
710  faceSets.set
711  (
712  seti,
713  new faceSet
714  (
715  mesh,
716  io.name(),
717  procSet.size()
718  )
719  );
720  }
721  faceSet& fSet = faceSets[seti];
722  fSet.instance() = runTime.timeName();
723 
724  for (const label facei : procSet)
725  {
726  fSet.insert(mag(faceMap[facei])-1);
727  }
728  }
729  // pointSets
730  const labelList& pointMap =
731  procMeshes.pointProcAddressing()[proci];
732 
733  for (const IOobject& io : objects.csorted<pointSet>())
734  {
735  // Load pointSet
736  const pointSet procSet(io);
737  const label seti = pSetNames[io.name()];
738  if (!pointSets.set(seti))
739  {
740  pointSets.set
741  (
742  seti,
743  new pointSet
744  (
745  mesh,
746  io.name(),
747  procSet.size()
748  )
749  );
750  }
751  pointSet& pSet = pointSets[seti];
752  pSet.instance() = runTime.timeName();
753 
754  for (const label pointi : procSet)
755  {
756  pSet.insert(pointMap[pointi]);
757  }
758  }
759  }
760 
761  // Write sets
762 
763  for (const auto& set : cellSets)
764  {
765  set.write();
766  }
767  for (const auto& set : faceSets)
768  {
769  set.write();
770  }
771  for (const auto& set : pointSets)
772  {
773  set.write();
774  }
775  }
776 
777 
778  // Reconstruct refinement data
779  {
780  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
781 
782  forAll(procMeshes.meshes(), procI)
783  {
784  const fvMesh& procMesh = procMeshes.meshes()[procI];
785 
786  procData.set
787  (
788  procI,
789  new hexRef8Data
790  (
791  IOobject
792  (
793  "dummy",
794  procMesh.time().timeName(),
795  polyMesh::meshSubDir,
796  procMesh,
797  IOobject::READ_IF_PRESENT,
798  IOobject::NO_WRITE,
799  IOobject::NO_REGISTER
800  )
801  )
802  );
803  }
804 
805  // Combine individual parts
806 
807  const PtrList<labelIOList>& cellAddr =
808  procMeshes.cellProcAddressing();
809 
810  UPtrList<const labelList> cellMaps(cellAddr.size());
811  forAll(cellAddr, i)
812  {
813  cellMaps.set(i, &cellAddr[i]);
814  }
815 
816  const PtrList<labelIOList>& pointAddr =
817  procMeshes.pointProcAddressing();
818 
819  UPtrList<const labelList> pointMaps(pointAddr.size());
820  forAll(pointAddr, i)
821  {
822  pointMaps.set(i, &pointAddr[i]);
823  }
824 
825  UPtrList<const hexRef8Data> procRefs(procData.size());
826  forAll(procData, i)
827  {
828  procRefs.set(i, &procData[i]);
829  }
830 
831  hexRef8Data
832  (
833  IOobject
834  (
835  "dummy",
836  mesh.time().timeName(),
837  polyMesh::meshSubDir,
838  mesh,
839  IOobject::NO_READ,
840  IOobject::NO_WRITE,
841  IOobject::NO_REGISTER
842  ),
843  cellMaps,
844  pointMaps,
845  procRefs
846  ).write();
847  }
848  }
849 
850 
851  // Reconstruct refinement data
852  {
853  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
854 
855  forAll(procMeshes.meshes(), procI)
856  {
857  const fvMesh& procMesh = procMeshes.meshes()[procI];
858 
859  procData.set
860  (
861  procI,
862  new hexRef8Data
863  (
864  IOobject
865  (
866  "dummy",
867  procMesh.time().timeName(),
868  polyMesh::meshSubDir,
869  procMesh,
870  IOobject::READ_IF_PRESENT,
871  IOobject::NO_WRITE,
872  IOobject::NO_REGISTER
873  )
874  )
875  );
876  }
877 
878  // Combine individual parts
879 
880  const PtrList<labelIOList>& cellAddr =
881  procMeshes.cellProcAddressing();
882 
883  UPtrList<const labelList> cellMaps(cellAddr.size());
884  forAll(cellAddr, i)
885  {
886  cellMaps.set(i, &cellAddr[i]);
887  }
888 
889  const PtrList<labelIOList>& pointAddr =
890  procMeshes.pointProcAddressing();
891 
892  UPtrList<const labelList> pointMaps(pointAddr.size());
893  forAll(pointAddr, i)
894  {
895  pointMaps.set(i, &pointAddr[i]);
896  }
897 
898  UPtrList<const hexRef8Data> procRefs(procData.size());
899  forAll(procData, i)
900  {
901  procRefs.set(i, &procData[i]);
902  }
903 
904  hexRef8Data
905  (
906  IOobject
907  (
908  "dummy",
909  mesh.time().timeName(),
910  polyMesh::meshSubDir,
911  mesh,
912  IOobject::NO_READ,
913  IOobject::NO_WRITE,
914  IOobject::NO_REGISTER
915  ),
916  cellMaps,
917  pointMaps,
918  procRefs
919  ).write();
920  }
921 
922  // If there is a "uniform" directory in the time region
923  // directory copy from the master processor
924  {
925  fileName uniformDir0
926  (
927  fileHandler().filePath
928  (
929  databases[0].timePath()/regionDir/"uniform"
930  )
931  );
932 
933  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
934  {
935  fileHandler().cp(uniformDir0, runTime.timePath()/regionDir);
936  }
937  }
938 
939  // For the first region of a multi-region case additionally
940  // copy the "uniform" directory in the time directory
941  if (regioni == 0 && !regionDir.empty())
942  {
943  fileName uniformDir0
944  (
945  fileHandler().filePath
946  (
947  databases[0].timePath()/"uniform"
948  )
949  );
950 
951  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
952  {
953  fileHandler().cp(uniformDir0, runTime.timePath());
954  }
955  }
956  }
957  }
958 
959  Info<< "\nEnd\n" << endl;
960 
961  return 0;
962 }
963 
964 
965 // ************************************************************************* //
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
List< bool > select(const label n, const labelUList &locations)
Construct a selection list of bools (all false) with the given pre-size, subsequently add specified l...
Definition: BitOps.C:134
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
GeometricField< symmTensor, faPatchField, areaMesh > areaSymmTensorField
Definition: areaFieldsFwd.H:84
wordList regionNames
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
bool allowFunctionObjects() const
The controlDict &#39;functions&#39; entry is allowed to be used.
Definition: argList.C:2159
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
StringType validate(const std::string &str, const UnaryPredicate &accept, const bool invert=false)
Return a copy of the input string with validated characters.
word timeName
Definition: getTimeIndex.H:3
bool allowLibs() const
The controlDict &#39;libs&#39; entry is allowed to be used. (eg, has not been disabled by the -no-libs option...
Definition: argList.C:2177
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const word cloudName(propsDict.get< word >("cloud"))
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:43
GeometricField< tensor, faPatchField, areaMesh > areaTensorField
Definition: areaFieldsFwd.H:88
const fileName & caseName() const noexcept
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:62
GeometricField< sphericalTensor, faPatchField, areaMesh > areaSphericalTensorField
Definition: areaFieldsFwd.H:80
const fileName & rootPath() const noexcept
Return root path.
Definition: argListI.H:56
bool readListIfPresent(const word &optName, List< T > &list) const
If named option is present, get a List of values treating a single entry like a list of size 1...
Definition: argListI.H:387
fileName path() const
Return the full path to the (processor local) case.
Definition: argListI.H:74
#define WarningInFunction
Report a warning using Foam::Warning.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
const word & regionDir
Required Classes.
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
List< instant > instantList
List of instants.
Definition: instantList.H:41
List< fileName > fileNameList
List of fileName.
Definition: fileNameList.H:32
Foam::argList args(argc, argv)
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:76
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:72