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