fvMeshTools.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) 2012-2016 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
27\*---------------------------------------------------------------------------*/
28
29#include "fvMeshTools.H"
30#include "pointSet.H"
31#include "faceSet.H"
32#include "cellSet.H"
33#include "IOobjectList.H"
35#include "processorPolyPatch.H"
38
39// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40
41// Adds patch if not yet there. Returns patchID.
43(
44 fvMesh& mesh,
45 const polyPatch& patch,
46 const dictionary& patchFieldDict,
47 const word& defaultPatchFieldType,
48 const bool validBoundary
49)
50{
51 auto& polyPatches = const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
52
53 label patchi = polyPatches.findPatchID(patch.name());
54 if (patchi != -1)
55 {
56 // Already there
57 return patchi;
58 }
59
60
61 // Append at end unless there are processor patches
62 label insertPatchi = polyPatches.size();
63 label startFacei = mesh.nFaces();
64
65 if (!isA<processorPolyPatch>(patch))
66 {
67 forAll(polyPatches, patchi)
68 {
69 const polyPatch& pp = polyPatches[patchi];
70
71 if (isA<processorPolyPatch>(pp))
72 {
73 insertPatchi = patchi;
74 startFacei = pp.start();
75 break;
76 }
77 }
78 }
79
80
81 // Below is all quite a hack. Feel free to change once there is a better
82 // mechanism to insert and reorder patches.
83
84 // Clear local fields and e.g. polyMesh parallelInfo.
85 mesh.clearOut();
86
87 label sz = polyPatches.size();
88
89 auto& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
90
91 // Add polyPatch at the end
92 polyPatches.resize(sz+1);
93 polyPatches.set
94 (
95 sz,
96 patch.clone
97 (
98 polyPatches,
99 insertPatchi, //index
100 0, //size
101 startFacei //start
102 )
103 );
104 fvPatches.resize(sz+1);
105 fvPatches.set
106 (
107 sz,
109 (
110 polyPatches[sz], // point to newly added polyPatch
111 mesh.boundary()
112 )
113 );
114
115
116 do
117 {
118 #undef doLocalCode
119 #define doLocalCode(FieldType) \
120 { \
121 addPatchFields<FieldType> \
122 ( \
123 mesh, patchFieldDict, defaultPatchFieldType, Zero \
124 ); \
125 }
126
127 // Volume fields
133
134 // Surface fields
140
141 #undef doLocalCode
142 }
143 while (false);
144
145
146 // Create reordering list
147 // - patches before insert position stay as is
148 // - patches after insert position move one up
149 // - appended patch gets moved to insert position
150
151 labelList oldToNew(identity(sz+1));
152 for (label i = insertPatchi; i < sz; ++i)
153 {
154 oldToNew[i] = i+1;
155 }
156 oldToNew[sz] = insertPatchi;
157
158
159 // Shuffle into place
160 polyPatches.reorder(oldToNew, validBoundary);
161 fvPatches.reorder(oldToNew);
162
163 do
164 {
165 #undef doLocalCode
166 #define doLocalCode(FieldType) \
167 { \
168 reorderPatchFields<FieldType>(mesh, oldToNew); \
169 }
170
171 // Volume fields
177
178 // Surface fields
184
185 #undef doLocalCode
186 }
187 while (false);
188
189 return insertPatchi;
190}
191
192
193void Foam::fvMeshTools::setPatchFields
194(
195 fvMesh& mesh,
196 const label patchi,
197 const dictionary& patchFieldDict
198)
199{
200 do
201 {
202 #undef doLocalCode
203 #define doLocalCode(FieldType) \
204 { \
205 setPatchFields<FieldType>(mesh, patchi, patchFieldDict); \
206 }
207
208 // Volume fields
214
215 // Surface fields
221
222 #undef doLocalCode
223 }
224 while (false);
225}
226
227
229{
230 do
231 {
232 #undef doLocalCode
233 #define doLocalCode(FieldType) \
234 { \
235 setPatchFields<FieldType>(mesh, patchi, Zero); \
236 }
237
238 // Volume fields
244
245 // Surface fields
251
252 #undef doLocalCode
253 }
254 while (false);
255}
256
257
258// Deletes last patch
259void Foam::fvMeshTools::trimPatches(fvMesh& mesh, const label nPatches)
260{
261 // Clear local fields and e.g. polyMesh globalMeshData.
262 mesh.clearOut();
263
264 auto& polyPatches = const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
265 auto& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
266
267 if (polyPatches.empty())
268 {
270 << "No patches in mesh"
271 << abort(FatalError);
272 }
273
274 label nFaces = 0;
275 for (label patchi = nPatches; patchi < polyPatches.size(); patchi++)
276 {
277 nFaces += polyPatches[patchi].size();
278 }
279 reduce(nFaces, sumOp<label>());
280
281 if (nFaces)
282 {
284 << "There are still " << nFaces
285 << " faces in " << polyPatches.size()-nPatches
286 << " patches to be deleted" << abort(FatalError);
287 }
288
289 // Remove actual patches
290 polyPatches.resize(nPatches);
291 fvPatches.resize(nPatches);
292
293 do
294 {
295 #undef doLocalCode
296 #define doLocalCode(FieldType) \
297 { \
298 trimPatchFields<FieldType>(mesh, nPatches); \
299 }
300
301 // Volume fields
307
308 // Surface fields
314
315 #undef doLocalCode
316 }
317 while (false);
318}
319
320
322(
323 fvMesh& mesh,
324 const labelList& oldToNew,
325 const label nNewPatches,
326 const bool validBoundary
327)
328{
329 auto& polyPatches = const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
330 auto& fvPatches = const_cast<fvBoundaryMesh&>(mesh.boundary());
331
332 // Shuffle into place
333 polyPatches.reorder(oldToNew, validBoundary);
334 fvPatches.reorder(oldToNew);
335
336 do
337 {
338 #undef doLocalCode
339 #define doLocalCode(FieldType) \
340 { \
341 reorderPatchFields<FieldType>(mesh, oldToNew); \
342 }
343
344 // Volume fields
350
351 // Surface fields
357
358 #undef doLocalCode
359 }
360 while (false);
361
362 // Remove last.
363 trimPatches(mesh, nNewPatches);
364}
365
366
368(
369 fvMesh& mesh,
370 const bool validBoundary
371)
372{
373 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
374
375 labelList newToOld(pbm.size());
376 labelList oldToNew(pbm.size(), -1);
377 label newI = 0;
378
379
380 // Assumes all non-coupled boundaries are on all processors!
381 forAll(pbm, patchI)
382 {
383 const polyPatch& pp = pbm[patchI];
384
385 if (!isA<processorPolyPatch>(pp))
386 {
387 label nFaces = pp.size();
388 if (validBoundary)
389 {
391 }
392
393 if (nFaces > 0)
394 {
395 newToOld[newI] = patchI;
396 oldToNew[patchI] = newI++;
397 }
398 }
399 }
400
401 // Same for processor patches (but need no reduction)
402 forAll(pbm, patchI)
403 {
404 const polyPatch& pp = pbm[patchI];
405
406 if (isA<processorPolyPatch>(pp) && pp.size())
407 {
408 newToOld[newI] = patchI;
409 oldToNew[patchI] = newI++;
410 }
411 }
412
413 newToOld.resize(newI);
414
415 // Move all deletable patches to the end
416 forAll(oldToNew, patchI)
417 {
418 if (oldToNew[patchI] == -1)
419 {
420 oldToNew[patchI] = newI++;
421 }
422 }
423
424 reorderPatches(mesh, oldToNew, newToOld.size(), validBoundary);
425
426 return newToOld;
427}
428
429
431{
432 // Set the fvGeometryScheme to basic since it does not require
433 // any parallel communication to construct the geometry. During
434 // redistributePar one might temporarily end up with processors
435 // with zero procBoundaries. Normally procBoundaries trigger geometry
436 // calculation (e.g. send over cellCentres) so on the processors without
437 // procBoundaries this will not happen. The call to the geometry calculation
438 // is not synchronised and this might lead to a hang for geometry schemes
439 // that do require synchronisation
440
441 tmp<fvGeometryScheme> basicGeometry
442 (
444 );
445 mesh.geometry(basicGeometry);
446}
447
448
451(
452 const IOobject& io,
453 const bool masterOnlyReading,
454 const bool verbose
455)
456{
457 // Region name
458 // ~~~~~~~~~~~
459
460 const fileName meshSubDir
461 (
463 );
464
465
466 fileName facesInstance;
467 fileName pointsInstance;
468
469 // Patch types
470 // ~~~~~~~~~~~
471 // Read and scatter master patches (without reading master mesh!)
472
473 PtrList<entry> patchEntries;
474 if (Pstream::master())
475 {
476 const bool oldParRun = Pstream::parRun(false);
477
478 facesInstance = io.time().findInstance
479 (
480 meshSubDir,
481 "faces",
483 );
484 pointsInstance = io.time().findInstance
485 (
486 meshSubDir,
487 "points",
489 );
490
491 patchEntries = polyBoundaryMeshEntries
492 (
494 (
495 "boundary",
496 facesInstance,
497 meshSubDir,
498 io.db(),
501 false
502 )
503 );
504
505 Pstream::parRun(oldParRun);
506 }
507
508 // Broadcast information to all
510 (
512 patchEntries,
513 facesInstance,
514 pointsInstance
515 );
516
517
518 // Dummy meshes
519 // ~~~~~~~~~~~~
520
521 // Set up to read-if-present.
522 // Note: does not search for mesh so set instance explicitly
523
524 IOobject meshIO(io);
525 meshIO.instance() = facesInstance;
527
528 // For mesh components (points, faces, ...)
529 IOobject cmptIO(meshIO, "points", meshSubDir);
532 cmptIO.registerObject(false);
533
534
535 // Check who has a mesh
536 const fileName meshDir = io.time().path()/facesInstance/meshSubDir;
537 bool haveMesh = isDir(meshDir);
538 if (masterOnlyReading && !Pstream::master())
539 {
540 haveMesh = false;
542 }
543
544 if (!haveMesh)
545 {
547 }
548
549
550 // Read mesh
551 // ~~~~~~~~~
552 // Now all processors use supplied points,faces etc
553 // Note: solution, schemes are also using the supplied IOobject so
554 // on slave will be NO_READ, on master READ_IF_PRESENT. This will
555 // conflict with e.g. timeStampMaster reading so switch off.
556 // Note: v2006 used the READ_IF_PRESENT flag in the meshIO.readOpt(). v2012
557 // (correctly) does no longer so below code explicitly addFvPatches
558 // using the separately read boundary file.
559
560 const auto oldCheckType = IOobject::fileModificationChecking;
562
563
564 // Points
565 cmptIO.instance() = pointsInstance;
566 cmptIO.rename("points");
567 pointIOField points(cmptIO);
568
569 // All other mesh components use facesInstance ...
570 cmptIO.instance() = facesInstance;
571
572 // Faces
573 cmptIO.rename("faces");
574 faceCompactIOList faces(cmptIO);
575
576 // Face owner
577 cmptIO.rename("owner");
578 labelIOList owner(cmptIO);
579
580 // Face neighbour
581 cmptIO.rename("neighbour");
582 labelIOList neighbour(cmptIO);
583
584
586 (
587 meshIO,
588 std::move(points),
589 std::move(faces),
590 std::move(owner),
591 std::move(neighbour)
592 );
593 fvMesh& mesh = *meshPtr;
594
596
597
598 // Some processors without patches? - add patches
599
601 {
602 DynamicList<polyPatch*> newPatches(patchEntries.size());
603
604 if (mesh.boundary().size() == patchEntries.size())
605 {
606 // Assumably we have correctly read the boundary and can clone.
607 // Note: for
608 // v2012 onwards this is probably never the case and this whole
609 // section can be removed.
610 forAll(mesh.boundary(), patchi)
611 {
612 newPatches.append
613 (
614 mesh.boundaryMesh()[patchi].clone(mesh.boundaryMesh()).ptr()
615 );
616 }
617 }
618 else
619 {
620 // Use patchEntries, which were read on master and broadcast
621 label nPatches = 0;
622
623 const bool isEmptyMesh = (mesh.nInternalFaces() == 0);
624
625 forAll(patchEntries, patchi)
626 {
627 const entry& e = patchEntries[patchi];
628 const word type(e.dict().get<word>("type"));
629 const word& name = e.keyword();
630
631 if
632 (
635 )
636 {
637 // Unlikely to work with inter-mixed proc-patches anyhow
638 // - could break out of loop here
639 }
640 else
641 {
642 dictionary patchDict(e.dict());
643
644 if (isEmptyMesh)
645 {
646 patchDict.set("nFaces", 0);
647 patchDict.set("startFace", 0);
648 }
649
650 newPatches.append
651 (
653 (
654 name,
655 patchDict,
656 nPatches++,
658 ).ptr()
659 );
660 }
661 }
662 }
663
665 mesh.addFvPatches(newPatches);
666 }
667
668
669 // Determine zones
670 // ~~~~~~~~~~~~~~~
671
672 wordList pointZoneNames(mesh.pointZones().names());
673 wordList faceZoneNames(mesh.faceZones().names());
674 wordList cellZoneNames(mesh.cellZones().names());
675
677 (
679 pointZoneNames,
680 faceZoneNames,
681 cellZoneNames
682 );
683
684 if (!haveMesh)
685 {
686 // Add the zones. Make sure to remove the old dummy ones first
688 mesh.faceZones().clear();
689 mesh.cellZones().clear();
690
691 List<pointZone*> pz(pointZoneNames.size());
692 forAll(pointZoneNames, i)
693 {
694 pz[i] = new pointZone
695 (
696 pointZoneNames[i],
697 i,
699 );
700 }
701 List<faceZone*> fz(faceZoneNames.size());
702 forAll(faceZoneNames, i)
703 {
704 fz[i] = new faceZone
705 (
706 faceZoneNames[i],
707 i,
709 );
710 }
711 List<cellZone*> cz(cellZoneNames.size());
712 forAll(cellZoneNames, i)
713 {
714 cz[i] = new cellZone
715 (
716 cellZoneNames[i],
717 i,
719 );
720 }
721
722 if (pz.size() || fz.size() || cz.size())
723 {
724 mesh.addZones(pz, fz, cz);
725 }
726 }
727
728 return meshPtr;
729}
730
731
734(
735 const IOobject& io,
736 const bool decompose,
737 const bool verbose
738)
739{
740
741 // Region name
742 // ~~~~~~~~~~~
743
744 const fileName meshSubDir
745 (
747 );
748
749
750 // Patch types
751 // ~~~~~~~~~~~
752 // Read and scatter master patches (without reading master mesh!)
753
754 PtrList<entry> patchEntries;
755 if (Pstream::master())
756 {
757 const bool oldParRun = Pstream::parRun(false);
758
759 patchEntries = polyBoundaryMeshEntries
760 (
762 (
763 "boundary",
764 io.instance(),
765 meshSubDir,
766 io.db(),
769 false
770 )
771 );
772
773 Pstream::parRun(oldParRun);
774 }
775
776 // Broadcast: send patches to all
777 Pstream::broadcast(patchEntries); // == worldComm;
778
779
780 // Dummy meshes
781 // ~~~~~~~~~~~~
782
783 // Check who has or needs a mesh.
784 // For 'decompose', only need mesh on master.
785 // Otherwise check for presence of the "faces" file
786
787 bool haveMesh =
788 (
789 decompose
792 (
793 fileHandler().filePath
794 (
795 io.time().path()/io.instance()/meshSubDir/"faces"
796 )
797 )
798 );
799
800
801 if (!haveMesh)
802 {
803 const bool oldParRun = Pstream::parRun(false);
804
805 // Create dummy mesh - on procs that don't already have a mesh
806 fvMesh dummyMesh
807 (
809 Foam::zero{},
810 false
811 );
812
813 // Add patches
814 polyPatchList patches(patchEntries.size());
815 label nPatches = 0;
816
817 forAll(patchEntries, patchi)
818 {
819 const entry& e = patchEntries[patchi];
820 const word type(e.dict().get<word>("type"));
821 const word& name = e.keyword();
822
823 if
824 (
827 )
828 {
829 // Stop at the first processor patch.
830 // - logic fails with inter-mixed proc-patches anyhow
831 break;
832 }
833 else
834 {
835 dictionary patchDict(e.dict());
836 patchDict.set("nFaces", 0);
837 patchDict.set("startFace", 0);
838
840 (
841 patchi,
843 (
844 name,
845 patchDict,
846 nPatches++,
847 dummyMesh.boundaryMesh()
848 )
849 );
850 }
851 }
853 dummyMesh.addFvPatches(patches, false); // No parallel comms
854
855
856 // Add some dummy zones so upon reading it does not read them
857 // from the undecomposed case. Should be done as extra argument to
858 // regIOobject::readStream?
859
861 (
862 1,
863 new pointZone("dummyPointZone", 0, dummyMesh.pointZones())
864 );
866 (
867 1,
868 new faceZone("dummyFaceZone", 0, dummyMesh.faceZones())
869 );
871 (
872 1,
873 new cellZone("dummyCellZone", 0, dummyMesh.cellZones())
874 );
875 dummyMesh.addZones(pz, fz, cz);
876 dummyMesh.pointZones().clear();
877 dummyMesh.faceZones().clear();
878 dummyMesh.cellZones().clear();
879 //Pout<< "Writing dummy mesh to " << dummyMesh.polyMesh::objectPath()
880 // << endl;
881 dummyMesh.write();
882
883 Pstream::parRun(oldParRun); // Restore parallel state
884 }
885
886
887 // Read mesh
888 // ~~~~~~~~~
889 // Now all processors have a (possibly zero size) mesh so read in
890 // parallel
891
892 //Pout<< "Reading mesh from " << io.objectRelPath() << endl;
894 fvMesh& mesh = *meshPtr;
895
896 // Make sure to use a non-parallel geometry calculation method
898
899
900 // Sync patches
901 // ~~~~~~~~~~~~
902
903 if (!Pstream::master() && haveMesh)
904 {
905 // Check master names against mine
906
908
909 forAll(patchEntries, patchi)
910 {
911 const entry& e = patchEntries[patchi];
912 const word type(e.dict().get<word>("type"));
913 const word& name = e.keyword();
914
915 if
916 (
919 )
920 {
921 break;
922 }
923
924 if (patchi >= patches.size())
925 {
927 << "Non-processor patches not synchronised."
928 << endl
929 << "Processor " << Pstream::myProcNo()
930 << " has only " << patches.size()
931 << " patches, master has "
932 << patchi
933 << exit(FatalError);
934 }
935
936 if
937 (
938 type != patches[patchi].type()
939 || name != patches[patchi].name()
940 )
941 {
943 << "Non-processor patches not synchronised."
944 << endl
945 << "Master patch " << patchi
946 << " name:" << type
947 << " type:" << type << endl
948 << "Processor " << Pstream::myProcNo()
949 << " patch " << patchi
950 << " has name:" << patches[patchi].name()
951 << " type:" << patches[patchi].type()
952 << exit(FatalError);
953 }
954 }
955 }
956
957
958 // Determine zones
959 // ~~~~~~~~~~~~~~~
960
961 wordList pointZoneNames(mesh.pointZones().names());
962 wordList faceZoneNames(mesh.faceZones().names());
963 wordList cellZoneNames(mesh.cellZones().names());
965 (
967 pointZoneNames,
968 faceZoneNames,
969 cellZoneNames
970 );
971
972 if (!haveMesh)
973 {
974 // Add the zones. Make sure to remove the old dummy ones first
976 mesh.faceZones().clear();
977 mesh.cellZones().clear();
978
979 List<pointZone*> pz(pointZoneNames.size());
980 forAll(pointZoneNames, i)
981 {
982 pz[i] = new pointZone(pointZoneNames[i], i, mesh.pointZones());
983 }
984 List<faceZone*> fz(faceZoneNames.size());
985 forAll(faceZoneNames, i)
986 {
987 fz[i] = new faceZone(faceZoneNames[i], i, mesh.faceZones());
988 }
989 List<cellZone*> cz(cellZoneNames.size());
990 forAll(cellZoneNames, i)
991 {
992 cz[i] = new cellZone(cellZoneNames[i], i, mesh.cellZones());
993 }
994 mesh.addZones(pz, fz, cz);
995 }
996
997
998 // Determine sets
999 // ~~~~~~~~~~~~~~
1000
1001 wordList pointSetNames;
1002 wordList faceSetNames;
1003 wordList cellSetNames;
1004 if (Pstream::master())
1005 {
1006 // Read sets
1007 const bool oldParRun = Pstream::parRun(false);
1008 IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1009 Pstream::parRun(oldParRun);
1010
1011 pointSetNames = objects.sortedNames<pointSet>();
1012 faceSetNames = objects.sortedNames<faceSet>();
1013 cellSetNames = objects.sortedNames<cellSet>();
1014 }
1016 (
1018 pointSetNames,
1019 faceSetNames,
1020 cellSetNames
1021 );
1022
1023 if (!haveMesh)
1024 {
1025 for (const word& setName : pointSetNames)
1026 {
1027 pointSet(mesh, setName, 0).write();
1028 }
1029 for (const word& setName : faceSetNames)
1030 {
1031 faceSet(mesh, setName, 0).write();
1032 }
1033 for (const word& setName : cellSetNames)
1034 {
1035 cellSet(mesh, setName, 0).write();
1036 }
1037 }
1038
1039
1040 // Force recreation of globalMeshData.
1041 mesh.globalData();
1042
1043
1044 // Do some checks.
1045
1046 // Check if the boundary definition is unique
1047 // and processor patches are correct
1050
1051 // Check names of zones are equal
1052 mesh.cellZones().checkDefinition(verbose);
1053 mesh.cellZones().checkParallelSync(verbose);
1054 mesh.faceZones().checkDefinition(verbose);
1055 mesh.faceZones().checkParallelSync(verbose);
1056 mesh.pointZones().checkDefinition(verbose);
1058
1059 return meshPtr;
1060}
1061
1062
1064(
1065 const objectRegistry& mesh,
1066 const word& regionName,
1067 const bool verbose
1068)
1069{
1070 // Create dummy system/fv*
1071 {
1072 IOobject io
1073 (
1074 "fvSchemes",
1075 mesh.time().system(),
1076 regionName,
1077 mesh.thisDb(),
1080 false
1081 );
1082
1083 if (!io.typeHeaderOk<IOdictionary>(false))
1084 {
1085 if (verbose)
1086 {
1087 Info<< "Writing dummy " << regionName/io.name() << endl;
1088 }
1090 dict.add("divSchemes", dictionary());
1091 dict.add("gradSchemes", dictionary());
1092 dict.add("laplacianSchemes", dictionary());
1093
1094 IOdictionary(io, dict).regIOobject::write();
1095 }
1096 }
1097 {
1098 IOobject io
1099 (
1100 "fvSolution",
1101 mesh.time().system(),
1102 regionName,
1103 mesh.thisDb(),
1106 false
1107 );
1108
1109 if (!io.typeHeaderOk<IOdictionary>(false))
1110 {
1111 if (verbose)
1112 {
1113 Info<< "Writing dummy " << regionName/io.name() << endl;
1114 }
1116 IOdictionary(io, dict).regIOobject::write();
1117 }
1118 }
1119}
1120
1121
1122// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
wordList sortedNames() const
The sorted names of the IOobjects.
Definition: IOobjectList.C:383
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:506
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
bool registerObject() const noexcept
Should object created with this IOobject be registered?
Definition: IOobjectI.H:107
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:303
virtual void rename(const word &newName)
Rename the object.
Definition: IOobject.H:497
writeOption writeOpt() const noexcept
The write option.
Definition: IOobjectI.H:179
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
static void broadcast(Type &value, const label comm=UPstream::worldComm)
static void broadcasts(const label comm, Type &arg1, Args &&... args)
Broadcast multiple items to all processes in communicator.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
PtrList< T > clone(Args &&... args) const
Make a copy by cloning each of the list elements.
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:103
const word & system() const
Return system name.
Definition: TimePathsI.H:102
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
fileName path() const
Return path.
Definition: Time.H:358
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Definition: Time.C:797
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type get(const label i) const
Definition: UList.H:528
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:293
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
void reorder(const labelUList &oldToNew, const bool check=false)
Definition: UPtrList.C:69
bool empty() const noexcept
True if the list is empty (ie, size() is zero)
Definition: UPtrListI.H:113
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
bool checkDefinition(const bool report=false) const
Check zone definition. Return true if in error.
Definition: ZoneMesh.C:739
void clear()
Clear the zones.
Definition: ZoneMesh.C:730
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:304
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order.
Definition: ZoneMesh.C:758
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Default geometry calculation scheme. Slight stabilisation for bad meshes.
A collection of cell labels.
Definition: cellSet.H:54
A subset of mesh cells.
Definition: cellZone.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:780
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
A list of face labels.
Definition: faceSet.H:54
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:67
A class for handling file names.
Definition: fileName.H:76
virtual bool isFile(const fileName &, const bool checkGzip=true, const bool followLink=true) const =0
Does the name exist as a FILE in the file system?
Foam::fvBoundaryMesh.
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
Definition: fvMeshTools.C:322
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Remove zero sized patches. All but processor patches are.
Definition: fvMeshTools.C:368
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
Definition: fvMeshTools.C:1064
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:43
static autoPtr< fvMesh > newMesh(const IOobject &io, const bool masterOnlyReading, const bool verbose=false)
Read mesh or create dummy mesh (0 cells, >0 patches).
Definition: fvMeshTools.C:451
static void setBasicGeometry(fvMesh &mesh)
Set the fvGeometryScheme to basic (to avoid parallel communication)
Definition: fvMeshTools.C:430
static autoPtr< fvMesh > loadOrCreateMesh(const IOobject &io, const bool decompose, const bool verbose=false)
Definition: fvMeshTools.C:734
static void zeroPatchFields(fvMesh &mesh, const label patchi)
Change patchField to zero on registered fields.
Definition: fvMeshTools.C:228
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:302
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:632
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
void removeFvBoundary()
Definition: fvMesh.C:662
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:239
Registry of regIOobjects.
A set of point labels.
Definition: pointSet.H:54
A subset of mesh points.
Definition: pointZone.H:68
Read and store dictionary entries for boundary patches.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
bool checkDefinition(const bool report=false) const
Check boundary definition.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:498
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:866
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:504
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:1013
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:492
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:324
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:854
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
label nInternalFaces() const noexcept
Number of internal faces.
label nFaces() const noexcept
Number of mesh faces.
int myProcNo() const noexcept
Return processor number.
virtual bool write(const bool valid=true) const
Write using setting from DB.
splitCell * master() const
Definition: splitCell.H:113
virtual const fvGeometryScheme & geometry() const
Return reference to geometry calculation scheme.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
Foam::autoPtr< Foam::dynamicFvMesh > meshPtr
Foam::word regionName(Foam::polyMesh::defaultRegion)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const label nPatches
const pointField & points
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
const fileOperation & fileHandler()
Get current file handler.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:87
GeometricField< tensor, fvsPatchField, surfaceMesh > surfaceTensorField
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:85
GeometricField< sphericalTensor, fvsPatchField, surfaceMesh > surfaceSphericalTensorField
errorManip< error > abort(error &err)
Definition: errorManip.H:144
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:86
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
GeometricField< symmTensor, fvsPatchField, surfaceMesh > surfaceSymmTensorField
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: MSwindows.C:651
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.
#define doLocalCode(GeoField)