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-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 "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 "stringListOps.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_.reset(new HashTable<labelList>(16));
81 auto& groupLookup = *groupIDsPtr_;
82
83 const polyPatchList& patches = *this;
84
85 forAll(patches, patchi)
86 {
87 const wordList& groups = patches[patchi].inGroups();
88
89 for (const word& groupName : groups)
90 {
91 groupLookup(groupName).append(patchi);
92 }
93 }
94
95 // Remove groups that clash with patch names
96 forAll(patches, patchi)
97 {
98 if (groupLookup.erase(patches[patchi].name()))
99 {
101 << "Removed group '" << patches[patchi].name()
102 << "' which clashes with patch " << patchi
103 << " of the same name."
104 << endl;
105 }
106 }
107}
108
109
110bool Foam::polyBoundaryMesh::readContents(const bool allowReadIfPresent)
111{
112 if
113 (
114 this->readOpt() == IOobject::MUST_READ
115 || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
116 ||
117 (
118 allowReadIfPresent
119 && (this->readOpt() == IOobject::READ_IF_PRESENT && this->headerOk())
120 )
121 )
122 {
123 // Warn for MUST_READ_IF_MODIFIED
124 warnNoRereading<polyBoundaryMesh>();
125
126 polyPatchList& patches = *this;
127
128 // Read polyPatchList
129 Istream& is = readStream(typeName);
130
131 // Read patches as entries
132 PtrList<entry> patchEntries(is);
133 patches.resize(patchEntries.size());
134
135 // Transcribe
136 forAll(patches, patchi)
137 {
139 (
140 patchi,
142 (
143 patchEntries[patchi].keyword(),
144 patchEntries[patchi].dict(),
145 patchi,
146 *this
147 )
148 );
149 }
150
151 is.check(FUNCTION_NAME);
152 close();
153 return true;
154 }
155
156 return false;
157}
158
159
160// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161
163(
164 const IOobject& io,
165 const polyMesh& mesh
166)
167:
170 mesh_(mesh)
171{
172 readContents(false); // READ_IF_PRESENT allowed: False
173}
174
175
177(
178 const IOobject& io,
179 const polyMesh& pm,
180 const label size
181)
182:
183 polyPatchList(size),
185 mesh_(pm)
186{}
187
188
190(
191 const IOobject& io,
192 const polyMesh& pm,
193 const polyPatchList& ppl
194)
195:
198 mesh_(pm)
199{
200 if (!readContents(true)) // READ_IF_PRESENT allowed: True
201 {
202 polyPatchList& patches = *this;
203 patches.resize(ppl.size());
204 forAll(patches, patchi)
205 {
206 patches.set(patchi, ppl[patchi].clone(*this));
207 }
208 }
209}
210
211
212// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
213
215{
216 polyPatchList& patches = *this;
217
218 for (polyPatch& p : patches)
219 {
220 p.clearGeom();
221 }
222}
223
224
226{
227 neighbourEdgesPtr_.clear();
228 patchIDPtr_.clear();
229 groupIDsPtr_.clear();
230
231 polyPatchList& patches = *this;
232
233 for (polyPatch& p : patches)
234 {
235 p.clearAddressing();
236 }
237}
238
239
240// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
241
242void Foam::polyBoundaryMesh::calcGeometry()
243{
245
246 if
247 (
248 pBufs.commsType() == Pstream::commsTypes::blocking
249 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
250 )
251 {
252 forAll(*this, patchi)
253 {
254 operator[](patchi).initGeometry(pBufs);
255 }
256
257 pBufs.finishedSends();
258
259 forAll(*this, patchi)
260 {
261 operator[](patchi).calcGeometry(pBufs);
262 }
263 }
264 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
265 {
266 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
267
268 // Dummy.
269 pBufs.finishedSends();
270
271 for (const auto& schedEval : patchSchedule)
272 {
273 const label patchi = schedEval.patch;
274
275 if (schedEval.init)
276 {
277 operator[](patchi).initGeometry(pBufs);
278 }
279 else
280 {
281 operator[](patchi).calcGeometry(pBufs);
282 }
283 }
284 }
285}
286
287
290{
291 const polyPatchList& patches = *this;
292
294
295 forAll(patches, patchi)
296 {
297 list.set(patchi, &patches[patchi].faceCells());
298 }
299
300 return list;
301}
302
303
306{
307 if (Pstream::parRun())
308 {
310 << "Neighbour edge addressing not correct across parallel"
311 << " boundaries." << endl;
312 }
313
314 if (!neighbourEdgesPtr_)
315 {
316 neighbourEdgesPtr_.reset(new List<labelPairList>(size()));
317 auto& neighbourEdges = *neighbourEdgesPtr_;
318
319 // Initialize.
320 label nEdgePairs = 0;
321 forAll(*this, patchi)
322 {
323 const polyPatch& pp = operator[](patchi);
324
325 neighbourEdges[patchi].setSize(pp.nEdges() - pp.nInternalEdges());
326
327 for (labelPair& edgeInfo : neighbourEdges[patchi])
328 {
329 edgeInfo[0] = -1;
330 edgeInfo[1] = -1;
331 }
332
333 nEdgePairs += pp.nEdges() - pp.nInternalEdges();
334 }
335
336 // From mesh edge (expressed as a point pair so as not to construct
337 // point addressing) to patch + relative edge index.
338 EdgeMap<labelPair> pointsToEdge(nEdgePairs);
339
340 forAll(*this, patchi)
341 {
342 const polyPatch& pp = operator[](patchi);
343
344 const edgeList& edges = pp.edges();
345
346 for
347 (
348 label edgei = pp.nInternalEdges();
349 edgei < edges.size();
350 edgei++
351 )
352 {
353 // Edge in patch local points
354 const edge& e = edges[edgei];
355
356 // Edge in mesh points.
357 edge meshEdge(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
358
359 auto fnd = pointsToEdge.find(meshEdge);
360
361 if (!fnd.found())
362 {
363 // First occurrence of mesh edge. Store patch and my
364 // local index.
365 pointsToEdge.insert
366 (
367 meshEdge,
369 (
370 patchi,
371 edgei - pp.nInternalEdges()
372 )
373 );
374 }
375 else
376 {
377 // Second occurrence. Store.
378 const labelPair& edgeInfo = fnd.val();
379
380 neighbourEdges[patchi][edgei - pp.nInternalEdges()] =
381 edgeInfo;
382
383 neighbourEdges[edgeInfo[0]][edgeInfo[1]]
384 = labelPair(patchi, edgei - pp.nInternalEdges());
385
386 // Found all two occurrences of this edge so remove from
387 // hash to save space. Note that this will give lots of
388 // problems if the polyBoundaryMesh is multiply connected.
389 pointsToEdge.erase(meshEdge);
390 }
391 }
392 }
393
394 if (pointsToEdge.size())
395 {
397 << "Not all boundary edges of patches match up." << nl
398 << "Is the outside of your mesh multiply connected?"
399 << abort(FatalError);
400 }
401
402 forAll(*this, patchi)
403 {
404 const polyPatch& pp = operator[](patchi);
405
406 const labelPairList& nbrEdges = neighbourEdges[patchi];
407
408 forAll(nbrEdges, i)
409 {
410 const labelPair& edgeInfo = nbrEdges[i];
411
412 if (edgeInfo[0] == -1 || edgeInfo[1] == -1)
413 {
414 const label edgei = pp.nInternalEdges() + i;
415 const edge& e = pp.edges()[edgei];
416
418 << "Not all boundary edges of patches match up." << nl
419 << "Edge " << edgei << " on patch " << pp.name()
420 << " end points " << pp.localPoints()[e[0]] << ' '
421 << pp.localPoints()[e[1]] << " is not matched to an"
422 << " edge on any other patch." << nl
423 << "Is the outside of your mesh multiply connected?"
424 << abort(FatalError);
425 }
426 }
427 }
428 }
429
430 return *neighbourEdgesPtr_;
431}
432
433
435{
436 if (!patchIDPtr_)
437 {
438 patchIDPtr_.reset(new labelList(mesh_.nBoundaryFaces()));
439 labelList& list = *patchIDPtr_;
440
441 const polyPatchList& patches = *this;
442
443 forAll(patches, patchi)
444 {
446 (
447 list,
448 patches[patchi].size(),
449 (patches[patchi].start() - mesh_.nInternalFaces())
450 ) = patchi;
451 }
452 }
453
454 return *patchIDPtr_;
455}
456
457
460{
461 if (!groupIDsPtr_)
462 {
463 calcGroupIDs();
464 }
465
466 return *groupIDsPtr_;
467}
468
469
471(
472 const word& groupName,
473 const labelUList& patchIDs
474)
475{
476 groupIDsPtr_.clear();
477
478 polyPatchList& patches = *this;
479
480 boolList donePatch(patches.size(), false);
481
482 // Add to specified patches
483 for (const label patchi : patchIDs)
484 {
485 patches[patchi].inGroups().appendUniq(groupName);
486 donePatch[patchi] = true;
487 }
488
489 // Remove from other patches
490 forAll(patches, patchi)
491 {
492 if (!donePatch[patchi])
493 {
494 wordList& groups = patches[patchi].inGroups();
495
496 if (groups.found(groupName))
497 {
498 label newi = 0;
499 forAll(groups, i)
500 {
501 if (groups[i] != groupName)
502 {
503 groups[newi++] = groups[i];
504 }
505 }
506 groups.resize(newi);
507 }
508 }
509 }
510}
511
512
514{
515 const polyPatchList& patches = *this;
516
517 label nonProc = 0;
518
519 for (const polyPatch& p : patches)
520 {
521 if (isA<processorPolyPatch>(p))
522 {
523 break;
524 }
525
526 ++nonProc;
527 }
528
529 return nonProc;
530}
531
532
534{
535 return PtrListOps::get<word>(*this, nameOp<polyPatch>());
536}
537
538
540{
541 return PtrListOps::get<word>(*this, typeOp<polyPatch>());
542}
543
544
546{
547 return
548 PtrListOps::get<word>
549 (
550 *this,
551 [](const polyPatch& p) { return p.physicalType(); }
552 );
553}
554
555
557{
558 return
559 PtrListOps::get<label>
560 (
561 *this,
562 [](const polyPatch& p) { return p.start(); }
563 );
564}
565
566
568{
569 return
570 PtrListOps::get<label>
571 (
572 *this,
573 [](const polyPatch& p) { return p.size(); }
574 );
575}
576
577
579{
580 return
581 PtrListOps::get<labelRange>
582 (
583 *this,
584 [](const polyPatch& p) { return p.range(); }
585 );
586}
587
588
590{
591 return mesh_.nInternalFaces();
592}
593
594
596{
597 return mesh_.nBoundaryFaces();
598}
599
600
602{
603 return labelRange(mesh_.nInternalFaces(), mesh_.nBoundaryFaces());
604}
605
606
608{
609 if (patchi < 0)
610 {
611 return labelRange(mesh_.nInternalFaces(), 0);
612 }
613
614 // Will fail if patchi >= size()
615 return (*this)[patchi].range();
616}
617
618
620(
621 const wordRe& matcher,
622 const bool useGroups
623) const
624{
625 if (matcher.empty())
626 {
627 return labelList();
628 }
629
630 // Only check groups if requested and they exist
631 const bool checkGroups = (useGroups && this->hasGroupIDs());
632
633 labelHashSet ids;
634
635 if (matcher.isPattern())
636 {
637 if (checkGroups)
638 {
639 const auto& groupLookup = groupPatchIDs();
640 forAllConstIters(groupLookup, iter)
641 {
642 if (matcher.match(iter.key()))
643 {
644 // Hash ids associated with the group
645 ids.insert(iter.val());
646 }
647 }
648 }
649
650 if (ids.empty())
651 {
652 return PtrListOps::findMatching(*this, matcher);
653 }
654 else
655 {
656 ids.insert(PtrListOps::findMatching(*this, matcher));
657 }
658 }
659 else
660 {
661 // Literal string.
662 // Special version of above for reduced memory footprint.
663
664 const label patchId = PtrListOps::firstMatching(*this, matcher);
665
666 if (patchId >= 0)
667 {
668 return labelList(one{}, patchId);
669 }
670 else if (checkGroups)
671 {
672 const auto iter = groupPatchIDs().cfind(matcher);
673
674 if (iter.found())
675 {
676 // Hash ids associated with the group
677 ids.insert(iter.val());
678 }
679 }
680 }
681
682 return ids.sortedToc();
683}
684
685
687(
688 const wordRes& matcher,
689 const bool useGroups
690) const
691{
692 if (matcher.empty())
693 {
694 return labelList();
695 }
696 else if (matcher.size() == 1)
697 {
698 return this->indices(matcher.first(), useGroups);
699 }
700
701 labelHashSet ids;
702
703 // Only check groups if requested and they exist
704 if (useGroups && this->hasGroupIDs())
705 {
706 ids.resize(2*this->size());
707
708 const auto& groupLookup = groupPatchIDs();
709 forAllConstIters(groupLookup, iter)
710 {
711 if (matcher.match(iter.key()))
712 {
713 // Hash ids associated with the group
714 ids.insert(iter.val());
715 }
716 }
717 }
718
719 if (ids.empty())
720 {
721 return PtrListOps::findMatching(*this, matcher);
722 }
723 else
724 {
725 ids.insert(PtrListOps::findMatching(*this, matcher));
726 }
727
728 return ids.sortedToc();
729}
730
731
732Foam::label Foam::polyBoundaryMesh::findIndex(const wordRe& key) const
733{
734 if (key.empty())
735 {
736 return -1;
737 }
738 return PtrListOps::firstMatching(*this, key);
739}
740
741
743(
744 const word& patchName,
745 bool allowNotFound
746) const
747{
748 if (patchName.empty())
749 {
750 return -1;
751 }
752
753 const label patchId = PtrListOps::firstMatching(*this, patchName);
754
755 if (patchId >= 0)
756 {
757 return patchId;
758 }
759
760 if (!allowNotFound)
761 {
763 << "Patch '" << patchName << "' not found. "
764 << "Available patch names";
765
766 if (polyMesh::defaultRegion != mesh_.name())
767 {
769 << " in region '" << mesh_.name() << "'";
770 }
771
773 << " include: " << names() << endl
774 << exit(FatalError);
775 }
776
777 // Patch not found
778 if (debug)
779 {
780 Pout<< "label polyBoundaryMesh::findPatchID(const word&) const"
781 << "Patch named " << patchName << " not found. "
782 << "Available patch names: " << names() << endl;
783 }
784
785 // Not found, return -1
786 return -1;
787}
788
789
790Foam::label Foam::polyBoundaryMesh::whichPatch(const label faceIndex) const
791{
792 // Find out which patch the current face belongs to by comparing label
793 // with patch start labels.
794 // If the face is internal, return -1;
795 // if it is off the end of the list, abort
796 if (faceIndex < mesh().nInternalFaces())
797 {
798 return -1;
799 }
800 else if (faceIndex >= mesh().nFaces())
801 {
803 << "Face " << faceIndex
804 << " out of bounds. Number of geometric faces " << mesh().nFaces()
805 << abort(FatalError);
806 }
807
808
809 // Patches are ordered, use binary search
810
811 const polyPatchList& patches = *this;
812
813 const label patchi =
815 (
816 patches,
817 faceIndex,
818 0,
819 // Must include the start in the comparison
820 [](const polyPatch& p, label val) { return (p.start() <= val); }
821 );
822
823 if (patchi < 0 || !patches[patchi].range().found(faceIndex))
824 {
825 // If not in any of above, it is trouble!
827 << "Face " << faceIndex << " not found in any of the patches "
828 << flatOutput(names()) << nl
829 << "The patches appear to be inconsistent with the mesh :"
830 << " internalFaces:" << mesh().nInternalFaces()
831 << " total number of faces:" << mesh().nFaces()
832 << abort(FatalError);
833
834 return -1;
835 }
836
837 return patchi;
838}
839
840
842(
844 const bool warnNotFound,
845 const bool useGroups
846) const
847{
848 const wordList allPatchNames(this->names());
849 labelHashSet ids(2*this->size());
850
851 // Only check groups if requested and they exist
852 const bool checkGroups = (useGroups && this->hasGroupIDs());
853
854 for (const wordRe& matcher : patchNames)
855 {
856 labelList matchIndices = findMatchingStrings(matcher, allPatchNames);
857 ids.insert(matchIndices);
858
859 bool missed = matchIndices.empty();
860
861 if (missed && checkGroups)
862 {
863 // Check group names
864 if (matcher.isPattern())
865 {
866 forAllConstIters(groupPatchIDs(), iter)
867 {
868 if (matcher.match(iter.key()))
869 {
870 // Hash ids associated with the group
871 ids.insert(iter.val());
872 missed = false;
873 }
874 }
875 }
876 else
877 {
878 const auto iter = groupPatchIDs().cfind(matcher);
879
880 if (iter.found())
881 {
882 // Hash ids associated with the group
883 ids.insert(iter.val());
884 missed = false;
885 }
886 }
887 }
888
889 if (missed && warnNotFound)
890 {
891 if (checkGroups)
892 {
894 << "Cannot find any patch or group names matching "
895 << matcher << endl;
896 }
897 else
898 {
900 << "Cannot find any patch names matching "
901 << matcher << endl;
902 }
903 }
904 }
905
906 return ids;
907}
908
909
911(
912 const labelUList& patchIDs,
913 wordList& groups,
914 labelHashSet& nonGroupPatches
915) const
916{
917 // Current matched groups
918 DynamicList<word> matchedGroups(1);
919
920 // Current set of unmatched patches
921 nonGroupPatches = labelHashSet(patchIDs);
922
923 const HashTable<labelList>& groupLookup = this->groupPatchIDs();
924 forAllConstIters(groupLookup, iter)
925 {
926 // Store currently unmatched patches so we can restore
927 labelHashSet oldNonGroupPatches(nonGroupPatches);
928
929 // Match by deleting patches in group from the current set and seeing
930 // if all have been deleted.
931 labelHashSet groupPatchSet(iter.val());
932
933 label nMatch = nonGroupPatches.erase(groupPatchSet);
934
935 if (nMatch == groupPatchSet.size())
936 {
937 matchedGroups.append(iter.key());
938 }
939 else if (nMatch != 0)
940 {
941 // No full match. Undo.
942 nonGroupPatches.transfer(oldNonGroupPatches);
943 }
944 }
945
946 groups.transfer(matchedGroups);
947}
948
949
950bool Foam::polyBoundaryMesh::checkParallelSync(const bool report) const
951{
952 if (!Pstream::parRun())
953 {
954 return false;
955 }
956
957 const polyBoundaryMesh& bm = *this;
958
959 bool hasError = false;
960
961 // Collect non-proc patches and check proc patches are last.
962 wordList localNames(bm.size());
963 wordList localTypes(bm.size());
964
965 label nonProci = 0;
966
967 forAll(bm, patchi)
968 {
969 if (!isA<processorPolyPatch>(bm[patchi]))
970 {
971 if (nonProci != patchi)
972 {
973 // A processor patch in between normal patches!
974 hasError = true;
975
976 if (debug || report)
977 {
978 Pout<< " ***Problem with boundary patch " << patchi
979 << " name:" << bm[patchi].name()
980 << " type:" << bm[patchi].type()
981 << " - seems to be preceeded by processor patches."
982 << " This is usually a problem." << endl;
983 }
984 }
985 else
986 {
987 localNames[nonProci] = bm[patchi].name();
988 localTypes[nonProci] = bm[patchi].type();
989 ++nonProci;
990 }
991 }
992 }
993 localNames.resize(nonProci);
994 localTypes.resize(nonProci);
995
996 // Check and report error(s) on master
997
999 (
1000 // Don't need to collect master itself
1001 (Pstream::master() ? 0 : nonProci),
1003 );
1004
1005 const wordList allNames(procAddr.gather(localNames));
1006 const wordList allTypes(procAddr.gather(localTypes));
1007
1008 // Automatically restricted to master
1009 for (const int proci : procAddr.subProcs())
1010 {
1011 const auto procNames(allNames.slice(procAddr.range(proci)));
1012 const auto procTypes(allTypes.slice(procAddr.range(proci)));
1013
1014 if (procNames != localNames || procTypes != localTypes)
1015 {
1016 hasError = true;
1017
1018 if (debug || report)
1019 {
1020 Info<< " ***Inconsistent patches across processors, "
1021 "processor0 has patch names:" << localNames
1022 << " patch types:" << localTypes
1023 << " processor" << proci
1024 << " has patch names:" << procNames
1025 << " patch types:" << procTypes
1026 << endl;
1027 }
1028 }
1029 }
1030
1031 // Reduce (not broadcast) to respect local out-of-order errors (first loop)
1032 reduce(hasError, orOp<bool>());
1033
1034 return hasError;
1035}
1036
1037
1038bool Foam::polyBoundaryMesh::checkDefinition(const bool report) const
1039{
1040 label nextPatchStart = mesh().nInternalFaces();
1041 const polyBoundaryMesh& bm = *this;
1042
1043 bool hasError = false;
1044
1045 wordHashSet patchNames(2*this->size());
1046
1047 forAll(bm, patchi)
1048 {
1049 if (bm[patchi].start() != nextPatchStart && !hasError)
1050 {
1051 hasError = true;
1052
1053 Info<< " ****Problem with boundary patch " << patchi
1054 << " named " << bm[patchi].name()
1055 << " of type " << bm[patchi].type()
1056 << ". The patch should start on face no " << nextPatchStart
1057 << " and the patch specifies " << bm[patchi].start()
1058 << "." << endl
1059 << "Possibly consecutive patches have this same problem."
1060 << " Suppressing future warnings." << endl;
1061 }
1062
1063 if (!patchNames.insert(bm[patchi].name()) && !hasError)
1064 {
1065 hasError = true;
1066
1067 Info<< " ****Duplicate boundary patch " << patchi
1068 << " named " << bm[patchi].name()
1069 << " of type " << bm[patchi].type()
1070 << "." << endl
1071 << "Suppressing future warnings." << endl;
1072 }
1073
1074 nextPatchStart += bm[patchi].size();
1075 }
1076
1077 reduce(hasError, orOp<bool>());
1078
1079 if (debug || report)
1080 {
1081 if (hasError)
1082 {
1083 Pout<< " ***Boundary definition is in error." << endl;
1084 }
1085 else
1086 {
1087 Info<< " Boundary definition OK." << endl;
1088 }
1089 }
1090
1091 return hasError;
1092}
1093
1094
1096{
1098
1099 if
1100 (
1103 )
1104 {
1105 forAll(*this, patchi)
1106 {
1107 operator[](patchi).initMovePoints(pBufs, p);
1108 }
1109
1110 pBufs.finishedSends();
1111
1112 forAll(*this, patchi)
1113 {
1114 operator[](patchi).movePoints(pBufs, p);
1115 }
1116 }
1117 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1118 {
1119 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1120
1121 // Dummy.
1122 pBufs.finishedSends();
1123
1124 for (const auto& schedEval : patchSchedule)
1125 {
1126 const label patchi = schedEval.patch;
1127
1128 if (schedEval.init)
1129 {
1130 operator[](patchi).initMovePoints(pBufs, p);
1131 }
1132 else
1133 {
1134 operator[](patchi).movePoints(pBufs, p);
1135 }
1136 }
1137 }
1138}
1139
1140
1142{
1143 neighbourEdgesPtr_.clear();
1144 patchIDPtr_.clear();
1145 groupIDsPtr_.clear();
1146
1148
1149 if
1150 (
1153 )
1154 {
1155 forAll(*this, patchi)
1156 {
1157 operator[](patchi).initUpdateMesh(pBufs);
1158 }
1159
1160 pBufs.finishedSends();
1161
1162 forAll(*this, patchi)
1163 {
1164 operator[](patchi).updateMesh(pBufs);
1165 }
1166 }
1167 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1168 {
1169 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1170
1171 // Dummy.
1172 pBufs.finishedSends();
1173
1174 for (const auto& schedEval : patchSchedule)
1175 {
1176 const label patchi = schedEval.patch;
1177
1178 if (schedEval.init)
1179 {
1180 operator[](patchi).initUpdateMesh(pBufs);
1181 }
1182 else
1183 {
1184 operator[](patchi).updateMesh(pBufs);
1185 }
1186 }
1187 }
1188}
1189
1190
1192(
1193 const labelUList& oldToNew,
1194 const bool validBoundary
1195)
1196{
1197 // Change order of patches
1198 polyPatchList::reorder(oldToNew);
1199
1200 // Adapt indices
1201 polyPatchList& patches = *this;
1202
1203 forAll(patches, patchi)
1204 {
1205 patches[patchi].index() = patchi;
1206 }
1207
1208 if (validBoundary)
1209 {
1210 updateMesh();
1211 }
1212}
1213
1214
1216{
1217 const polyPatchList& patches = *this;
1218
1219 os << patches.size() << nl << token::BEGIN_LIST << incrIndent << nl;
1220
1221 for (const polyPatch& pp : patches)
1222 {
1223 os.beginBlock(pp.name());
1224 os << pp;
1225 os.endBlock();
1226 }
1227
1229
1231 return os.good();
1232}
1233
1234
1236(
1237 IOstreamOption streamOpt,
1238 const bool valid
1239) const
1240{
1242 return regIOobject::writeObject(streamOpt, valid);
1243}
1244
1245
1246// * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
1247
1249(
1250 const word& patchName
1251) const
1252{
1253 const label patchi = findPatchID(patchName);
1254
1255 if (patchi < 0)
1256 {
1258 << "Patch named " << patchName << " not found." << nl
1259 << "Available patch names: " << names() << endl
1260 << abort(FatalError);
1261 }
1262
1263 return operator[](patchi);
1264}
1265
1266
1268(
1269 const word& patchName
1270)
1271{
1272 const label patchi = findPatchID(patchName);
1273
1274 if (patchi < 0)
1275 {
1277 << "Patch named " << patchName << " not found." << nl
1278 << "Available patch names: " << names() << endl
1279 << abort(FatalError);
1280 }
1281
1282 return operator[](patchi);
1283}
1284
1285
1286// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1287
1289{
1290 pbm.writeData(os);
1291 return os;
1292}
1293
1294
1295// ************************************************************************* //
scalar range
Functions to operate on Pointer Lists.
bool found
globalIndex procAddr(aMesh.nFaces())
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
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:54
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTableI.H:59
void resize(const label sz)
Resize the hash table for efficiency.
Definition: HashTable.C:601
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
Definition: HashTable.C:719
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition: HashTableI.H:114
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:440
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
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:459
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:180
The IOstreamOption is a simple container for options an IOstream can normally have.
compressionType compression() const noexcept
Get the stream compression.
@ UNCOMPRESSED
compression = false
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void transfer(List< T > &list)
Definition: List.C:447
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:105
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:87
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 sends as done.
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
void append(T *ptr)
Append an element to the end of the list.
Definition: PtrListI.H:113
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:103
const labelUList & indices() const noexcept
Return the list of sorted indices (updated every sort).
Definition: SortListI.H:54
A List obtained as a section of another List.
Definition: SubList.H:70
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:265
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
Definition: SubList.H:165
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
const T * set(const label i) const
Definition: UPtrList.H:248
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
labelRange range() const
Return start/size range of local processor data.
Definition: globalIndexI.H:232
static void gather(const labelUList &offsets, const label comm, const ProcIDsContainer &procIDs, const UList< Type > &fld, List< Type > &allFld, const int tag=UPstream::msgType(), const UPstream::commsTypes=UPstream::commsTypes::nonBlocking)
Collect data in processor order on master (== procIDs[0]).
labelRange subProcs() const noexcept
Range of process indices for addressed sub-offsets (processes)
Definition: globalIndexI.H:159
const lduSchedule & patchSchedule() const noexcept
const labelHashSet & patchSet() const
Wall patches to process.
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:58
void movePoints()
Update for new mesh geometry.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
bool checkDefinition(const bool report=false) const
Check boundary definition.
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.
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
writeData member function required by regIOobject
void setGroup(const word &groupName, const labelUList &patchIDs)
Set/add group with patches.
labelRange range() const noexcept
The face range for all boundary faces.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
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.
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
const labelList & patchID() const
Per boundary face label the patch index.
labelList patchSizes() const
Return a list of patch sizes.
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.
label findIndex(const wordRe &key) const
Return patch index for the first match, return -1 if not found.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:321
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
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:76
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
splitCell * master() const
Definition: splitCell.H:113
@ BEGIN_LIST
Begin list [isseparator].
Definition: token.H:155
@ END_LIST
End list [isseparator].
Definition: token.H:156
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:83
bool isPattern() const noexcept
The wordRe is a pattern, not a literal string.
Definition: wordReI.H:111
bool match(const std::string &text, bool literal=false) const
Smart match as regular expression or as a string.
Definition: wordReI.H:200
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
bool match(const std::string &text, bool literal=false) const
Smart match as literal or regex, stopping on the first match.
Definition: wordResI.H:91
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
label patchId(-1)
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
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 OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:57
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:63
List< word > wordList
A List of words.
Definition: fileName.H:63
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
labelList findMatchingStrings(const UnaryMatchPredicate &matcher, const UList< StringType > &input, const bool invert=false)
Extract list indices for all matches.
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:50
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:349
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
errorManip< error > abort(error &err)
Definition: errorManip.H:144
const direction noexcept
Definition: Scalar.H:223
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:356
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
wordList patchNames(nPatches)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Operations on lists of strings.
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:238
Extract type (as a word) from an object, typically using its type() method.
Definition: word.H:255