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