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-2021 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::addOption
96  (
97  "fields",
98  "wordRes",
99  "Specify single or multiple fields to reconstruct (all by default)."
100  " Eg, 'T' or '(p T U \"alpha.*\")'"
101  );
102  argList::addBoolOption
103  (
104  "no-fields", // noFields
105  "Skip reconstructing fields"
106  );
107  argList::addOptionCompat("no-fields", {"noFields", 2106});
108  argList::addOption
109  (
110  "lagrangianFields",
111  "wordRes",
112  "Specify single or multiple lagrangian fields to reconstruct"
113  " (all by default)."
114  " Eg, '(U d)'"
115  " - Positions are always included."
116  );
117  argList::addBoolOption
118  (
119  "no-lagrangian", // noLagrangian
120  "Skip reconstructing lagrangian positions and fields"
121  );
122  argList::addOptionCompat("no-lagrangian", {"noLagrangian", 2106});
123 
124  argList::addBoolOption
125  (
126  "no-sets",
127  "Skip reconstructing cellSets, faceSets, pointSets"
128  );
129  argList::addOptionCompat("no-sets", {"noSets", 2106});
130 
131  argList::addBoolOption
132  (
133  "newTimes",
134  "Only reconstruct new times (i.e. that do not exist already)"
135  );
136 
137  #include "setRootCase.H"
138  #include "createTime.H"
139 
140 
141  const bool doFields = !args.found("no-fields");
142  wordRes selectedFields;
143 
144  if (doFields)
145  {
146  args.readListIfPresent<wordRe>("fields", selectedFields);
147  }
148  else
149  {
150  Info<< "Skipping reconstructing fields";
151  if (args.found("fields"))
152  {
153  Info<< ". Ignore -fields option";
154  }
155  Info<< nl << endl;
156  }
157 
158 
159  const bool doFiniteArea = !args.found("no-finite-area");
160  if (!doFiniteArea)
161  {
162  Info<< "Skipping reconstructing finiteArea mesh/fields"
163  << nl << endl;
164  }
165 
166 
167  const bool doLagrangian = !args.found("no-lagrangian");
168  wordRes selectedLagrangianFields;
169 
170  if (doLagrangian)
171  {
172  args.readListIfPresent<wordRe>
173  (
174  "lagrangianFields", selectedLagrangianFields
175  );
176  }
177  else
178  {
179  Info<< "Skipping reconstructing lagrangian positions/fields";
180  if (args.found("lagrangianFields"))
181  {
182  Info<< ". Ignore -lagrangianFields option";
183  }
184  Info<< nl << endl;
185  }
186 
187 
188  const bool doReconstructSets = !args.found("no-sets");
189 
190  if (!doReconstructSets)
191  {
192  Info<< "Skipping reconstructing cellSets, faceSets and pointSets"
193  << nl << endl;
194  }
195 
196  const bool newTimes = args.found("newTimes");
197 
198  // Get region names
199  #include "getAllRegionOptions.H"
200 
201  // Determine the processor count
202  label nProcs{0};
203 
204  if (regionNames.empty())
205  {
207  << "No regions specified or detected."
208  << exit(FatalError);
209  }
210  else if (regionNames[0] == polyMesh::defaultRegion)
211  {
212  nProcs = fileHandler().nProcs(args.path());
213  }
214  else
215  {
216  nProcs = fileHandler().nProcs(args.path(), regionNames[0]);
217 
218  if (regionNames.size() == 1)
219  {
220  Info<< "Using region: " << regionNames[0] << nl << endl;
221  }
222  }
223 
224  if (!nProcs)
225  {
227  << "No processor* directories found"
228  << exit(FatalError);
229  }
230 
231  // Warn fileHandler of number of processors
232  const_cast<fileOperation&>(fileHandler()).setNProcs(nProcs);
233 
234  // Create the processor databases
235  PtrList<Time> databases(nProcs);
236 
237  forAll(databases, proci)
238  {
239  databases.set
240  (
241  proci,
242  new Time
243  (
244  Time::controlDictName,
245  args.rootPath(),
246  args.caseName()/("processor" + Foam::name(proci))
247  )
248  );
249  }
250 
251  // Use the times list from the master processor
252  // and select a subset based on the command-line options
253  instantList timeDirs = timeSelector::select
254  (
255  databases[0].times(),
256  args
257  );
258 
259  // Note that we do not set the runTime time so it is still the
260  // one set through the controlDict. The -time option
261  // only affects the selected set of times from processor0.
262  // - can be illogical
263  // + any point motion handled through mesh.readUpdate
264 
265 
266  if (timeDirs.empty())
267  {
268  WarningInFunction << "No times selected";
269  exit(1);
270  }
271 
272 
273  // Get current times if -newTimes
274  instantList masterTimeDirs;
275  if (newTimes)
276  {
277  masterTimeDirs = runTime.times();
278  }
279  wordHashSet masterTimeDirSet(2*masterTimeDirs.size());
280  for (const instant& t : masterTimeDirs)
281  {
282  masterTimeDirSet.insert(t.name());
283  }
284 
285 
286  // Set all times on processor meshes equal to reconstructed mesh
287  forAll(databases, proci)
288  {
289  databases[proci].setTime(runTime);
290  }
291 
292 
293  forAll(regionNames, regioni)
294  {
295  const word& regionName = regionNames[regioni];
296  const word& regionDir =
297  (
298  regionName != polyMesh::defaultRegion
299  ? regionName
300  : word::null
301  );
302 
303  Info<< "\n\nReconstructing fields" << nl
304  << "region=" << regionName << nl << endl;
305 
306  if
307  (
308  newTimes
309  && regionNames.size() == 1
310  && regionDir.empty()
311  && haveAllTimes(masterTimeDirSet, timeDirs)
312  )
313  {
314  Info<< "Skipping region " << regionName
315  << " since already have all times"
316  << endl << endl;
317  continue;
318  }
319 
320 
321  fvMesh mesh
322  (
323  IOobject
324  (
325  regionName,
326  runTime.timeName(),
327  runTime,
329  )
330  );
331 
332 
333  // Read all meshes and addressing to reconstructed mesh
334  processorMeshes procMeshes(databases, regionName);
335 
336 
337  // Check face addressing for meshes that have been decomposed
338  // with a very old foam version
339  #include "checkFaceAddressingComp.H"
340 
341  // Loop over all times
342  forAll(timeDirs, timei)
343  {
344  if (newTimes && masterTimeDirSet.found(timeDirs[timei].name()))
345  {
346  Info<< "Skipping time " << timeDirs[timei].name()
347  << endl << endl;
348  continue;
349  }
350 
351 
352  // Set time for global database
353  runTime.setTime(timeDirs[timei], timei);
354 
355  Info<< "Time = " << runTime.timeName() << endl << endl;
356 
357  // Set time for all databases
358  forAll(databases, proci)
359  {
360  databases[proci].setTime(timeDirs[timei], timei);
361  }
362 
363  // Check if any new meshes need to be read.
364  fvMesh::readUpdateState meshStat = mesh.readUpdate();
365 
366  fvMesh::readUpdateState procStat = procMeshes.readUpdate();
367 
368  if (procStat == fvMesh::POINTS_MOVED)
369  {
370  // Reconstruct the points for moving mesh cases and write
371  // them out
372  procMeshes.reconstructPoints(mesh);
373  }
374  else if (meshStat != procStat)
375  {
377  << "readUpdate for the reconstructed mesh:"
378  << meshStat << nl
379  << "readUpdate for the processor meshes :"
380  << procStat << nl
381  << "These should be equal or your addressing"
382  << " might be incorrect."
383  << " Please check your time directories for any "
384  << "mesh directories." << endl;
385  }
386 
387 
388  // Get list of objects from processor0 database
389  IOobjectList objects
390  (
391  procMeshes.meshes()[0],
392  databases[0].timeName()
393  );
394 
395  if (doFields)
396  {
397  // If there are any FV fields, reconstruct them
398  Info<< "Reconstructing FV fields" << nl << endl;
399 
400  fvFieldReconstructor reconstructor
401  (
402  mesh,
403  procMeshes.meshes(),
404  procMeshes.faceProcAddressing(),
405  procMeshes.cellProcAddressing(),
406  procMeshes.boundaryProcAddressing()
407  );
408 
409  reconstructor.reconstructFvVolumeInternalFields<scalar>
410  (
411  objects,
412  selectedFields
413  );
414  reconstructor.reconstructFvVolumeInternalFields<vector>
415  (
416  objects,
417  selectedFields
418  );
419  reconstructor.reconstructFvVolumeInternalFields<sphericalTensor>
420  (
421  objects,
422  selectedFields
423  );
424  reconstructor.reconstructFvVolumeInternalFields<symmTensor>
425  (
426  objects,
427  selectedFields
428  );
429  reconstructor.reconstructFvVolumeInternalFields<tensor>
430  (
431  objects,
432  selectedFields
433  );
434 
435  reconstructor.reconstructFvVolumeFields<scalar>
436  (
437  objects,
438  selectedFields
439  );
440  reconstructor.reconstructFvVolumeFields<vector>
441  (
442  objects,
443  selectedFields
444  );
445  reconstructor.reconstructFvVolumeFields<sphericalTensor>
446  (
447  objects,
448  selectedFields
449  );
450  reconstructor.reconstructFvVolumeFields<symmTensor>
451  (
452  objects,
453  selectedFields
454  );
455  reconstructor.reconstructFvVolumeFields<tensor>
456  (
457  objects,
458  selectedFields
459  );
460 
461  reconstructor.reconstructFvSurfaceFields<scalar>
462  (
463  objects,
464  selectedFields
465  );
466  reconstructor.reconstructFvSurfaceFields<vector>
467  (
468  objects,
469  selectedFields
470  );
471  reconstructor.reconstructFvSurfaceFields<sphericalTensor>
472  (
473  objects,
474  selectedFields
475  );
476  reconstructor.reconstructFvSurfaceFields<symmTensor>
477  (
478  objects,
479  selectedFields
480  );
481  reconstructor.reconstructFvSurfaceFields<tensor>
482  (
483  objects,
484  selectedFields
485  );
486 
487  if (reconstructor.nReconstructed() == 0)
488  {
489  Info<< "No FV fields" << nl << endl;
490  }
491  }
492 
493  if (doFields)
494  {
495  Info<< "Reconstructing point fields" << nl << endl;
496 
497  const pointMesh& pMesh = pointMesh::New(mesh);
498  PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
499 
500  forAll(pMeshes, proci)
501  {
502  pMeshes.set
503  (
504  proci,
505  new pointMesh(procMeshes.meshes()[proci])
506  );
507  }
508 
509  pointFieldReconstructor reconstructor
510  (
511  pMesh,
512  pMeshes,
513  procMeshes.pointProcAddressing(),
514  procMeshes.boundaryProcAddressing()
515  );
516 
517  reconstructor.reconstructFields<scalar>
518  (
519  objects,
520  selectedFields
521  );
522  reconstructor.reconstructFields<vector>
523  (
524  objects,
525  selectedFields
526  );
527  reconstructor.reconstructFields<sphericalTensor>
528  (
529  objects,
530  selectedFields
531  );
532  reconstructor.reconstructFields<symmTensor>
533  (
534  objects,
535  selectedFields
536  );
537  reconstructor.reconstructFields<tensor>
538  (
539  objects,
540  selectedFields
541  );
542 
543  if (reconstructor.nReconstructed() == 0)
544  {
545  Info<< "No point fields" << nl << endl;
546  }
547  }
548 
549 
550  // If there are any clouds, reconstruct them.
551  // The problem is that a cloud of size zero will not get written so
552  // in pass 1 we determine the cloud names and per cloud name the
553  // fields. Note that the fields are stored as IOobjectList from
554  // the first processor that has them. They are in pass2 only used
555  // for name and type (scalar, vector etc).
556 
557  if (doLagrangian)
558  {
559  HashTable<IOobjectList> allCloudObjects;
560 
561  forAll(databases, proci)
562  {
563  fileName lagrangianDir
564  (
565  fileHandler().filePath
566  (
567  databases[proci].timePath()
568  / regionDir
569  / cloud::prefix
570  )
571  );
572 
573  fileNameList cloudDirs;
574  if (!lagrangianDir.empty())
575  {
576  cloudDirs = fileHandler().readDir
577  (
578  lagrangianDir,
579  fileName::DIRECTORY
580  );
581  }
582 
583  for (const fileName& cloudDir : cloudDirs)
584  {
585  // Check if we already have cloud objects for this
586  // cloudname
587  if (!allCloudObjects.found(cloudDir))
588  {
589  // Do local scan for valid cloud objects
590  IOobjectList localObjs
591  (
592  procMeshes.meshes()[proci],
593  databases[proci].timeName(),
594  cloud::prefix/cloudDir
595  );
596 
597  if
598  (
599  localObjs.found("coordinates")
600  || localObjs.found("positions")
601  )
602  {
603  allCloudObjects.insert(cloudDir, localObjs);
604  }
605  }
606  }
607  }
608 
609 
610  if (allCloudObjects.size())
611  {
612  lagrangianReconstructor reconstructor
613  (
614  mesh,
615  procMeshes.meshes(),
616  procMeshes.faceProcAddressing(),
617  procMeshes.cellProcAddressing()
618  );
619 
620  // Pass2: reconstruct the cloud
621  forAllConstIters(allCloudObjects, iter)
622  {
623  const word cloudName = word::validate(iter.key());
624 
625  // Objects (on arbitrary processor)
626  const IOobjectList& cloudObjs = iter.val();
627 
628  Info<< "Reconstructing lagrangian fields for cloud "
629  << cloudName << nl << endl;
630 
631  reconstructor.reconstructPositions(cloudName);
632 
633  reconstructor.reconstructFields<label>
634  (
635  cloudName,
636  cloudObjs,
637  selectedLagrangianFields
638  );
639  reconstructor.reconstructFieldFields<label>
640  (
641  cloudName,
642  cloudObjs,
643  selectedLagrangianFields
644  );
645 
646  reconstructor.reconstructFields<scalar>
647  (
648  cloudName,
649  cloudObjs,
650  selectedLagrangianFields
651  );
652  reconstructor.reconstructFieldFields<scalar>
653  (
654  cloudName,
655  cloudObjs,
656  selectedLagrangianFields
657  );
658 
659  reconstructor.reconstructFields<vector>
660  (
661  cloudName,
662  cloudObjs,
663  selectedLagrangianFields
664  );
665  reconstructor.reconstructFieldFields<vector>
666  (
667  cloudName,
668  cloudObjs,
669  selectedLagrangianFields
670  );
671 
672  reconstructor.reconstructFields<sphericalTensor>
673  (
674  cloudName,
675  cloudObjs,
676  selectedLagrangianFields
677  );
678  reconstructor.reconstructFieldFields<sphericalTensor>
679  (
680  cloudName,
681  cloudObjs,
682  selectedLagrangianFields
683  );
684 
685  reconstructor.reconstructFields<symmTensor>
686  (
687  cloudName,
688  cloudObjs,
689  selectedLagrangianFields
690  );
691  reconstructor.reconstructFieldFields<symmTensor>
692  (
693  cloudName,
694  cloudObjs,
695  selectedLagrangianFields
696  );
697 
698  reconstructor.reconstructFields<tensor>
699  (
700  cloudName,
701  cloudObjs,
702  selectedLagrangianFields
703  );
704  reconstructor.reconstructFieldFields<tensor>
705  (
706  cloudName,
707  cloudObjs,
708  selectedLagrangianFields
709  );
710  }
711  }
712  else
713  {
714  Info<< "No lagrangian fields" << nl << endl;
715  }
716  }
717 
718 
719  // If there are any FA fields, reconstruct them
720 
721  if (!doFiniteArea)
722  {
723  }
724  else if
725  (
726  objects.lookupClass(areaScalarField::typeName).size()
727  || objects.lookupClass(areaVectorField::typeName).size()
728  || objects.lookupClass(areaSphericalTensorField::typeName).size()
729  || objects.lookupClass(areaSymmTensorField::typeName).size()
730  || objects.lookupClass(areaTensorField::typeName).size()
731  || objects.lookupClass(edgeScalarField::typeName).size()
732  )
733  {
734  Info << "Reconstructing FA fields" << nl << endl;
735 
736  faMesh aMesh(mesh);
737 
738  processorFaMeshes procFaMeshes(procMeshes.meshes());
739 
740  faFieldReconstructor reconstructor
741  (
742  aMesh,
743  procFaMeshes.meshes(),
744  procFaMeshes.edgeProcAddressing(),
745  procFaMeshes.faceProcAddressing(),
746  procFaMeshes.boundaryProcAddressing()
747  );
748 
749  reconstructor.reconstructFaAreaFields<scalar>(objects);
750  reconstructor.reconstructFaAreaFields<vector>(objects);
751  reconstructor.reconstructFaAreaFields<sphericalTensor>(objects);
752  reconstructor.reconstructFaAreaFields<symmTensor>(objects);
753  reconstructor.reconstructFaAreaFields<tensor>(objects);
754 
755  reconstructor.reconstructFaEdgeFields<scalar>(objects);
756  }
757  else
758  {
759  Info << "No FA fields" << nl << endl;
760  }
761 
762  if (doReconstructSets)
763  {
764  // Scan to find all sets
765  HashTable<label> cSetNames;
766  HashTable<label> fSetNames;
767  HashTable<label> pSetNames;
768 
769  forAll(procMeshes.meshes(), proci)
770  {
771  const fvMesh& procMesh = procMeshes.meshes()[proci];
772 
773  // Note: look at sets in current time only or between
774  // mesh and current time?. For now current time. This will
775  // miss out on sets in intermediate times that have not
776  // been reconstructed.
777  IOobjectList objects
778  (
779  procMesh,
780  databases[0].timeName(), //procMesh.facesInstance()
781  polyMesh::meshSubDir/"sets"
782  );
783 
784  IOobjectList cSets(objects.lookupClass(cellSet::typeName));
785  forAllConstIters(cSets, iter)
786  {
787  cSetNames.insert(iter.key(), cSetNames.size());
788  }
789 
790  IOobjectList fSets(objects.lookupClass(faceSet::typeName));
791  forAllConstIters(fSets, iter)
792  {
793  fSetNames.insert(iter.key(), fSetNames.size());
794  }
795  IOobjectList pSets(objects.lookupClass(pointSet::typeName));
796  forAllConstIters(pSets, iter)
797  {
798  pSetNames.insert(iter.key(), pSetNames.size());
799  }
800  }
801 
802  if (cSetNames.size() || fSetNames.size() || pSetNames.size())
803  {
804  // Construct all sets
805  PtrList<cellSet> cellSets(cSetNames.size());
806  PtrList<faceSet> faceSets(fSetNames.size());
807  PtrList<pointSet> pointSets(pSetNames.size());
808 
809  Info<< "Reconstructing sets:" << endl;
810  if (cSetNames.size())
811  {
812  Info<< " cellSets "
813  << cSetNames.sortedToc() << endl;
814  }
815  if (fSetNames.size())
816  {
817  Info<< " faceSets "
818  << fSetNames.sortedToc() << endl;
819  }
820  if (pSetNames.size())
821  {
822  Info<< " pointSets "
823  << pSetNames.sortedToc() << endl;
824  }
825 
826  // Load sets
827  forAll(procMeshes.meshes(), proci)
828  {
829  const fvMesh& procMesh = procMeshes.meshes()[proci];
830 
831  IOobjectList objects
832  (
833  procMesh,
834  databases[0].timeName(),
835  polyMesh::meshSubDir/"sets"
836  );
837 
838  // cellSets
839  const labelList& cellMap =
840  procMeshes.cellProcAddressing()[proci];
841 
842  IOobjectList cSets
843  (
844  objects.lookupClass(cellSet::typeName)
845  );
846 
847  forAllConstIters(cSets, iter)
848  {
849  // Load cellSet
850  const cellSet procSet(*iter());
851  label setI = cSetNames[iter.key()];
852  if (!cellSets.set(setI))
853  {
854  cellSets.set
855  (
856  setI,
857  new cellSet
858  (
859  mesh,
860  iter.key(),
861  procSet.size()
862  )
863  );
864  }
865  cellSet& cSet = cellSets[setI];
866  cSet.instance() = runTime.timeName();
867 
868  for (const label celli : procSet)
869  {
870  cSet.insert(cellMap[celli]);
871  }
872  }
873 
874  // faceSets
875  const labelList& faceMap =
876  procMeshes.faceProcAddressing()[proci];
877 
878  IOobjectList fSets
879  (
880  objects.lookupClass(faceSet::typeName)
881  );
882 
883  forAllConstIters(fSets, iter)
884  {
885  // Load faceSet
886  const faceSet procSet(*iter());
887  label setI = fSetNames[iter.key()];
888  if (!faceSets.set(setI))
889  {
890  faceSets.set
891  (
892  setI,
893  new faceSet
894  (
895  mesh,
896  iter.key(),
897  procSet.size()
898  )
899  );
900  }
901  faceSet& fSet = faceSets[setI];
902  fSet.instance() = runTime.timeName();
903 
904  for (const label facei : procSet)
905  {
906  fSet.insert(mag(faceMap[facei])-1);
907  }
908  }
909  // pointSets
910  const labelList& pointMap =
911  procMeshes.pointProcAddressing()[proci];
912 
913  IOobjectList pSets
914  (
915  objects.lookupClass(pointSet::typeName)
916  );
917  forAllConstIters(pSets, iter)
918  {
919  // Load pointSet
920  const pointSet propSet(*iter());
921  label setI = pSetNames[iter.key()];
922  if (!pointSets.set(setI))
923  {
924  pointSets.set
925  (
926  setI,
927  new pointSet
928  (
929  mesh,
930  iter.key(),
931  propSet.size()
932  )
933  );
934  }
935  pointSet& pSet = pointSets[setI];
936  pSet.instance() = runTime.timeName();
937 
938  for (const label pointi : propSet)
939  {
940  pSet.insert(pointMap[pointi]);
941  }
942  }
943  }
944 
945  // Write sets
946 
947  for (const auto& set : cellSets)
948  {
949  set.write();
950  }
951  for (const auto& set : faceSets)
952  {
953  set.write();
954  }
955  for (const auto& set : pointSets)
956  {
957  set.write();
958  }
959  }
960 
961 
962  // Reconstruct refinement data
963  {
964  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
965 
966  forAll(procMeshes.meshes(), procI)
967  {
968  const fvMesh& procMesh = procMeshes.meshes()[procI];
969 
970  procData.set
971  (
972  procI,
973  new hexRef8Data
974  (
975  IOobject
976  (
977  "dummy",
978  procMesh.time().timeName(),
979  polyMesh::meshSubDir,
980  procMesh,
981  IOobject::READ_IF_PRESENT,
982  IOobject::NO_WRITE,
983  false
984  )
985  )
986  );
987  }
988 
989  // Combine individual parts
990 
991  const PtrList<labelIOList>& cellAddr =
992  procMeshes.cellProcAddressing();
993 
994  UPtrList<const labelList> cellMaps(cellAddr.size());
995  forAll(cellAddr, i)
996  {
997  cellMaps.set(i, &cellAddr[i]);
998  }
999 
1000  const PtrList<labelIOList>& pointAddr =
1001  procMeshes.pointProcAddressing();
1002 
1003  UPtrList<const labelList> pointMaps(pointAddr.size());
1004  forAll(pointAddr, i)
1005  {
1006  pointMaps.set(i, &pointAddr[i]);
1007  }
1008 
1009  UPtrList<const hexRef8Data> procRefs(procData.size());
1010  forAll(procData, i)
1011  {
1012  procRefs.set(i, &procData[i]);
1013  }
1014 
1015  hexRef8Data
1016  (
1017  IOobject
1018  (
1019  "dummy",
1020  mesh.time().timeName(),
1021  polyMesh::meshSubDir,
1022  mesh,
1023  IOobject::NO_READ,
1024  IOobject::NO_WRITE,
1025  false
1026  ),
1027  cellMaps,
1028  pointMaps,
1029  procRefs
1030  ).write();
1031  }
1032  }
1033 
1034 
1035  // Reconstruct refinement data
1036  {
1037  PtrList<hexRef8Data> procData(procMeshes.meshes().size());
1038 
1039  forAll(procMeshes.meshes(), procI)
1040  {
1041  const fvMesh& procMesh = procMeshes.meshes()[procI];
1042 
1043  procData.set
1044  (
1045  procI,
1046  new hexRef8Data
1047  (
1048  IOobject
1049  (
1050  "dummy",
1051  procMesh.time().timeName(),
1052  polyMesh::meshSubDir,
1053  procMesh,
1054  IOobject::READ_IF_PRESENT,
1055  IOobject::NO_WRITE,
1056  false
1057  )
1058  )
1059  );
1060  }
1061 
1062  // Combine individual parts
1063 
1064  const PtrList<labelIOList>& cellAddr =
1065  procMeshes.cellProcAddressing();
1066 
1067  UPtrList<const labelList> cellMaps(cellAddr.size());
1068  forAll(cellAddr, i)
1069  {
1070  cellMaps.set(i, &cellAddr[i]);
1071  }
1072 
1073  const PtrList<labelIOList>& pointAddr =
1074  procMeshes.pointProcAddressing();
1075 
1076  UPtrList<const labelList> pointMaps(pointAddr.size());
1077  forAll(pointAddr, i)
1078  {
1079  pointMaps.set(i, &pointAddr[i]);
1080  }
1081 
1082  UPtrList<const hexRef8Data> procRefs(procData.size());
1083  forAll(procData, i)
1084  {
1085  procRefs.set(i, &procData[i]);
1086  }
1087 
1088  hexRef8Data
1089  (
1090  IOobject
1091  (
1092  "dummy",
1093  mesh.time().timeName(),
1094  polyMesh::meshSubDir,
1095  mesh,
1096  IOobject::NO_READ,
1097  IOobject::NO_WRITE,
1098  false
1099  ),
1100  cellMaps,
1101  pointMaps,
1102  procRefs
1103  ).write();
1104  }
1105 
1106  // If there is a "uniform" directory in the time region
1107  // directory copy from the master processor
1108  {
1109  fileName uniformDir0
1110  (
1111  fileHandler().filePath
1112  (
1113  databases[0].timePath()/regionDir/"uniform"
1114  )
1115  );
1116 
1117  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
1118  {
1119  fileHandler().cp(uniformDir0, runTime.timePath()/regionDir);
1120  }
1121  }
1122 
1123  // For the first region of a multi-region case additionally
1124  // copy the "uniform" directory in the time directory
1125  if (regioni == 0 && regionDir != word::null)
1126  {
1127  fileName uniformDir0
1128  (
1129  fileHandler().filePath
1130  (
1131  databases[0].timePath()/"uniform"
1132  )
1133  );
1134 
1135  if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
1136  {
1137  fileHandler().cp(uniformDir0, runTime.timePath());
1138  }
1139  }
1140  }
1141  }
1142 
1143  Info<< "\nEnd\n" << endl;
1144 
1145  return 0;
1146 }
1147 
1148 
1149 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::sphericalTensor
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
Definition: sphericalTensor.H:54
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
cloudName
const word cloudName(propsDict.get< word >("cloud"))
faCFD.H
regionProperties.H
fvFieldReconstructor.H
lagrangianReconstructor.H
processorFaMeshes.H
Foam::argList::caseName
const fileName & caseName() const noexcept
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:69
Foam::argList::readListIfPresent
bool readListIfPresent(const word &optName, List< T > &list) const
Definition: argListI.H:394
Foam::fileHandler
const fileOperation & fileHandler()
Get current file handler.
Definition: fileOperation.C:1485
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
faMesh.H
processorMeshes.H
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:44
regionNames
wordList regionNames
Definition: getAllRegionOptions.H:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
aMesh
faMesh aMesh(mesh)
checkFaceAddressingComp.H
Foam::fileOperation::readDir
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.
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
pointFieldReconstructor.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fileOperation::nProcs
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
Definition: fileOperation.C:1193
argList.H
faceSet.H
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::argList::path
fileName path() const
Return the full path to the (processor local) case.
Definition: argListI.H:81
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
hexRef8Data.H
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
getAllRegionOptions.H
Priority.
Foam::fileNameList
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:58
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
setRootCase.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
addAllRegionOptions.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::fileOperation::isDir
virtual bool isDir(const fileName &, const bool followLink=true) const =0
Does the name exist as a DIRECTORY in the file system?
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
timeSelector.H
createTime.H
faFieldReconstructor.H
Foam::wordHashSet
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:77
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fileOperation::cp
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
fvCFD.H
cellSet.H
regionDir
const word & regionDir
Definition: findMeshDefinitionDict.H:34
Foam::stringOps::validate
StringType validate(const std::string &str, const UnaryPredicate &accept, const bool invert=false)
Return a copy of the input string with validated characters.
Definition: stringOpsTemplates.C:71
args
Foam::argList args(argc, argv)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::argList::rootPath
const fileName & rootPath() const noexcept
Return root path.
Definition: argListI.H:63
Foam::argList::found
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
pointSet.H