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