Loading...
Searching...
No Matches
faBoundaryMesh.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) 2016-2017 Wikki Ltd
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 "faBoundaryMesh.H"
30#include "faMesh.H"
31#include "globalIndex.H"
32#include "processorFaPatch.H"
33#include "wordRes.H"
34#include "PtrListOps.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
46bool Foam::faBoundaryMesh::hasGroupIDs() const
47{
48 if (groupIDsPtr_)
49 {
50 // Use existing cache
51 return !groupIDsPtr_->empty();
52 }
53
54 const faPatchList& patches = *this;
55
56 for (const faPatch& p : patches)
57 {
58 if (!p.inGroups().empty())
59 {
60 return true;
61 }
62 }
63
64 return false;
65}
66
67
68void Foam::faBoundaryMesh::calcGroupIDs() const
69{
70 if (groupIDsPtr_)
71 {
72 return; // Or FatalError
73 }
74
75 groupIDsPtr_.emplace(16);
76 auto& groupLookup = *groupIDsPtr_;
77
78 const faPatchList& patches = *this;
79
80 forAll(patches, patchi)
81 {
82 for (const word& groupName : patches[patchi].inGroups())
83 {
84 groupLookup(groupName).push_back(patchi);
85 }
86 }
87
88 // Remove groups that clash with patch names
89 forAll(patches, patchi)
90 {
91 if (groupLookup.empty())
92 {
93 break; // Early termination
94 }
95 else if (groupLookup.erase(patches[patchi].name()))
96 {
98 << "Removed group '" << patches[patchi].name()
99 << "' which clashes with patch " << patchi
100 << " of the same name."
101 << endl;
102 }
103 }
104}
105
106
107void Foam::faBoundaryMesh::populate(PtrList<entry>&& entries)
108{
109 clearLocalAddressing();
110
111 faPatchList& patches = *this;
112
113 patches.resize_null(entries.size());
114
115 // Transcribe.
116 // Does not handle nullptr at all (what could possibly be done?)
117 forAll(patches, patchi)
118 {
120 (
121 patchi,
123 (
124 entries[patchi].keyword(),
125 entries[patchi].dict(),
126 patchi,
127 *this
128 )
129 );
130 }
131
132 entries.clear();
133}
134
135
136bool Foam::faBoundaryMesh::readIOcontents(const bool allowOptionalRead)
137{
138 bool updated = false;
139 PtrList<entry> entries;
140
141 if
142 (
143 this->isReadRequired()
144 || (allowOptionalRead && this->isReadOptional() && this->headerOk())
145 )
146 {
147 // Warn for MUST_READ_IF_MODIFIED
148 warnNoRereading<faBoundaryMesh>();
149
150 // Read entries
151 Istream& is = readStream(typeName);
152
153 is >> entries;
154
155 is.check(FUNCTION_NAME);
156 close();
157 updated = true;
158 }
159 // Future: support master-only and broadcast?
160
161 if (updated)
162 {
163 populate(std::move(entries));
164 }
166 return updated;
167}
168
169
170// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
171
173(
174 const IOobject& io,
175 const faMesh& mesh
176)
177:
178 faPatchList(),
180 mesh_(mesh)
181{
182 readIOcontents(false); // allowOptionalRead = false
183}
184
185
187(
188 const IOobject& io,
189 const faMesh& pm,
191)
193 faPatchList(),
195 mesh_(pm)
196{}
197
198
200(
201 const IOobject& io,
202 const faMesh& pm,
203 const label size
204)
208 mesh_(pm)
209{}
210
211
213(
214 const IOobject& io,
215 const faMesh& fam,
216 const faPatchList& list
217)
218:
219 faPatchList(),
221 mesh_(fam)
222{
223 if (!readIOcontents(true)) // allowOptionalRead = true
224 {
225 // Nothing read. Use supplied patches
226 faPatchList& patches = *this;
227 patches.resize(list.size());
228
229 forAll(patches, patchi)
231 patches.set(patchi, list[patchi].clone(*this));
232 }
233 }
234}
235
236
238(
239 const IOobject& io,
240 const faMesh& fam,
241 PtrList<entry>&& entries
242)
243:
244 faPatchList(),
246 mesh_(fam)
247{
248 if (!readIOcontents(true)) // allowOptionalRead = true
249 {
250 populate(std::move(entries));
252 entries.clear();
253}
254
255
256// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
257
259{
260 clearLocalAddressing();
262}
263
264
265void Foam::faBoundaryMesh::clearLocalAddressing()
267 groupIDsPtr_.reset(nullptr);
268}
269
270
271// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
272
274{
275 // processor initGeometry send/recv the following:
276 // - edgeCentres() : faMesh::edgeCentres()
277 // - edgeLengths() : faMesh::Le()
278 // - edgeFaceCentres() : faMesh::areaCentres()
279 //
280 // faMesh::Le() has its own point-to-point communication (OK) but
281 // triggers either/or edgeAreaNormals(), pointAreaNormals()
282 // with their own communication that can block.
283
284 // This uses parallel comms and hence will not be trigggered
285 // on processors that do not have a processorFaPatch so instead
286 // force construction.
287
288 (void)mesh_.edgeAreaNormals();
289 (void)mesh_.pointAreaNormals();
290
291 (void)mesh_.areaCentres();
292 (void)mesh_.faceAreaNormals();
293
294
296
297 if
298 (
299 pBufs.commsType() == Pstream::commsTypes::buffered
300 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
301 )
302 {
303 forAll(*this, patchi)
304 {
305 operator[](patchi).initGeometry(pBufs);
306 }
307
308 pBufs.finishedSends();
309
310 forAll(*this, patchi)
311 {
312 operator[](patchi).calcGeometry(pBufs);
313 }
314 }
315 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
316 {
317 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
318
319 // Dummy.
320 pBufs.finishedSends();
321
322 for (const auto& patchEval : patchSchedule)
323 {
324 const label patchi = patchEval.patch;
325
326 if (patchEval.init)
327 {
328 operator[](patchi).initGeometry(pBufs);
329 }
330 else
331 {
332 operator[](patchi).calcGeometry(pBufs);
334 }
335 }
336}
337
338
339Foam::UPtrList<const Foam::labelUList>
341{
342 const faPatchList& patches = *this;
343
345
346 forAll(list, patchi)
347 {
348 list.set(patchi, &patches[patchi].edgeLabels());
350
351 return list;
352}
353
354
357{
358 const faPatchList& patches = *this;
359
360 UPtrList<const labelUList> list(patches.size());
361
362 forAll(list, patchi)
363 {
364 list.set(patchi, &patches[patchi].edgeFaces());
365 }
366
367 return list;
368}
369
370
372{
373 const faPatchList& patches = *this;
374
375 lduInterfacePtrsList list(patches.size());
376
377 forAll(list, patchi)
378 {
379 const lduInterface* lduPtr = isA<lduInterface>(patches[patchi]);
380
381 if (lduPtr)
382 {
383 list.set(patchi, lduPtr);
385 }
386
387 return list;
388}
389
390
391Foam::label Foam::faBoundaryMesh::nNonProcessor() const
392{
393 const faPatchList& patches = *this;
394
395 label count = 0;
396
397 for (const faPatch& p : patches)
398 {
400 {
401 break;
402 }
403
405 }
406
407 return count;
408}
409
410
412{
413 const faPatchList& patches = *this;
414
415 label count = 0;
416
417 for (const faPatch& p : patches)
418 {
420 {
421 ++count;
423 }
424
425 return count;
426}
427
428
430{
431 const faPatchList& patches = *this;
432
433 label count = 0;
434
435 for (const faPatch& p : patches)
436 {
438 {
439 break;
440 }
441
442 count += p.nEdges();
444
445 return count;
446}
447
448
451{
452 if (!groupIDsPtr_)
453 {
454 calcGroupIDs();
455 }
456
457 return *groupIDsPtr_;
458}
459
460
462(
463 const word& groupName,
464 const labelUList& patchIDs
465)
466{
467 groupIDsPtr_.reset(nullptr);
468
469 faPatchList& patches = *this;
470
471 boolList pending(patches.size(), true);
472
473 // Add to specified patches
474 for (const label patchi : patchIDs)
475 {
476 if (pending.test(patchi))
477 {
478 pending.unset(patchi);
479 patches[patchi].addGroup(groupName);
480 }
481 }
482
483 // Remove from other patches
484 forAll(patches, patchi)
485 {
486 if (pending.test(patchi))
488 patches[patchi].removeGroup(groupName);
489 }
490 }
491}
492
497}
498
501{
502 return PtrListOps::get<word>(*this, typeOp<faPatch>());
503}
504
505
507{
508 // Manually: faPatch does not have independent start() information
509
510 const faPatchList& patches = *this;
511
512 labelList list(patches.size());
513
514 label beg = mesh_.nInternalEdges();
515 forAll(patches, patchi)
516 {
517 const label len = patches[patchi].nEdges();
518 list[patchi] = beg;
519 beg += len;
520 }
521 return list;
522}
523
524
526{
527 return
530 *this,
531 [](const faPatch& p) { return p.nEdges(); } // avoid virtual
532 );
533}
534
535
537{
538 const faPatchList& patches = *this;
539
540 List<labelRange> list(patches.size());
541
542 label beg = mesh_.nInternalEdges();
543 forAll(patches, patchi)
544 {
545 const label len = patches[patchi].nEdges();
546 list[patchi].reset(beg, len);
547 beg += len;
548 }
549 return list;
550}
551
554{
555 return this->groupPatchIDs().sortedToc();
556}
557
559Foam::label Foam::faBoundaryMesh::start() const
560{
561 return mesh_.nInternalEdges();
562}
563
565Foam::label Foam::faBoundaryMesh::nEdges() const
566{
567 return mesh_.nBoundaryEdges();
568}
569
572{
573 return labelRange(mesh_.nInternalEdges(), mesh_.nBoundaryEdges());
574}
575
576
578(
579 const wordRe& matcher,
580 const bool useGroups
581) const
582{
583 if (matcher.empty())
584 {
585 return labelList();
586 }
587
588 // Only check groups if requested and they exist
589 const bool checkGroups = (useGroups && this->hasGroupIDs());
590
591 labelHashSet ids;
592
593 if (matcher.isPattern())
594 {
595 if (checkGroups)
596 {
597 ids.reserve(this->size());
598
599 const auto& groupLookup = groupPatchIDs();
600 forAllConstIters(groupLookup, iter)
601 {
602 if (matcher(iter.key()))
603 {
604 // Add ids associated with the group
605 ids.insert(iter.val());
606 }
607 }
608 }
609
610 if (ids.empty())
611 {
612 return PtrListOps::findMatching(*this, matcher);
613 }
614 else
615 {
616 ids.insert(PtrListOps::findMatching(*this, matcher));
617 }
618 }
619 else
620 {
621 // Literal string.
622 // Special version of above for reduced memory footprint
623
624 const label patchId = PtrListOps::firstMatching(*this, matcher);
625
626 if (patchId >= 0)
627 {
628 return labelList(one{}, patchId);
629 }
630 else if (checkGroups)
631 {
632 const auto iter = groupPatchIDs().cfind(matcher);
633
634 if (iter.good())
635 {
636 // Add ids associated with the group
637 ids.insert(iter.val());
638 }
640 }
641
642 return ids.sortedToc();
643}
644
645
647(
648 const wordRes& matcher,
649 const bool useGroups
650) const
651{
652 if (matcher.empty())
653 {
654 return labelList();
655 }
656 else if (matcher.size() == 1)
657 {
658 return this->indices(matcher.front(), useGroups);
659 }
660
661 labelHashSet ids;
662
663 // Only check groups if requested and they exist
664 if (useGroups && this->hasGroupIDs())
665 {
666 ids.reserve(this->size());
667
668 const auto& groupLookup = groupPatchIDs();
669 forAllConstIters(groupLookup, iter)
670 {
671 if (matcher(iter.key()))
672 {
673 // Add ids associated with the group
674 ids.insert(iter.val());
675 }
676 }
677 }
678
679 if (ids.empty())
680 {
681 return PtrListOps::findMatching(*this, matcher);
682 }
683 else
684 {
685 ids.insert(PtrListOps::findMatching(*this, matcher));
686 }
687
688 return ids.sortedToc();
689}
690
691
693(
694 const wordRes& allow,
695 const wordRes& deny,
696 const bool useGroups
697) const
698{
699 if (allow.empty() && deny.empty())
700 {
701 // Fast-path: select all
702 return identity(this->size());
703 }
704
705 const wordRes::filter matcher(allow, deny);
706
707 labelHashSet ids;
708
709 // Only check groups if requested and they exist
710 if (useGroups && this->hasGroupIDs())
711 {
712 ids.reserve(this->size());
713
714 const auto& groupLookup = groupPatchIDs();
715 forAllConstIters(groupLookup, iter)
716 {
717 if (matcher(iter.key()))
718 {
719 // Add ids associated with the group
720 ids.insert(iter.val());
721 }
722 }
723 }
724
725 if (ids.empty())
726 {
727 return PtrListOps::findMatching(*this, matcher);
728 }
729 else
730 {
731 ids.insert(PtrListOps::findMatching(*this, matcher));
732 }
733
734 return ids.sortedToc();
735}
736
737
738Foam::label Foam::faBoundaryMesh::findIndex(const wordRe& key) const
739{
740 if (key.empty())
742 return -1;
743 }
744 return PtrListOps::firstMatching(*this, key);
745}
746
747
749(
750 const word& patchName,
751 bool allowNotFound
752) const
753{
754 if (patchName.empty())
755 {
756 return -1;
757 }
758
759 const label patchId = PtrListOps::firstMatching(*this, patchName);
760
761 if (patchId >= 0)
762 {
763 return patchId;
764 }
765
766 if (!allowNotFound)
767 {
769 << "Patch '" << patchName << "' not found. "
770 << "Available patch names: " << names() << endl
771 << exit(FatalError);
772 }
773
774 // Patch not found
775 if (debug)
776 {
777 Pout<< "label faBoundaryMesh::findPatchID(const word&) const"
778 << "Patch named " << patchName << " not found. "
779 << "Available patch names: " << names() << endl;
781
782 // Not found, return -1
783 return -1;
784}
785
786
788(
789 const word& patchName
790) const
791{
792 const faPatchList& patches = *this;
793
794 if (!patchName.empty())
795 {
796 // Note: get() handles out-of-range access properly
798 }
799
800 return nullptr;
801}
802
803
804Foam::label Foam::faBoundaryMesh::whichPatch(const label edgeIndex) const
805{
806 if (edgeIndex < mesh().nInternalEdges())
807 {
808 // Internal edge
809 return -1;
810 }
811 else if (edgeIndex >= mesh().nEdges())
812 {
813 // Bounds error: abort
815 << "Edge " << edgeIndex
816 << " out of bounds. Number of geometric edges " << mesh().nEdges()
817 << abort(FatalError);
818
819 return -1;
820 }
821
822 // Find patch that the edgeIndex belongs to.
823
824 forAll(*this, patchi)
825 {
826 label start = mesh_.patchStarts()[patchi];
827 label size = operator[](patchi).faPatch::size();
828
829 if (edgeIndex >= start && edgeIndex < start + size)
830 {
831 return patchi;
832 }
833 }
834
835 // If not in any of above, it's trouble!
837 << "Error in patch search algorithm"
838 << abort(FatalError);
839
840 return -1;
841}
842
843
844bool Foam::faBoundaryMesh::checkParallelSync(const bool report) const
845{
846 if (!Pstream::parRun())
847 {
848 return false;
849 }
850
851 const faBoundaryMesh& bm = *this;
852
853 bool hasError = false;
854
855 // Collect non-proc patches and check proc patches are last.
856 wordList localNames(bm.size());
857 wordList localTypes(bm.size());
858
859 label nonProci = 0;
860
861 forAll(bm, patchi)
862 {
863 if (!isA<processorFaPatch>(bm[patchi]))
864 {
865 if (nonProci != patchi)
866 {
867 // A processor patch in between normal patches!
868 hasError = true;
869
870 if (debug || report)
871 {
872 Pout<< " ***Problem with boundary patch " << patchi
873 << " name:" << bm[patchi].name()
874 << " type:" << bm[patchi].type()
875 << " - seems to be preceeded by processor patches."
876 << " This is usually a problem." << endl;
877 }
878 }
879 else
880 {
881 localNames[nonProci] = bm[patchi].name();
882 localTypes[nonProci] = bm[patchi].type();
883 ++nonProci;
884 }
885 }
886 }
887 localNames.resize(nonProci);
888 localTypes.resize(nonProci);
889
890 // Check and report error(s) on master
891 // - don't need indexing on master itself
892
893 const globalIndex procAddr(globalIndex::gatherNonLocal{}, nonProci);
894
895 const wordList allNames(procAddr.gather(localNames));
896 const wordList allTypes(procAddr.gather(localTypes));
897
898 // Automatically restricted to master
899 for (const int proci : procAddr.subProcs())
900 {
901 const auto procNames(allNames.slice(procAddr.range(proci)));
902 const auto procTypes(allTypes.slice(procAddr.range(proci)));
903
904 if (procNames != localNames || procTypes != localTypes)
905 {
906 hasError = true;
907
908 if (debug || report)
909 {
910 Info<< " ***Inconsistent patches across processors, "
911 "processor0 has patch names:" << localNames
912 << " patch types:" << localTypes
913 << " processor" << proci
914 << " has patch names:" << procNames
915 << " patch types:" << procTypes
916 << endl;
917 }
918 }
920
921 // Reduce (not broadcast) to respect local out-of-order errors (first loop)
922 return returnReduceOr(hasError);
923}
924
925
926bool Foam::faBoundaryMesh::checkDefinition(const bool report) const
927{
928 label nextPatchStart = mesh().nInternalEdges();
929 const faBoundaryMesh& bm = *this;
930
931 bool hasError = false;
932
933 forAll(bm, patchi)
934 {
935 if (bm[patchi].start() != nextPatchStart && !hasError)
936 {
937 hasError = true;
938
940 << " ****Problem with boundary patch " << patchi
941 << " named " << bm[patchi].name()
942 << " of type " << bm[patchi].type()
943 << ". The patch should start on face no " << nextPatchStart
944 << " and the patch specifies " << bm[patchi].start()
945 << "." << endl
946 << "Possibly consecutive patches have this same problem."
947 << " Suppressing future warnings." << endl;
948 }
949
950 // Warn about duplicate boundary patches?
951
952 nextPatchStart += bm[patchi].faPatch::size();
953 }
954
955 if (hasError)
956 {
958 << "This mesh is not valid: boundary definition is in error."
959 << endl;
960 }
961 else
962 {
963 if (debug || report)
964 {
965 Info << "Boundary definition OK." << endl;
967 }
968
969 return hasError;
970}
971
972
974{
975 // See comments in calcGeometry()
976
977 (void)mesh_.edgeAreaNormals();
978 (void)mesh_.pointAreaNormals();
979
980 (void)mesh_.areaCentres();
981 (void)mesh_.faceAreaNormals();
982
983
985
986 if
987 (
988 pBufs.commsType() == Pstream::commsTypes::buffered
989 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
990 )
991 {
992 forAll(*this, patchi)
993 {
994 operator[](patchi).initMovePoints(pBufs, p);
995 }
996
997 pBufs.finishedSends();
998
999 forAll(*this, patchi)
1000 {
1001 operator[](patchi).movePoints(pBufs, p);
1002 }
1003 }
1004 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1005 {
1006 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1007
1008 // Dummy.
1009 pBufs.finishedSends();
1010
1011 for (const auto& schedEval : patchSchedule)
1012 {
1013 const label patchi = schedEval.patch;
1014
1015 if (schedEval.init)
1016 {
1017 operator[](patchi).initMovePoints(pBufs, p);
1018 }
1019 else
1020 {
1021 operator[](patchi).movePoints(pBufs, p);
1022 }
1023 }
1024 }
1025}
1026
1027
1029{
1031
1032 if
1033 (
1034 pBufs.commsType() == Pstream::commsTypes::buffered
1035 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
1036 )
1037 {
1038 forAll(*this, patchi)
1039 {
1040 operator[](patchi).initUpdateMesh(pBufs);
1041 }
1042
1043 pBufs.finishedSends();
1044
1045 forAll(*this, patchi)
1046 {
1047 operator[](patchi).updateMesh(pBufs);
1048 }
1049 }
1050 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1051 {
1052 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1053
1054 // Dummy.
1055 pBufs.finishedSends();
1056
1057 for (const auto& schedEval : patchSchedule)
1058 {
1059 const label patchi = schedEval.patch;
1060
1061 if (schedEval.init)
1062 {
1063 operator[](patchi).initUpdateMesh(pBufs);
1064 }
1065 else
1066 {
1067 operator[](patchi).updateMesh(pBufs);
1068 }
1069 }
1070 }
1071}
1072
1073
1075{
1076 const faPatchList& entries = *this;
1077
1078 os << entries.size();
1079
1080 if (entries.empty())
1081 {
1082 // 0-sized : can write with less vertical space
1084 }
1085 else
1086 {
1087 os << nl << token::BEGIN_LIST << incrIndent << nl;
1088
1089 for (const auto& pp : entries)
1090 {
1091 os.beginBlock(pp.name());
1092 os << pp;
1093 os.endBlock();
1096 }
1097 os.check(FUNCTION_NAME);
1098}
1099
1100
1102(
1103 const word& keyword,
1104 Ostream& os
1105) const
1106{
1107 const faPatchList& entries = *this;
1108
1109 if (!keyword.empty())
1110 {
1111 os.write(keyword);
1112 os << (entries.empty() ? token::SPACE : token::NL);
1113 }
1115 writeEntry(os);
1116
1117 if (!keyword.empty()) os.endEntry();
1118}
1119
1120
1122{
1123 writeEntry(os);
1124 return os.good();
1125}
1126
1127
1129(
1130 IOstreamOption streamOpt,
1131 const bool writeOnProc
1132) const
1133{
1135 return regIOobject::writeObject(streamOpt, writeOnProc);
1136}
1137
1138
1139// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1140
1141Foam::Ostream& Foam::operator<<(Ostream& os, const faBoundaryMesh& bm)
1142{
1143 bm.writeData(os);
1144 return os;
1145}
1146
1147
1148// ************************************************************************* //
Functions to operate on Pointer Lists.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
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
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
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 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
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
Definition OSstream.H:134
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
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
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 commsTypes defaultCommsType
Default commsType.
Definition UPstream.H:1045
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
const T & operator[](const label i) const
Return const reference to the element at given position. FatalError for bounds problem or nullptr....
Definition UPtrListI.H:289
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
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
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
UPtrList< const labelUList > edgeFaces() const
Return a list of edgeFaces for each patch.
bool checkDefinition(const bool report=false) const
Check boundary definition.
label nNonProcessor() const
The number of patches before the first processor patch.
void calcGeometry()
Calculate the geometry for the patches.
const faPatch * cfindPatch(const word &patchName) const
Find patch by name and return const pointer.
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
labelList patchStarts() const
Return a list of patch start indices.
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch with only those pointing to interfaces being set.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write using stream options, but always UNCOMPRESSED.
const faMesh & mesh() const noexcept
Return the mesh reference.
labelRange range() const
The edge range for all boundary edges.
wordList types() const
Return a list of patch types.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
virtual bool writeData(Ostream &os) const
The writeData member function required by regIOobject.
void setGroup(const word &groupName, const labelUList &patchIDs)
Set/add group with patches.
List< labelRange > patchRanges() const
Return a list of patch ranges.
faBoundaryMesh(const faBoundaryMesh &)=delete
No copy construct.
label whichPatch(const label edgeIndex) const
Return patch index for a given edge label.
void movePoints(const pointField &)
Correct faBoundaryMesh after moving points.
void writeEntry(Ostream &os) const
Write as a plain list of entries.
labelList patchSizes() const
Return a list of patch sizes (number of edges in each patch).
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 faBoundaryMesh after topology update.
label nProcessorPatches() const
The number of processorFaPatch patches.
UPtrList< const labelUList > edgeLabels() const
Return a list of edgeLabels for each patch.
label nNonProcessorEdges() const
The number of boundary edges before the first processor patch.
labelList indices(const wordRe &matcher, const bool useGroups=true) const
The (sorted) patch indices for all matches, optionally matching patch groups.
label nEdges() const
The number of boundary edges for the underlying mesh.
label findIndex(const wordRe &key) const
Return patch index for the first match, return -1 if not found.
label start() const
The start label of the edges in the faMesh edges list.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition faPatch.H:76
static autoPtr< faPatch > New(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm)
Return pointer to a new patch created on freestore from dictionary.
Definition faPatchNew.C:28
Type type(bool followLink=true, bool checkGzip=false) const
Return the directory entry type: UNDEFINED, FILE, DIRECTORY (or SYMLINK).
Definition fileName.C:353
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
An abstract base class for implicitly-coupled interfaces e.g. processor and cyclic patches.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
label nInternalEdges() const
Internal edges using 0,1 or 2 boundary points.
label nEdges() const
Number of mesh edges.
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 & names
label patchId(-1)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
#define InfoInFunction
Report an information message using Foam::Info.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
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< 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.
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
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
messageStream Info
Information stream (stdout output on master, null elsewhere).
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
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
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.
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition faPatch.H:64
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
Calculate the matrix for the second temporal derivative.
dictionary dict
#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