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