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-2018 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 
36 #include "fvCFD.H"
37 #include "IOobjectList.H"
38 #include "processorMeshes.H"
39 #include "regionProperties.H"
40 #include "fvFieldReconstructor.H"
42 #include "reconstructLagrangian.H"
43 
44 #include "cellSet.H"
45 #include "faceSet.H"
46 #include "pointSet.H"
47 
48 #include "hexRef8Data.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54  bool haveAllTimes
55  (
56  const HashSet<word>& masterTimeDirSet,
57  const instantList& timeDirs
58  )
59  {
60  // Loop over all times
61  forAll(timeDirs, timei)
62  {
63  if (!masterTimeDirSet.found(timeDirs[timei].name()))
64  {
65  return false;
66  }
67  }
68  return true;
69  }
70 }
71 
72 
73 int main(int argc, char *argv[])
74 {
76  (
77  "Reconstruct fields of a parallel case"
78  );
79 
80  // Enable -constant ... if someone really wants it
81  // Enable -withZero to prevent accidentally trashing the initial fields
82  timeSelector::addOptions(true, true);
84  #include "addRegionOption.H"
85  #include "addAllRegionsOption.H"
87  (
88  "fields",
89  "list",
90  "specify a list of fields to be reconstructed. Eg, '(U T p)' - "
91  "regular expressions not currently supported"
92  );
94  (
95  "noFields",
96  "skip reconstructing fields"
97  );
99  (
100  "lagrangianFields",
101  "list",
102  "specify a list of lagrangian fields to be reconstructed. Eg, '(U d)' -"
103  "regular expressions not currently supported, "
104  "positions always included."
105  );
107  (
108  "noLagrangian",
109  "skip reconstructing lagrangian positions and fields"
110  );
112  (
113  "noSets",
114  "skip reconstructing cellSets, faceSets, pointSets"
115  );
117  (
118  "newTimes",
119  "only reconstruct new times (i.e. that do not exist already)"
120  );
121 
122  #include "setRootCase.H"
123  #include "createTime.H"
124 
125  HashSet<word> selectedFields;
126  if (args.optionFound("fields"))
127  {
128  args.optionLookup("fields")() >> selectedFields;
129  }
130 
131  const bool noFields = args.optionFound("noFields");
132 
133  if (noFields)
134  {
135  Info<< "Skipping reconstructing fields"
136  << nl << endl;
137  }
138 
139  const bool noLagrangian = args.optionFound("noLagrangian");
140 
141  if (noLagrangian)
142  {
143  Info<< "Skipping reconstructing lagrangian positions and fields"
144  << nl << endl;
145  }
146 
147 
148  const bool noReconstructSets = args.optionFound("noSets");
149 
150  if (noReconstructSets)
151  {
152  Info<< "Skipping reconstructing cellSets, faceSets and pointSets"
153  << nl << endl;
154  }
155 
156 
157  HashSet<word> selectedLagrangianFields;
158  if (args.optionFound("lagrangianFields"))
159  {
160  if (noLagrangian)
161  {
163  << "Cannot specify noLagrangian and lagrangianFields "
164  << "options together."
165  << exit(FatalError);
166  }
167 
168  args.optionLookup("lagrangianFields")() >> selectedLagrangianFields;
169  }
170 
171  const wordList regionNames(selectRegionNames(args, runTime));
172 
173  // Determine the processor count
174  const label nProcs =
175  fileHandler().nProcs(args.path(), regionDir(regionNames[0]));
176 
177  if (!nProcs)
178  {
180  << "No processor* directories found"
181  << exit(FatalError);
182  }
183 
184  // Warn fileHandler of number of processors
185  const_cast<fileOperation&>(fileHandler()).setNProcs(nProcs);
186 
187  // Create the processor databases
188  PtrList<Time> databases(nProcs);
189 
190  forAll(databases, proci)
191  {
192  databases.set
193  (
194  proci,
195  new Time
196  (
198  args.rootPath(),
199  args.caseName()/fileName(word("processor") + name(proci))
200  )
201  );
202  }
203 
204  // Use the times list from the master processor
205  // and select a subset based on the command-line options
207  (
208  databases[0].times(),
209  args
210  );
211 
212  // Note that we do not set the runTime time so it is still the
213  // one set through the controlDict. The -time option
214  // only affects the selected set of times from processor0.
215  // - can be illogical
216  // + any point motion handled through mesh.readUpdate
217 
218  if (timeDirs.empty())
219  {
220  WarningInFunction << "No times selected" << endl;
221  exit(1);
222  }
223 
224  // Get current times if -newTimes
225  const bool newTimes = args.optionFound("newTimes");
226  instantList masterTimeDirs;
227  if (newTimes)
228  {
229  masterTimeDirs = runTime.times();
230  }
231 
232  HashSet<word> masterTimeDirSet(2*masterTimeDirs.size());
233  forAll(masterTimeDirs, i)
234  {
235  masterTimeDirSet.insert(masterTimeDirs[i].name());
236  }
237 
238  if
239  (
240  newTimes
241  && regionNames.size() == 1
242  && regionNames[0] == fvMesh::defaultRegion
243  && haveAllTimes(masterTimeDirSet, timeDirs)
244  )
245  {
246  Info<< "All times already reconstructed.\n\nEnd\n" << endl;
247  return 0;
248  }
249 
250  // Set all times on processor meshes equal to reconstructed mesh
251  forAll(databases, proci)
252  {
253  databases[proci].setTime(runTime);
254  }
255 
256  forAll(regionNames, regioni)
257  {
258  const word& regionName = regionNames[regioni];
259  const word& regionDir = Foam::regionDir(regionName);
260 
261  Info<< "\n\nReconstructing fields for mesh " << regionName << nl
262  << endl;
263 
264  fvMesh mesh
265  (
266  IOobject
267  (
268  regionName,
269  runTime.timeName(),
270  runTime,
272  )
273  );
274 
275 
276  // Read all meshes and addressing to reconstructed mesh
277  processorMeshes procMeshes(databases, regionName);
278 
279 
280  // Check face addressing for meshes that have been decomposed
281  // with a very old foam version
282  #include "checkFaceAddressingComp.H"
283 
284  // Loop over all times
285  forAll(timeDirs, timei)
286  {
287  if (newTimes && masterTimeDirSet.found(timeDirs[timei].name()))
288  {
289  Info<< "Skipping time " << timeDirs[timei].name()
290  << endl << endl;
291  continue;
292  }
293 
294 
295  // Set time for global database
296  runTime.setTime(timeDirs[timei], timei);
297 
298  Info<< "Time = " << runTime.timeName() << endl << endl;
299 
300  // Set time for all databases
301  forAll(databases, proci)
302  {
303  databases[proci].setTime(timeDirs[timei], timei);
304  }
305 
306  // Check if any new meshes need to be read.
307  fvMesh::readUpdateState meshStat = mesh.readUpdate();
308 
309  fvMesh::readUpdateState procStat = procMeshes.readUpdate();
310 
311  if (procStat == fvMesh::POINTS_MOVED)
312  {
313  // Reconstruct the points for moving mesh cases and write
314  // them out
315  procMeshes.reconstructPoints(mesh);
316  }
317  else if (meshStat != procStat)
318  {
320  << "readUpdate for the reconstructed mesh:"
321  << meshStat << nl
322  << "readUpdate for the processor meshes :"
323  << procStat << nl
324  << "These should be equal or your addressing"
325  << " might be incorrect."
326  << " Please check your time directories for any "
327  << "mesh directories." << endl;
328  }
329 
330 
331  // Get list of objects from processor0 database
332  IOobjectList objects
333  (
334  procMeshes.meshes()[0],
335  databases[0].timeName()
336  );
337 
338  if (!noFields)
339  {
340  // If there are any FV fields, reconstruct them
341  Info<< "Reconstructing FV fields" << nl << endl;
342 
343  fvFieldReconstructor fvReconstructor
344  (
345  mesh,
346  procMeshes.meshes(),
347  procMeshes.faceProcAddressing(),
348  procMeshes.cellProcAddressing(),
349  procMeshes.boundaryProcAddressing()
350  );
351 
352  fvReconstructor.reconstructFvVolumeInternalFields<scalar>
353  (
354  objects,
355  selectedFields
356  );
357  fvReconstructor.reconstructFvVolumeInternalFields<vector>
358  (
359  objects,
360  selectedFields
361  );
362  fvReconstructor.reconstructFvVolumeInternalFields
364  (
365  objects,
366  selectedFields
367  );
368  fvReconstructor.reconstructFvVolumeInternalFields<symmTensor>
369  (
370  objects,
371  selectedFields
372  );
373  fvReconstructor.reconstructFvVolumeInternalFields<tensor>
374  (
375  objects,
376  selectedFields
377  );
378 
379  fvReconstructor.reconstructFvVolumeFields<scalar>
380  (
381  objects,
382  selectedFields
383  );
384  fvReconstructor.reconstructFvVolumeFields<vector>
385  (
386  objects,
387  selectedFields
388  );
389  fvReconstructor.reconstructFvVolumeFields<sphericalTensor>
390  (
391  objects,
392  selectedFields
393  );
394  fvReconstructor.reconstructFvVolumeFields<symmTensor>
395  (
396  objects,
397  selectedFields
398  );
399  fvReconstructor.reconstructFvVolumeFields<tensor>
400  (
401  objects,
402  selectedFields
403  );
404 
405  fvReconstructor.reconstructFvSurfaceFields<scalar>
406  (
407  objects,
408  selectedFields
409  );
410  fvReconstructor.reconstructFvSurfaceFields<vector>
411  (
412  objects,
413  selectedFields
414  );
415  fvReconstructor.reconstructFvSurfaceFields<sphericalTensor>
416  (
417  objects,
418  selectedFields
419  );
420  fvReconstructor.reconstructFvSurfaceFields<symmTensor>
421  (
422  objects,
423  selectedFields
424  );
425  fvReconstructor.reconstructFvSurfaceFields<tensor>
426  (
427  objects,
428  selectedFields
429  );
430 
431  if (fvReconstructor.nReconstructed() == 0)
432  {
433  Info<< "No FV fields" << nl << endl;
434  }
435  }
436 
437  if (!noFields)
438  {
439  Info<< "Reconstructing point fields" << nl << endl;
440 
441  const pointMesh& pMesh = pointMesh::New(mesh);
442  PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
443 
444  forAll(pMeshes, proci)
445  {
446  pMeshes.set
447  (
448  proci,
449  new pointMesh(procMeshes.meshes()[proci])
450  );
451  }
452 
453  pointFieldReconstructor pointReconstructor
454  (
455  pMesh,
456  pMeshes,
457  procMeshes.pointProcAddressing(),
458  procMeshes.boundaryProcAddressing()
459  );
460 
461  pointReconstructor.reconstructFields<scalar>
462  (
463  objects,
464  selectedFields
465  );
466  pointReconstructor.reconstructFields<vector>
467  (
468  objects,
469  selectedFields
470  );
471  pointReconstructor.reconstructFields<sphericalTensor>
472  (
473  objects,
474  selectedFields
475  );
476  pointReconstructor.reconstructFields<symmTensor>
477  (
478  objects,
479  selectedFields
480  );
481  pointReconstructor.reconstructFields<tensor>
482  (
483  objects,
484  selectedFields
485  );
486 
487  if (pointReconstructor.nReconstructed() == 0)
488  {
489  Info<< "No point fields" << nl << endl;
490  }
491  }
492 
493 
494  // If there are any clouds, reconstruct them.
495  // The problem is that a cloud of size zero will not get written so
496  // in pass 1 we determine the cloud names and per cloud name the
497  // fields. Note that the fields are stored as IOobjectList from
498  // the first processor that has them. They are in pass2 only used
499  // for name and type (scalar, vector etc).
500 
501  if (!noLagrangian)
502  {
503  HashTable<IOobjectList> cloudObjects;
504 
505  forAll(databases, proci)
506  {
507  fileName lagrangianDir
508  (
509  fileHandler().filePath
510  (
511  databases[proci].timePath()
512  / regionDir
513  / cloud::prefix
514  )
515  );
516 
517  fileNameList cloudDirs;
518  if (!lagrangianDir.empty())
519  {
520  cloudDirs = fileHandler().readDir
521  (
522  lagrangianDir,
524  );
525  }
526 
527  forAll(cloudDirs, i)
528  {
529  // Check if we already have cloud objects for this
530  // cloudname
532  cloudObjects.find(cloudDirs[i]);
533 
534  if (iter == cloudObjects.end())
535  {
536  // Do local scan for valid cloud objects
537  IOobjectList sprayObjs
538  (
539  procMeshes.meshes()[proci],
540  databases[proci].timeName(),
541  cloud::prefix/cloudDirs[i]
542  );
543 
544  IOobject* positionsPtr =
545  sprayObjs.lookup(word("positions"));
546 
547  if (positionsPtr)
548  {
549  cloudObjects.insert(cloudDirs[i], sprayObjs);
550  }
551  }
552  }
553  }
554 
555 
556  if (cloudObjects.size())
557  {
558  // Pass2: reconstruct the cloud
559  forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
560  {
561  const word cloudName = string::validate<word>
562  (
563  iter.key()
564  );
565 
566  // Objects (on arbitrary processor)
567  const IOobjectList& sprayObjs = iter();
568 
569  Info<< "Reconstructing lagrangian fields for cloud "
570  << cloudName << nl << endl;
571 
573  (
574  mesh,
575  cloudName,
576  procMeshes.meshes(),
577  procMeshes.faceProcAddressing(),
578  procMeshes.cellProcAddressing()
579  );
580  reconstructLagrangianFields<label>
581  (
582  cloudName,
583  mesh,
584  procMeshes.meshes(),
585  sprayObjs,
586  selectedLagrangianFields
587  );
588  reconstructLagrangianFieldFields<label>
589  (
590  cloudName,
591  mesh,
592  procMeshes.meshes(),
593  sprayObjs,
594  selectedLagrangianFields
595  );
596  reconstructLagrangianFields<scalar>
597  (
598  cloudName,
599  mesh,
600  procMeshes.meshes(),
601  sprayObjs,
602  selectedLagrangianFields
603  );
604  reconstructLagrangianFieldFields<scalar>
605  (
606  cloudName,
607  mesh,
608  procMeshes.meshes(),
609  sprayObjs,
610  selectedLagrangianFields
611  );
612  reconstructLagrangianFields<vector>
613  (
614  cloudName,
615  mesh,
616  procMeshes.meshes(),
617  sprayObjs,
618  selectedLagrangianFields
619  );
620  reconstructLagrangianFieldFields<vector>
621  (
622  cloudName,
623  mesh,
624  procMeshes.meshes(),
625  sprayObjs,
626  selectedLagrangianFields
627  );
628  reconstructLagrangianFields<sphericalTensor>
629  (
630  cloudName,
631  mesh,
632  procMeshes.meshes(),
633  sprayObjs,
634  selectedLagrangianFields
635  );
636  reconstructLagrangianFieldFields<sphericalTensor>
637  (
638  cloudName,
639  mesh,
640  procMeshes.meshes(),
641  sprayObjs,
642  selectedLagrangianFields
643  );
644  reconstructLagrangianFields<symmTensor>
645  (
646  cloudName,
647  mesh,
648  procMeshes.meshes(),
649  sprayObjs,
650  selectedLagrangianFields
651  );
652  reconstructLagrangianFieldFields<symmTensor>
653  (
654  cloudName,
655  mesh,
656  procMeshes.meshes(),
657  sprayObjs,
658  selectedLagrangianFields
659  );
660  reconstructLagrangianFields<tensor>
661  (
662  cloudName,
663  mesh,
664  procMeshes.meshes(),
665  sprayObjs,
666  selectedLagrangianFields
667  );
668  reconstructLagrangianFieldFields<tensor>
669  (
670  cloudName,
671  mesh,
672  procMeshes.meshes(),
673  sprayObjs,
674  selectedLagrangianFields
675  );
676  }
677  }
678  else
679  {
680  Info<< "No lagrangian fields" << nl << endl;
681  }
682  }
683 
684 
685  if (!noReconstructSets)
686  {
687  // Scan to find all sets
688  HashTable<label> cSetNames;
689  HashTable<label> fSetNames;
690  HashTable<label> pSetNames;
691 
692  forAll(procMeshes.meshes(), proci)
693  {
694  const fvMesh& procMesh = procMeshes.meshes()[proci];
695 
696  // Note: look at sets in current time only or between
697  // mesh and current time?. For now current time. This will
698  // miss out on sets in intermediate times that have not
699  // been reconstructed.
700  IOobjectList objects
701  (
702  procMesh,
703  databases[0].timeName(), // procMesh.facesInstance()
704  polyMesh::meshSubDir/"sets"
705  );
706 
707  IOobjectList cSets(objects.lookupClass(cellSet::typeName));
708  forAllConstIter(IOobjectList, cSets, iter)
709  {
710  cSetNames.insert(iter.key(), cSetNames.size());
711  }
712 
713  IOobjectList fSets(objects.lookupClass(faceSet::typeName));
714  forAllConstIter(IOobjectList, fSets, iter)
715  {
716  fSetNames.insert(iter.key(), fSetNames.size());
717  }
718  IOobjectList pSets(objects.lookupClass(pointSet::typeName));
719  forAllConstIter(IOobjectList, pSets, iter)
720  {
721  pSetNames.insert(iter.key(), pSetNames.size());
722  }
723  }
724 
725  if (cSetNames.size() || fSetNames.size() || pSetNames.size())
726  {
727  // Construct all sets
728  PtrList<cellSet> cellSets(cSetNames.size());
729  PtrList<faceSet> faceSets(fSetNames.size());
730  PtrList<pointSet> pointSets(pSetNames.size());
731 
732  Info<< "Reconstructing sets:" << endl;
733  if (cSetNames.size())
734  {
735  Info<< " cellSets "
736  << cSetNames.sortedToc() << endl;
737  }
738  if (fSetNames.size())
739  {
740  Info<< " faceSets "
741  << fSetNames.sortedToc() << endl;
742  }
743  if (pSetNames.size())
744  {
745  Info<< " pointSets "
746  << pSetNames.sortedToc() << endl;
747  }
748 
749  // Load sets
750  forAll(procMeshes.meshes(), proci)
751  {
752  const fvMesh& procMesh = procMeshes.meshes()[proci];
753 
754  IOobjectList objects
755  (
756  procMesh,
757  databases[0].timeName(),
758  polyMesh::meshSubDir/"sets"
759  );
760 
761  // cellSets
762  const labelList& cellMap =
763  procMeshes.cellProcAddressing()[proci];
764 
765  IOobjectList cSets
766  (
767  objects.lookupClass(cellSet::typeName)
768  );
769 
770  forAllConstIter(IOobjectList, cSets, iter)
771  {
772  // Load cellSet
773  const cellSet procSet(*iter());
774  label setI = cSetNames[iter.key()];
775  if (!cellSets.set(setI))
776  {
777  cellSets.set
778  (
779  setI,
780  new cellSet
781  (
782  mesh,
783  iter.key(),
784  procSet.size()
785  )
786  );
787  }
788  cellSet& cSet = cellSets[setI];
789  cSet.instance() = runTime.timeName();
790 
791  forAllConstIter(cellSet, procSet, iter)
792  {
793  cSet.insert(cellMap[iter.key()]);
794  }
795  }
796 
797  // faceSets
798  const labelList& faceMap =
799  procMeshes.faceProcAddressing()[proci];
800 
801  IOobjectList fSets
802  (
803  objects.lookupClass(faceSet::typeName)
804  );
805 
806  forAllConstIter(IOobjectList, fSets, iter)
807  {
808  // Load faceSet
809  const faceSet procSet(*iter());
810  label setI = fSetNames[iter.key()];
811  if (!faceSets.set(setI))
812  {
813  faceSets.set
814  (
815  setI,
816  new faceSet
817  (
818  mesh,
819  iter.key(),
820  procSet.size()
821  )
822  );
823  }
824  faceSet& fSet = faceSets[setI];
825  fSet.instance() = runTime.timeName();
826 
827  forAllConstIter(faceSet, procSet, iter)
828  {
829  fSet.insert(mag(faceMap[iter.key()])-1);
830  }
831  }
832  // pointSets
833  const labelList& pointMap =
834  procMeshes.pointProcAddressing()[proci];
835 
836  IOobjectList pSets
837  (
838  objects.lookupClass(pointSet::typeName)
839  );
840  forAllConstIter(IOobjectList, pSets, iter)
841  {
842  // Load pointSet
843  const pointSet propSet(*iter());
844  label setI = pSetNames[iter.key()];
845  if (!pointSets.set(setI))
846  {
847  pointSets.set
848  (
849  setI,
850  new pointSet
851  (
852  mesh,
853  iter.key(),
854  propSet.size()
855  )
856  );
857  }
858  pointSet& pSet = pointSets[setI];
859  pSet.instance() = runTime.timeName();
860 
861  forAllConstIter(pointSet, propSet, iter)
862  {
863  pSet.insert(pointMap[iter.key()]);
864  }
865  }
866  }
867 
868  // Write sets
869  forAll(cellSets, i)
870  {
871  cellSets[i].write();
872  }
873  forAll(faceSets, i)
874  {
875  faceSets[i].write();
876  }
877  forAll(pointSets, i)
878  {
879  pointSets[i].write();
880  }
881  }
882  }
883 
884 
885  // Reconstruct refinement data
886  {
887  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
888 
889  forAll(procMeshes.meshes(), procI)
890  {
891  const fvMesh& procMesh = procMeshes.meshes()[procI];
892 
893  procData.set
894  (
895  procI,
896  new hexRef8Data
897  (
898  IOobject
899  (
900  "dummy",
901  procMesh.time().timeName(),
903  procMesh,
906  false
907  )
908  )
909  );
910  }
911 
912  // Combine individual parts
913 
914  const PtrList<labelIOList>& cellAddr =
915  procMeshes.cellProcAddressing();
916 
917  UPtrList<const labelList> cellMaps(cellAddr.size());
918  forAll(cellAddr, i)
919  {
920  cellMaps.set(i, &cellAddr[i]);
921  }
922 
923  const PtrList<labelIOList>& pointAddr =
924  procMeshes.pointProcAddressing();
925 
926  UPtrList<const labelList> pointMaps(pointAddr.size());
927  forAll(pointAddr, i)
928  {
929  pointMaps.set(i, &pointAddr[i]);
930  }
931 
932  UPtrList<const hexRef8Data> procRefs(procData.size());
933  forAll(procData, i)
934  {
935  procRefs.set(i, &procData[i]);
936  }
937 
938  hexRef8Data
939  (
940  IOobject
941  (
942  "dummy",
943  mesh.time().timeName(),
945  mesh,
948  false
949  ),
950  cellMaps,
951  pointMaps,
952  procRefs
953  ).write();
954  }
955 
956  // If there is a "uniform" directory in the time region
957  // directory copy from the master processor
958  {
959  fileName uniformDir0
960  (
961  fileHandler().filePath
962  (
963  databases[0].timePath()/regionDir/"uniform"
964  )
965  );
966 
967  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
968  {
969  fileHandler().cp(uniformDir0, runTime.timePath()/regionDir);
970  }
971  }
972 
973  // For the first region of a multi-region case additionally
974  // copy the "uniform" directory in the time directory
975  if (regioni == 0 && regionDir != word::null)
976  {
977  fileName uniformDir0
978  (
979  fileHandler().filePath
980  (
981  databases[0].timePath()/"uniform"
982  )
983  );
984 
985  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
986  {
987  fileHandler().cp(uniformDir0, runTime.timePath());
988  }
989  }
990  }
991  }
992 
993  Info<< "\nEnd\n" << endl;
994 
995  return 0;
996 }
997 
998 
999 // ************************************************************************* //
List< instant > instantList
List of instants.
Definition: instantList.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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:319
fileName timePath() const
Return current time path.
Definition: Time.H:282
const fileName & rootPath() const
Return root path.
Definition: argListI.H:36
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static void noParallel()
Remove the parallel options.
Definition: argList.C:174
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reconstructLagrangianPositions(const polyMesh &mesh, const word &cloudName, const PtrList< fvMesh > &meshes, const PtrList< labelIOList > &faceProcAddressing, const PtrList< labelIOList > &cellProcAddressing)
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
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)
dynamicFvMesh & mesh
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:127
static const word null
An empty word.
Definition: word.H:77
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:873
const fileOperation & fileHandler()
Get current file handler.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:198
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
instantList select(const instantList &) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:100
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:42
static const char nl
Definition: Ostream.H:265
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:60
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?
instantList times() const
Search the case for valid time directories.
Definition: Time.C:642
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
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:117
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:158
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
Foam::argList args(argc, argv)
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
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:114