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