Loading...
Searching...
No Matches
polyBoundaryMesh.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) 2018-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
27\*---------------------------------------------------------------------------*/
28
29#include "polyBoundaryMesh.H"
30#include "polyMesh.H"
31#include "primitiveMesh.H"
32#include "processorPolyPatch.H"
33#include "PstreamBuffers.H"
34#include "lduSchedule.H"
35#include "globalMeshData.H"
36#include "wordRes.H"
37#include "DynamicList.H"
38#include "PtrListOps.H"
39#include "edgeHashes.H"
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
46}
47
48
49// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50
51bool Foam::polyBoundaryMesh::hasGroupIDs() const
52{
53 if (groupIDsPtr_)
54 {
55 // Use existing cache
56 return !groupIDsPtr_->empty();
57 }
58
59 const polyPatchList& patches = *this;
60
61 for (const polyPatch& p : patches)
62 {
63 if (!p.inGroups().empty())
64 {
65 return true;
66 }
67 }
68
69 return false;
70}
71
72
73void Foam::polyBoundaryMesh::calcGroupIDs() const
74{
75 if (groupIDsPtr_)
76 {
77 return; // Or FatalError
78 }
79
80 groupIDsPtr_.emplace(16);
81 auto& groupLookup = *groupIDsPtr_;
82
83 const polyPatchList& patches = *this;
84
85 forAll(patches, patchi)
86 {
87 for (const word& groupName : patches[patchi].inGroups())
88 {
89 groupLookup(groupName).push_back(patchi);
90 }
91 }
92
93 // Remove groups that clash with patch names
94 forAll(patches, patchi)
95 {
96 if (groupLookup.empty())
97 {
98 break; // Early termination
99 }
100 else if (groupLookup.erase(patches[patchi].name()))
101 {
103 << "Removed group '" << patches[patchi].name()
104 << "' which clashes with patch " << patchi
105 << " of the same name."
106 << endl;
107 }
108 }
109}
110
111
112void Foam::polyBoundaryMesh::populate(PtrList<entry>&& entries)
113{
114 clearLocalAddressing();
115
116 polyPatchList& patches = *this;
117
118 patches.resize_null(entries.size());
119
120 // Transcribe.
121 // Does not handle nullptr at all (what could possibly be done?)
122 forAll(patches, patchi)
123 {
125 (
126 patchi,
128 (
129 entries[patchi].keyword(),
130 entries[patchi].dict(),
131 patchi,
132 *this
133 )
134 );
135 }
136
137 entries.clear();
138}
139
140
141bool Foam::polyBoundaryMesh::readIOcontents(const bool allowOptionalRead)
142{
143 bool updated = false;
144 PtrList<entry> entries;
145
146 if
147 (
148 this->isReadRequired()
149 || (allowOptionalRead && this->isReadOptional() && this->headerOk())
150 )
151 {
152 // Warn for MUST_READ_IF_MODIFIED
153 warnNoRereading<polyBoundaryMesh>();
154
155 // Read entries
156 Istream& is = readStream(typeName);
157
158 is >> entries;
159
160 is.check(FUNCTION_NAME);
161 close();
162 updated = true;
163 }
164 // Future: support master-only and broadcast?
165
166 if (updated)
167 {
168 populate(std::move(entries));
169 }
171 return updated;
172}
173
174
175// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
176
178(
179 const IOobject& io,
180 const polyMesh& mesh
181)
182:
185 mesh_(mesh)
186{
187 readIOcontents(false); // allowOptionalRead = false
188}
189
190
192(
193 const IOobject& io,
194 const polyMesh& pm,
196)
200 mesh_(pm)
201{}
202
203
205(
206 const IOobject& io,
207 const polyMesh& pm,
208 const label size
209)
213 mesh_(pm)
214{}
215
216
218(
219 const IOobject& io,
220 const polyMesh& pm,
221 const polyPatchList& list
222)
223:
226 mesh_(pm)
227{
228 if (!readIOcontents(true)) // allowOptionalRead = true
229 {
230 // Nothing read. Use supplied patches
231 polyPatchList& patches = *this;
232 patches.resize(list.size());
233
234 forAll(patches, patchi)
236 patches.set(patchi, list[patchi].clone(*this));
237 }
238 }
239}
240
241
243(
244 const IOobject& io,
245 const polyMesh& pm,
246 PtrList<entry>&& entries
247)
248:
251 mesh_(pm)
252{
253 if (!readIOcontents(true)) // allowOptionalRead = true
254 {
255 populate(std::move(entries));
257 entries.clear();
258}
259
260
261// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
262
272 polyPatchList& patches = *this;
273
274 for (polyPatch& p : patches)
275 {
276 p.clearGeom();
277 }
278}
279
280
281// Private until it is more generally required (and gets a better name?)
282void Foam::polyBoundaryMesh::clearLocalAddressing()
284 neighbourEdgesPtr_.reset(nullptr);
285 patchIDPtr_.reset(nullptr);
286 groupIDsPtr_.reset(nullptr);
287}
288
289
291{
292 clearLocalAddressing();
293
294 polyPatchList& patches = *this;
295
296 for (polyPatch& p : patches)
297 {
298 p.clearAddressing();
299 }
300}
301
302
303// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
304
305void Foam::polyBoundaryMesh::calcGeometry()
306{
307 // Make sure messages don't interact by having unique tag
308 const int oldTag = UPstream::incrMsgType();
309
311
312 if
313 (
314 pBufs.commsType() == Pstream::commsTypes::buffered
315 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
316 )
317 {
318 forAll(*this, patchi)
319 {
320 operator[](patchi).initGeometry(pBufs);
321 }
322
323 pBufs.finishedSends();
324
325 forAll(*this, patchi)
326 {
327 operator[](patchi).calcGeometry(pBufs);
328 }
329 }
330 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
331 {
332 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
333
334 // Dummy.
335 pBufs.finishedSends();
336
337 for (const auto& schedEval : patchSchedule)
338 {
339 const label patchi = schedEval.patch;
340
341 if (schedEval.init)
342 {
343 operator[](patchi).initGeometry(pBufs);
344 }
345 else
346 {
347 operator[](patchi).calcGeometry(pBufs);
348 }
349 }
351
352 // Reset tag
353 UPstream::msgType(oldTag);
354}
355
356
358{
359 return faceList::subList
360 (
361 mesh_.faces(),
362 mesh_.nBoundaryFaces(),
363 mesh_.nInternalFaces()
364 );
365}
366
367
369{
370 return labelList::subList
371 (
372 mesh_.faceOwner(),
373 mesh_.nBoundaryFaces(),
374 mesh_.nInternalFaces()
375 );
376}
377
378// Potentially useful to simplify logic elsewhere?
379// const Foam::labelList::subList Foam::polyBoundaryMesh::faceNeighbour() const
380// {
381// return labelList::subList();
382// }
383
384
387{
388 const polyPatchList& patches = *this;
389
391
392 forAll(patches, patchi)
393 {
394 list.set(patchi, &patches[patchi].faceCells());
396
397 return list;
398}
399
400
403{
404 if (Pstream::parRun())
405 {
407 << "Neighbour edge addressing not correct across parallel"
408 << " boundaries." << endl;
409 }
410
411 if (!neighbourEdgesPtr_)
412 {
413 neighbourEdgesPtr_.emplace(size());
414 auto& neighbourEdges = *neighbourEdgesPtr_;
415
416 // Initialize.
417 label nEdgePairs = 0;
418 forAll(*this, patchi)
419 {
420 const polyPatch& pp = operator[](patchi);
421
422 neighbourEdges[patchi].setSize(pp.nEdges() - pp.nInternalEdges());
423
424 for (labelPair& edgeInfo : neighbourEdges[patchi])
425 {
426 edgeInfo[0] = -1;
427 edgeInfo[1] = -1;
428 }
429
430 nEdgePairs += pp.nEdges() - pp.nInternalEdges();
431 }
432
433 // From mesh edge (expressed as a point pair so as not to construct
434 // point addressing) to patch + relative edge index.
435 EdgeMap<labelPair> pointsToEdge(nEdgePairs);
436
437 forAll(*this, patchi)
438 {
439 const polyPatch& pp = operator[](patchi);
440
441 const edgeList& edges = pp.edges();
442
443 for
444 (
445 label edgei = pp.nInternalEdges();
446 edgei < edges.size();
447 edgei++
448 )
449 {
450 // Edge in patch local points
451 const edge& e = edges[edgei];
452
453 // Edge in mesh points.
454 edge meshEdge(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
455
456 auto fnd = pointsToEdge.find(meshEdge);
457
458 if (!fnd.good())
459 {
460 // First occurrence of mesh edge. Store patch and my
461 // local index.
462 pointsToEdge.insert
463 (
464 meshEdge,
466 (
467 patchi,
468 edgei - pp.nInternalEdges()
469 )
470 );
471 }
472 else
473 {
474 // Second occurrence. Store.
475 const labelPair& edgeInfo = fnd.val();
476
477 neighbourEdges[patchi][edgei - pp.nInternalEdges()] =
478 edgeInfo;
479
480 neighbourEdges[edgeInfo[0]][edgeInfo[1]]
481 = labelPair(patchi, edgei - pp.nInternalEdges());
482
483 // Found all two occurrences of this edge so remove from
484 // hash to save space. Note that this will give lots of
485 // problems if the polyBoundaryMesh is multiply connected.
486 pointsToEdge.erase(meshEdge);
487 }
488 }
489 }
490
491 if (pointsToEdge.size())
492 {
494 << "Not all boundary edges of patches match up." << nl
495 << "Is the outside of your mesh multiply connected?"
496 << abort(FatalError);
497 }
498
499 forAll(*this, patchi)
500 {
501 const polyPatch& pp = operator[](patchi);
502
503 const labelPairList& nbrEdges = neighbourEdges[patchi];
504
505 forAll(nbrEdges, i)
506 {
507 const labelPair& edgeInfo = nbrEdges[i];
508
509 if (edgeInfo[0] == -1 || edgeInfo[1] == -1)
510 {
511 const label edgei = pp.nInternalEdges() + i;
512 const edge& e = pp.edges()[edgei];
513
515 << "Not all boundary edges of patches match up." << nl
516 << "Edge " << edgei << " on patch " << pp.name()
517 << " end points " << pp.localPoints()[e[0]] << ' '
518 << pp.localPoints()[e[1]] << " is not matched to an"
519 << " edge on any other patch." << nl
520 << "Is the outside of your mesh multiply connected?"
521 << abort(FatalError);
522 }
523 }
525 }
526
527 return *neighbourEdgesPtr_;
528}
529
530
532{
533 if (!patchIDPtr_)
534 {
535 patchIDPtr_.emplace(mesh_.nBoundaryFaces());
536 auto& list = *patchIDPtr_;
537
538 const polyPatchList& patches = *this;
539
540 forAll(patches, patchi)
541 {
543 (
544 list,
545 patches[patchi].size(),
546 (patches[patchi].start() - mesh_.nInternalFaces())
547 ) = patchi;
549 }
550
551 return *patchIDPtr_;
552}
553
554
555Foam::label Foam::polyBoundaryMesh::patchID(const label meshFacei) const
556{
557 const label bndFacei = (meshFacei - mesh_.nInternalFaces());
558
559 return
560 (
561 (bndFacei >= 0 && bndFacei < mesh_.nBoundaryFaces())
562 ? this->patchID()[bndFacei]
563 : -1
564 );
565}
566
567
569Foam::polyBoundaryMesh::patchID(const labelUList& meshFaceIndices) const
570{
571 labelList output(meshFaceIndices.size());
572 forAll(meshFaceIndices, i)
573 {
574 output[i] = patchID(meshFaceIndices[i]);
575 }
576 return output;
577}
578
579
582{
583 if (!groupIDsPtr_)
584 {
585 calcGroupIDs();
586 }
587
588 return *groupIDsPtr_;
589}
590
591
593(
594 const word& groupName,
595 const labelUList& patchIDs
596)
597{
598 groupIDsPtr_.reset(nullptr);
599
600 polyPatchList& patches = *this;
601
602 boolList pending(patches.size(), true);
603
604 // Add to specified patches
605 for (const label patchi : patchIDs)
606 {
607 if (pending.test(patchi))
608 {
609 pending.unset(patchi);
610 patches[patchi].addGroup(groupName);
611 }
612 }
613
614 // Remove from other patches
615 forAll(patches, patchi)
616 {
617 if (pending.test(patchi))
619 patches[patchi].removeGroup(groupName);
620 }
621 }
622}
623
624
626{
627 const polyPatchList& patches = *this;
628
629 label count = 0;
630
631 for (const polyPatch& p : patches)
632 {
634 {
635 break;
636 }
637
639 }
640
641 return count;
642}
643
644
646{
647 const polyPatchList& patches = *this;
648
649 label count = 0;
650
651 for (const polyPatch& p : patches)
652 {
654 {
655 ++count;
657 }
658
659 return count;
660}
661
662
664{
665 const polyPatchList& patches = *this;
666
667 label count = 0;
668
669 for (const polyPatch& p : patches)
670 {
672 {
673 break;
674 }
675
676 count += p.nFaces();
677 }
678
679 return count;
680}
681
686}
687
697 return
700 *this,
701 [](const polyPatch& p) { return p.physicalType(); }
702 );
703}
704
705
707{
708 return
711 *this,
712 [](const polyPatch& p) { return p.start(); }
713 );
714}
715
716
718{
719 return
722 *this,
723 [](const polyPatch& p) { return p.size(); }
724 );
725}
726
727
729{
730 return
733 *this,
734 [](const polyPatch& p) { return p.range(); }
735 );
736}
737
740{
741 return this->groupPatchIDs().sortedToc();
742}
743
745Foam::label Foam::polyBoundaryMesh::start() const noexcept
746{
747 return mesh_.nInternalFaces();
748}
749
752{
753 return mesh_.nBoundaryFaces();
754}
755
758{
759 return labelRange(mesh_.nInternalFaces(), mesh_.nBoundaryFaces());
760}
761
762
763Foam::labelRange Foam::polyBoundaryMesh::range(const label patchi) const
764{
765 if (patchi < 0)
766 {
767 return labelRange(mesh_.nInternalFaces(), 0);
769
770 // Will fail if patchi >= size()
771 return (*this)[patchi].range();
772}
773
774
776(
777 const wordRe& matcher,
778 const bool useGroups
779) const
780{
781 if (matcher.empty())
782 {
783 return labelList();
784 }
785
786 // Only check groups if requested and they exist
787 const bool checkGroups = (useGroups && this->hasGroupIDs());
788
789 labelHashSet ids;
790
791 if (matcher.isPattern())
792 {
793 if (checkGroups)
794 {
795 const auto& groupLookup = groupPatchIDs();
796 forAllConstIters(groupLookup, iter)
797 {
798 if (matcher(iter.key()))
799 {
800 // Add patch ids associated with the group
801 ids.insert(iter.val());
802 }
803 }
804 }
805
806 if (ids.empty())
807 {
808 return PtrListOps::findMatching(*this, matcher);
809 }
810 else
811 {
812 ids.insert(PtrListOps::findMatching(*this, matcher));
813 }
814 }
815 else
816 {
817 // Literal string.
818 // Special version of above for reduced memory footprint.
819
820 const label patchId = PtrListOps::firstMatching(*this, matcher);
821
822 if (patchId >= 0)
823 {
824 return labelList(one{}, patchId);
825 }
826 else if (checkGroups)
827 {
828 const auto iter = groupPatchIDs().cfind(matcher);
829
830 if (iter.good())
831 {
832 // Hash ids associated with the group
833 ids.insert(iter.val());
834 }
836 }
837
838 return ids.sortedToc();
839}
840
841
843(
844 const wordRes& matcher,
845 const bool useGroups
846) const
847{
848 if (matcher.empty())
849 {
850 return labelList();
851 }
852 else if (matcher.size() == 1)
853 {
854 return this->indices(matcher.front(), useGroups);
855 }
856
857 labelHashSet ids;
858
859 // Only check groups if requested and they exist
860 if (useGroups && this->hasGroupIDs())
861 {
862 ids.reserve(this->size());
863
864 const auto& groupLookup = groupPatchIDs();
865 forAllConstIters(groupLookup, iter)
866 {
867 if (matcher(iter.key()))
868 {
869 // Add patch ids associated with the group
870 ids.insert(iter.val());
871 }
872 }
873 }
874
875 if (ids.empty())
876 {
877 return PtrListOps::findMatching(*this, matcher);
878 }
879 else
880 {
881 ids.insert(PtrListOps::findMatching(*this, matcher));
882 }
883
884 return ids.sortedToc();
885}
886
887
889(
890 const wordRes& allow,
891 const wordRes& deny,
892 const bool useGroups
893) const
894{
895 if (allow.empty() && deny.empty())
896 {
897 // Fast-path: select all
898 return identity(this->size());
899 }
900
901 const wordRes::filter matcher(allow, deny);
902
903 labelHashSet ids;
904
905 // Only check groups if requested and they exist
906 if (useGroups && this->hasGroupIDs())
907 {
908 ids.reserve(this->size());
909
910 const auto& groupLookup = groupPatchIDs();
911 forAllConstIters(groupLookup, iter)
912 {
913 if (matcher(iter.key()))
914 {
915 // Add patch ids associated with the group
916 ids.insert(iter.val());
917 }
918 }
919 }
920
921 if (ids.empty())
922 {
923 return PtrListOps::findMatching(*this, matcher);
924 }
925 else
926 {
927 ids.insert(PtrListOps::findMatching(*this, matcher));
928 }
929
930 return ids.sortedToc();
931}
932
933
934Foam::label Foam::polyBoundaryMesh::findIndex(const wordRe& key) const
935{
936 if (key.empty())
938 return -1;
939 }
940 return PtrListOps::firstMatching(*this, key);
941}
942
943
945(
946 const word& patchName,
947 bool allowNotFound
948) const
949{
950 if (patchName.empty())
951 {
952 return -1;
953 }
954
955 const label patchId = PtrListOps::firstMatching(*this, patchName);
956
957 if (patchId >= 0)
958 {
959 return patchId;
960 }
961
962 if (!allowNotFound)
963 {
965 << "Patch '" << patchName << "' not found. "
966 << "Available patch names";
967
968 if (polyMesh::defaultRegion != mesh_.name())
969 {
971 << " in region '" << mesh_.name() << "'";
972 }
973
975 << " include: " << names() << endl
976 << exit(FatalError);
977 }
978
979 // Patch not found
980 if (debug)
981 {
982 Pout<< "label polyBoundaryMesh::findPatchID(const word&) const"
983 << "Patch named " << patchName << " not found. "
984 << "Available patch names: " << names() << endl;
985 }
987 // Not found, return -1
988 return -1;
989}
990
991
992const Foam::polyPatch*
993Foam::polyBoundaryMesh::cfindPatch(const word& patchName) const
994{
995 const polyPatchList& patches = *this;
996
997 if (!patchName.empty())
998 {
999 // Note: get() handles out-of-range access properly
1000 return patches.get(PtrListOps::firstMatching(patches, patchName));
1002
1003 return nullptr;
1004}
1005
1006
1008Foam::polyBoundaryMesh::whichPatchFace(const label meshFacei) const
1009{
1010 if (meshFacei < mesh().nInternalFaces())
1011 {
1012 // Internal face: return (-1, meshFace)
1013 return labelPair(-1, meshFacei);
1014 }
1015 else if (meshFacei >= mesh().nFaces())
1016 {
1017 // Bounds error: abort
1019 << "Face " << meshFacei
1020 << " out of bounds. Number of geometric faces " << mesh().nFaces()
1021 << abort(FatalError);
1022
1023 return labelPair(-1, meshFacei);
1024 }
1025
1026
1027 const polyPatchList& patches = *this;
1028
1029 // Do we have cached patch info?
1030 if (patchIDPtr_)
1031 {
1032 const label patchi =
1033 this->patchID()[meshFacei - mesh().nInternalFaces()];
1034
1035 // (patch, local face index)
1036 return labelPair(patchi, meshFacei - patches[patchi].start());
1037 }
1038
1039
1040 // Patches are ordered, use binary search
1041 // Find out which patch index and local patch face the specified
1042 // mesh face belongs to by comparing label with patch start labels.
1043
1044 const label patchi =
1045 findLower
1046 (
1047 patches,
1048 meshFacei,
1049 0,
1050 // Must include the start in the comparison
1051 [](const polyPatch& p, label val) { return (p.start() <= val); }
1052 );
1053
1054 if (patchi < 0 || !patches[patchi].range().contains(meshFacei))
1055 {
1056 // If not in any of above, it is trouble!
1058 << "Face " << meshFacei << " not found in any of the patches "
1059 << flatOutput(names()) << nl
1060 << "The patches appear to be inconsistent with the mesh :"
1061 << " internalFaces:" << mesh().nInternalFaces()
1062 << " total number of faces:" << mesh().nFaces()
1063 << abort(FatalError);
1064
1065 return labelPair(-1, meshFacei);
1066 }
1068 // (patch, local face index)
1069 return labelPair(patchi, meshFacei - patches[patchi].start());
1070}
1071
1072
1074Foam::polyBoundaryMesh::whichPatchFace(const labelUList& meshFaceIndices) const
1075{
1076 labelPairList output(meshFaceIndices.size());
1077 forAll(meshFaceIndices, i)
1079 output[i] = whichPatchFace(meshFaceIndices[i]);
1080 }
1081 return output;
1082}
1083
1084
1086(
1087 const UList<wordRe>& select,
1088 const bool warnNotFound,
1089 const bool useGroups
1090) const
1091{
1092 labelHashSet ids;
1093 if (select.empty())
1094 {
1095 return ids;
1096 }
1097
1098 const polyPatchList& patches = *this;
1099
1100 const label len = patches.size();
1101
1102 ids.reserve(len);
1103
1104 // Only check groups if requested and they exist
1105 const bool checkGroups = (useGroups && this->hasGroupIDs());
1106
1107 for (const wordRe& matcher : select)
1108 {
1109 bool missed = true;
1110
1111 for (label i = 0; i < len; ++i)
1112 {
1113 if (matcher(patches[i].name()))
1114 {
1115 ids.insert(i);
1116 missed = false;
1117 }
1118 }
1119
1120 if (missed && checkGroups)
1121 {
1122 // Check group names
1123 if (matcher.isPattern())
1124 {
1125 forAllConstIters(groupPatchIDs(), iter)
1126 {
1127 if (matcher.match(iter.key()))
1128 {
1129 // Hash ids associated with the group
1130 ids.insert(iter.val());
1131 missed = false;
1132 }
1133 }
1134 }
1135 else
1136 {
1137 const auto iter = groupPatchIDs().cfind(matcher);
1138
1139 if (iter.good())
1140 {
1141 // Hash ids associated with the group
1142 ids.insert(iter.val());
1143 missed = false;
1144 }
1145 }
1146 }
1147
1148 if (missed && warnNotFound)
1149 {
1150 if (checkGroups)
1151 {
1153 << "Cannot find any patch or group names matching "
1154 << matcher << endl;
1155 }
1156 else
1157 {
1159 << "Cannot find any patch names matching "
1160 << matcher << endl;
1161 }
1163 }
1164
1165 return ids;
1166}
1167
1168
1170(
1171 const labelUList& patchIDs,
1172 wordList& groups,
1173 labelHashSet& nonGroupPatches
1174) const
1175{
1176 // Current matched groups
1177 DynamicList<word> matchedGroups(1);
1178
1179 // Current set of unmatched patches
1180 nonGroupPatches = labelHashSet(patchIDs);
1181
1182 const HashTable<labelList>& groupLookup = this->groupPatchIDs();
1183 forAllConstIters(groupLookup, iter)
1184 {
1185 // Store currently unmatched patches so we can restore
1186 labelHashSet oldNonGroupPatches(nonGroupPatches);
1187
1188 // Match by deleting patches in group from the current set and seeing
1189 // if all have been deleted.
1190 labelHashSet groupPatchSet(iter.val());
1191
1192 label nMatch = nonGroupPatches.erase(groupPatchSet);
1193
1194 if (nMatch == groupPatchSet.size())
1195 {
1196 matchedGroups.push_back(iter.key());
1197 }
1198 else if (nMatch != 0)
1199 {
1200 // No full match. Undo.
1201 nonGroupPatches.transfer(oldNonGroupPatches);
1203 }
1204
1205 groups.transfer(matchedGroups);
1206}
1207
1208
1209bool Foam::polyBoundaryMesh::checkParallelSync(const bool report) const
1210{
1211 if (!Pstream::parRun())
1212 {
1213 return false;
1214 }
1215
1216 const polyBoundaryMesh& bm = *this;
1217
1218 bool hasError = false;
1219
1220 // Collect non-proc patches and check proc patches are last.
1221 wordList localNames(bm.size());
1222 wordList localTypes(bm.size());
1223
1224 label nonProci = 0;
1225
1226 forAll(bm, patchi)
1227 {
1228 if (!isA<processorPolyPatch>(bm[patchi]))
1229 {
1230 if (nonProci != patchi)
1231 {
1232 // A processor patch in between normal patches!
1233 hasError = true;
1234
1235 if (debug || report)
1236 {
1237 Pout<< " ***Problem with boundary patch " << patchi
1238 << " name:" << bm[patchi].name()
1239 << " type:" << bm[patchi].type()
1240 << " - seems to be preceeded by processor patches."
1241 << " This is usually a problem." << endl;
1242 }
1243 }
1244 else
1245 {
1246 localNames[nonProci] = bm[patchi].name();
1247 localTypes[nonProci] = bm[patchi].type();
1248 ++nonProci;
1249 }
1250 }
1251 }
1252 localNames.resize(nonProci);
1253 localTypes.resize(nonProci);
1254
1255 // Check and report error(s) on master
1256 // - don't need indexing on master itself
1257
1258 const globalIndex procAddr(globalIndex::gatherNonLocal{}, nonProci);
1259
1260 const wordList allNames(procAddr.gather(localNames));
1261 const wordList allTypes(procAddr.gather(localTypes));
1262
1263 // Automatically restricted to master
1264 for (const int proci : procAddr.subProcs())
1265 {
1266 const auto procNames(allNames.slice(procAddr.range(proci)));
1267 const auto procTypes(allTypes.slice(procAddr.range(proci)));
1268
1269 if (procNames != localNames || procTypes != localTypes)
1270 {
1271 hasError = true;
1272
1273 if (debug || report)
1274 {
1275 Info<< " ***Inconsistent patches across processors, "
1276 "processor0 has patch names:" << localNames
1277 << " patch types:" << localTypes
1278 << " processor" << proci
1279 << " has patch names:" << procNames
1280 << " patch types:" << procTypes
1281 << endl;
1282 }
1283 }
1285
1286 // Reduce (not broadcast) to respect local out-of-order errors (first loop)
1287 return returnReduceOr(hasError);
1288}
1289
1290
1291bool Foam::polyBoundaryMesh::checkDefinition(const bool report) const
1292{
1293 label nextPatchStart = mesh().nInternalFaces();
1294 const polyBoundaryMesh& bm = *this;
1295
1296 bool hasError = false;
1297
1298 wordHashSet patchNames(2*this->size());
1299
1300 forAll(bm, patchi)
1301 {
1302 if (bm[patchi].start() != nextPatchStart && !hasError)
1303 {
1304 hasError = true;
1305
1306 Info<< " ****Problem with boundary patch " << patchi
1307 << " named " << bm[patchi].name()
1308 << " of type " << bm[patchi].type()
1309 << ". The patch should start on face no " << nextPatchStart
1310 << " and the patch specifies " << bm[patchi].start()
1311 << "." << endl
1312 << "Possibly consecutive patches have this same problem."
1313 << " Suppressing future warnings." << endl;
1314 }
1315
1316 if (!patchNames.insert(bm[patchi].name()) && !hasError)
1317 {
1318 hasError = true;
1319
1320 Info<< " ****Duplicate boundary patch " << patchi
1321 << " named " << bm[patchi].name()
1322 << " of type " << bm[patchi].type()
1323 << "." << endl
1324 << "Suppressing future warnings." << endl;
1325 }
1326
1327 nextPatchStart += bm[patchi].size();
1328 }
1329
1330 Pstream::reduceOr(hasError);
1331
1332 if (debug || report)
1333 {
1334 if (hasError)
1335 {
1336 Pout<< " ***Boundary definition is in error." << endl;
1337 }
1338 else
1339 {
1340 Info<< " Boundary definition OK." << endl;
1342 }
1343
1344 return hasError;
1345}
1346
1347
1349{
1351
1352 if
1353 (
1354 pBufs.commsType() == Pstream::commsTypes::buffered
1355 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
1356 )
1357 {
1358 forAll(*this, patchi)
1359 {
1360 operator[](patchi).initMovePoints(pBufs, p);
1361 }
1362
1363 pBufs.finishedSends();
1364
1365 forAll(*this, patchi)
1366 {
1367 operator[](patchi).movePoints(pBufs, p);
1368 }
1369 }
1370 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1371 {
1372 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1373
1374 // Dummy.
1375 pBufs.finishedSends();
1376
1377 for (const auto& schedEval : patchSchedule)
1378 {
1379 const label patchi = schedEval.patch;
1380
1381 if (schedEval.init)
1382 {
1383 operator[](patchi).initMovePoints(pBufs, p);
1384 }
1385 else
1386 {
1387 operator[](patchi).movePoints(pBufs, p);
1388 }
1389 }
1390 }
1391}
1392
1393
1395{
1396 neighbourEdgesPtr_.reset(nullptr);
1397 patchIDPtr_.reset(nullptr);
1398 groupIDsPtr_.reset(nullptr);
1399
1401
1402 if
1403 (
1404 pBufs.commsType() == Pstream::commsTypes::buffered
1405 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
1406 )
1407 {
1408 forAll(*this, patchi)
1409 {
1410 operator[](patchi).initUpdateMesh(pBufs);
1411 }
1412
1413 pBufs.finishedSends();
1414
1415 forAll(*this, patchi)
1416 {
1417 operator[](patchi).updateMesh(pBufs);
1418 }
1419 }
1420 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1421 {
1422 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1423
1424 // Dummy.
1425 pBufs.finishedSends();
1426
1427 for (const auto& schedEval : patchSchedule)
1428 {
1429 const label patchi = schedEval.patch;
1430
1431 if (schedEval.init)
1432 {
1433 operator[](patchi).initUpdateMesh(pBufs);
1434 }
1435 else
1436 {
1437 operator[](patchi).updateMesh(pBufs);
1438 }
1439 }
1440 }
1441}
1442
1443
1445(
1446 const labelUList& oldToNew,
1447 const bool validBoundary
1448)
1449{
1450 // Change order of patches
1451 polyPatchList::reorder(oldToNew);
1452
1453 // Adapt indices
1454 polyPatchList& patches = *this;
1455
1456 forAll(patches, patchi)
1457 {
1458 patches[patchi].index() = patchi;
1459 }
1460
1461 // Clear group-to-patch addressing. Note: could re-calculate
1462 groupIDsPtr_.reset(nullptr);
1463
1464 if (validBoundary)
1465 {
1466 updateMesh();
1467 }
1468}
1469
1470
1471void Foam::polyBoundaryMesh::writeEntry(Ostream& os) const
1472{
1473 const polyPatchList& entries = *this;
1474
1475 os << entries.size();
1476
1477 if (entries.empty())
1478 {
1479 // 0-sized : can write with less vertical space
1481 }
1482 else
1483 {
1484 os << nl << token::BEGIN_LIST << incrIndent << nl;
1485
1486 for (const auto& pp : entries)
1487 {
1488 os.beginBlock(pp.name());
1489 os << pp;
1490 os.endBlock();
1493 }
1494 os.check(FUNCTION_NAME);
1495}
1496
1497
1499(
1500 const word& keyword,
1501 Ostream& os
1502) const
1503{
1504 const polyPatchList& entries = *this;
1505
1506 if (!keyword.empty())
1507 {
1508 os.write(keyword);
1509 os << (entries.empty() ? token::SPACE : token::NL);
1510 }
1512 writeEntry(os);
1513
1514 if (!keyword.empty()) os.endEntry();
1515}
1516
1517
1519{
1520 writeEntry(os);
1521 return os.good();
1522}
1523
1524
1526(
1527 IOstreamOption streamOpt,
1528 const bool writeOnProc
1529) const
1530{
1532 return regIOobject::writeObject(streamOpt, writeOnProc);
1533}
1534
1535
1536// * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
1537
1538const Foam::polyPatch& Foam::polyBoundaryMesh::operator[]
1539(
1540 const word& patchName
1541) const
1542{
1543 const label patchi = findPatchID(patchName);
1544
1545 if (patchi < 0)
1546 {
1548 << "Patch named " << patchName << " not found." << nl
1549 << "Available patch names: " << names() << endl
1551 }
1552
1553 return operator[](patchi);
1554}
1555
1556
1557Foam::polyPatch& Foam::polyBoundaryMesh::operator[]
1558(
1559 const word& patchName
1560)
1561{
1562 const label patchi = findPatchID(patchName);
1563
1564 if (patchi < 0)
1565 {
1567 << "Patch named " << patchName << " not found." << nl
1568 << "Available patch names: " << names() << endl
1569 << abort(FatalError);
1570 }
1572 return operator[](patchi);
1573}
1574
1575
1576// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1577
1578Foam::Ostream& Foam::operator<<(Ostream& os, const polyBoundaryMesh& pbm)
1579{
1580 pbm.writeData(os);
1581 return os;
1582}
1583
1584
1585// ************************************************************************* //
scalar range
Functions to operate on Pointer Lists.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void push_back(const T &val)
Copy append an element to the end of this list.
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
Foam::label erase(InputIter first, InputIter last)
Definition HashTable.C:515
bool empty() const noexcept
True if the hash table is empty.
Definition HashTable.H:353
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
Definition HashTable.C:794
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition HashTable.C:489
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
autoPtr< IOobject > clone() const
Clone.
Definition IOobject.H:641
A simple container for options an IOstream can normally have.
compressionType compression() const noexcept
Get the stream compression.
@ UNCOMPRESSED
compression = false
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
SubList< face > subList
Definition List.H:129
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual Ostream & endBlock()
Write end block group.
Definition Ostream.C:108
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition Ostream.C:90
label nEdges() const
Number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
UPstream::commsTypes commsType() const noexcept
The communications type of the stream.
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void resize_null(const label newLen)
Set the addressed list to the given size, deleting all existing entries. Afterwards the list contains...
Definition PtrListI.H:113
constexpr PtrList() noexcept
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool test(const label i) const
Test bool value at specified position, always false for out-of-range access.
Definition UList.H:852
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
Definition SubList.H:258
T & front()
Access first element of the list, position [0].
Definition UListI.H:239
bool unset(const label i)
Unset the bool entry at specified position, always false for out-of-range access.
Definition UList.H:885
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
@ scheduled
"scheduled" (MPI standard) : (MPI_Send, MPI_Recv)
Definition UPstream.H:83
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static void reduceOr(bool &value, const int communicator=worldComm)
Logical (or) reduction (MPI_AllReduce).
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
static commsTypes defaultCommsType
Default commsType.
Definition UPstream.H:1045
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
void reorder(const labelUList &oldToNew, const bool check=false)
Reorder elements. Reordering must be unique (ie, shuffle).
Definition UPtrList.C:62
bool empty() const noexcept
True if the list is empty (ie, size() is zero).
Definition UPtrListI.H:99
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
label count() const noexcept
The number of non-nullptr entries in the list.
Definition UPtrList.H:890
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
const lduSchedule & patchSchedule() const noexcept
Order in which the patches should be initialised/evaluated corresponding to the schedule.
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
bool checkDefinition(const bool report=false) const
Check boundary definition.
label nNonProcessorFaces() const
The number of boundary faces before the first processor patch.
label nNonProcessor() const
The number of patches before the first processor patch.
wordList physicalTypes() const
Return a list of physical types.
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
labelList patchStarts() const
Return a list of patch start face indices.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write using stream options, but always UNCOMPRESSED.
const List< labelPairList > & neighbourEdges() const
Per patch the edges on the neighbouring patch.
void reorder(const labelUList &oldToNew, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
wordList types() const
Return a list of patch types.
void clearGeom()
Clear geometry at this level and at patches.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
virtual bool writeData(Ostream &os) const
The writeData member function required by regIOobject.
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
const polyPatch & operator[](const word &patchName) const
Return const reference to polyPatch by name.
const faceList::subList faces() const
Return mesh faces for the entire boundary.
void setGroup(const word &groupName, const labelUList &patchIDs)
Set/add group with patches.
labelRange range() const noexcept
The face range for all boundary faces.
void movePoints(const pointField &p)
Correct polyBoundaryMesh after moving points.
List< labelRange > patchRanges() const
Return a list of patch ranges.
void matchGroups(const labelUList &patchIDs, wordList &groups, labelHashSet &nonGroupPatches) const
Match the patches to groups.
void clearAddressing()
Clear addressing at this level and at patches.
label nFaces() const noexcept
The number of boundary faces in the underlying mesh.
friend class polyMesh
Declare friendship with polyMesh.
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
void writeEntry(Ostream &os) const
Write as a plain list of entries.
const labelList & patchID() const
Per boundary face label the patch index.
const polyMesh & mesh() const noexcept
Return the mesh reference.
labelList patchSizes() const
Return a list of patch sizes.
void clear()
Clear the patch list and all demand-driven data.
wordList groupNames() const
A list of the group names (if any).
wordList names() const
Return a list of patch names.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
void updateMesh()
Correct polyBoundaryMesh after topology update.
const polyPatch * cfindPatch(const word &patchName) const
Find patch by name and return const pointer.
const labelList::subList faceOwner() const
Return face owner for the entire boundary.
label nProcessorPatches() const
The number of processorPolyPatch patches.
polyBoundaryMesh(const polyBoundaryMesh &)=delete
No copy construct.
labelList indices(const wordRe &matcher, const bool useGroups=true) const
The (sorted) patch indices for all matches, optionally matching patch groups.
label findIndex(const wordRe &key) const
Return patch index for the first match, return -1 if not found.
labelPair whichPatchFace(const label meshFacei) const
Lookup mesh face index and return (patchi, patchFacei) tuple or (-1, meshFacei) for internal faces.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
label nInternalFaces() const noexcept
Number of internal faces.
label nFaces() const noexcept
Number of mesh faces.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
@ BEGIN_LIST
Begin list [isseparator].
Definition token.H:174
@ END_LIST
End list [isseparator].
Definition token.H:175
@ SPACE
Space [isspace].
Definition token.H:144
@ NL
Newline [isspace].
Definition token.H:143
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
bool isPattern() const noexcept
The wordRe is a pattern, not a literal string.
Definition wordReI.H:104
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & name
auto & names
label patchId(-1)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
labelList findMatching(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Extract list indices for all items with 'name()' that matches.
label firstMatching(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Find first list item with 'name()' that matches, -1 on failure.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
List< word > wordList
List of word.
Definition fileName.H:60
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< labelPair > labelPairList
List of labelPair.
Definition labelPair.H:33
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition lduSchedule.H:46
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
List< bool > boolList
A List of bools.
Definition List.H:60
const direction noexcept
Definition scalarImpl.H:265
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
dictionary dict
volScalarField & e
#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
Dispatch tag: Construct 'one-sided' from the non-master local sizes using gather but no broadcast.
Extract name (as a word) from an object, typically using its name() method.
Definition word.H:341
Extract type (as a word) from an object, typically using its type() method.
Definition word.H:362
Functor wrapper of allow/deny lists of wordRe for filtering.
Definition wordRes.H:275