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-2023 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 "wordRes.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_.emplace(16);
81  auto& groupLookup = *groupIDsPtr_;
82 
83  const polyPatchList& patches = *this;
84 
85  forAll(patches, patchi)
86  {
87  for (const word& groupName : patches[patchi].inGroups())
88  {
89  groupLookup(groupName).push_back(patchi);
90  }
91  }
92 
93  // Remove groups that clash with patch names
94  forAll(patches, patchi)
95  {
96  if (groupLookup.erase(patches[patchi].name()))
97  {
99  << "Removed group '" << patches[patchi].name()
100  << "' which clashes with patch " << patchi
101  << " of the same name."
102  << endl;
103  }
104  }
105 }
106 
107 
108 void Foam::polyBoundaryMesh::populate(PtrList<entry>&& entries)
109 {
110  clearLocalAddressing();
111 
112  polyPatchList& patches = *this;
113 
114  patches.resize_null(entries.size());
115 
116  // Transcribe.
117  // Does not handle nullptr at all (what could possibly be done?)
118  forAll(patches, patchi)
119  {
120  patches.set
121  (
122  patchi,
124  (
125  entries[patchi].keyword(),
126  entries[patchi].dict(),
127  patchi,
128  *this
129  )
130  );
131  }
132 
133  entries.clear();
134 }
135 
136 
137 bool Foam::polyBoundaryMesh::readContents(const bool allowOptionalRead)
138 {
139  bool updated = false;
140  PtrList<entry> entries;
141 
142  if
143  (
144  this->isReadRequired()
145  || (allowOptionalRead && this->isReadOptional() && this->headerOk())
146  )
147  {
148  // Warn for MUST_READ_IF_MODIFIED
149  warnNoRereading<polyBoundaryMesh>();
150 
151  // Read entries
152  Istream& is = readStream(typeName);
153 
154  is >> entries;
155 
156  is.check(FUNCTION_NAME);
157  close();
158  updated = true;
159  }
160  // Future: support master-only and broadcast?
161 
162  if (updated)
163  {
164  populate(std::move(entries));
165  }
166 
167  return updated;
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
172 
174 (
175  const IOobject& io,
176  const polyMesh& mesh
177 )
178 :
179  polyPatchList(),
180  regIOobject(io),
181  mesh_(mesh)
182 {
183  readContents(false); // allowOptionalRead = false
184 }
185 
186 
188 (
189  const IOobject& io,
190  const polyMesh& pm,
191  Foam::zero
192 )
193 :
195  regIOobject(io),
196  mesh_(pm)
197 {}
198 
199 
201 (
202  const IOobject& io,
203  const polyMesh& pm,
204  const label size
205 )
206 :
208  regIOobject(io),
209  mesh_(pm)
210 {}
211 
212 
214 (
215  const IOobject& io,
216  const polyMesh& pm,
217  const polyPatchList& list
218 )
219 :
220  polyPatchList(),
221  regIOobject(io),
222  mesh_(pm)
223 {
224  if (!readContents(true)) // allowOptionalRead = true
225  {
226  // Nothing read. Use supplied patches
227  polyPatchList& patches = *this;
228  patches.resize(list.size());
229 
230  forAll(patches, patchi)
231  {
232  patches.set(patchi, list[patchi].clone(*this));
233  }
234  }
235 }
236 
237 
239 (
240  const IOobject& io,
241  const polyMesh& pm,
242  PtrList<entry>&& entries
243 )
244 :
245  polyPatchList(),
246  regIOobject(io),
247  mesh_(pm)
248 {
249  if (!readContents(true)) // allowOptionalRead = true
250  {
251  populate(std::move(entries));
252  }
253  entries.clear();
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
258 
260 {
262  clearAddressing();
263 }
264 
265 
267 {
268  polyPatchList& patches = *this;
269 
270  for (polyPatch& p : patches)
271  {
272  p.clearGeom();
273  }
274 }
275 
276 
277 // Private until it is more generally required (and gets a better name?)
278 void Foam::polyBoundaryMesh::clearLocalAddressing()
279 {
280  neighbourEdgesPtr_.reset(nullptr);
281  patchIDPtr_.reset(nullptr);
282  groupIDsPtr_.reset(nullptr);
283 }
284 
285 
287 {
288  clearLocalAddressing();
289 
290  polyPatchList& patches = *this;
291 
292  for (polyPatch& p : patches)
293  {
294  p.clearAddressing();
295  }
296 }
297 
298 
299 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
300 
301 void Foam::polyBoundaryMesh::calcGeometry()
302 {
303  PstreamBuffers pBufs(Pstream::defaultCommsType);
304 
305  if
306  (
307  pBufs.commsType() == Pstream::commsTypes::blocking
308  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
309  )
310  {
311  forAll(*this, patchi)
312  {
313  operator[](patchi).initGeometry(pBufs);
314  }
315 
316  pBufs.finishedSends();
317 
318  forAll(*this, patchi)
319  {
320  operator[](patchi).calcGeometry(pBufs);
321  }
322  }
323  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
324  {
325  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
326 
327  // Dummy.
328  pBufs.finishedSends();
329 
330  for (const auto& schedEval : patchSchedule)
331  {
332  const label patchi = schedEval.patch;
333 
334  if (schedEval.init)
335  {
336  operator[](patchi).initGeometry(pBufs);
337  }
338  else
339  {
340  operator[](patchi).calcGeometry(pBufs);
341  }
342  }
343  }
344 }
345 
346 
348 {
349  return faceList::subList
350  (
351  mesh_.faces(),
352  mesh_.nBoundaryFaces(),
353  mesh_.nInternalFaces()
354  );
355 }
356 
357 
359 {
360  return labelList::subList
361  (
362  mesh_.faceOwner(),
363  mesh_.nBoundaryFaces(),
364  mesh_.nInternalFaces()
365  );
366 }
367 
368 // Potentially useful to simplify logic elsewhere?
369 // const Foam::labelList::subList Foam::polyBoundaryMesh::faceNeighbour() const
370 // {
371 // return labelList::subList();
372 // }
373 
374 
377 {
378  const polyPatchList& patches = *this;
379 
381 
382  forAll(patches, patchi)
383  {
384  list.set(patchi, &patches[patchi].faceCells());
385  }
386 
387  return list;
388 }
389 
390 
393 {
394  if (Pstream::parRun())
395  {
397  << "Neighbour edge addressing not correct across parallel"
398  << " boundaries." << endl;
399  }
400 
401  if (!neighbourEdgesPtr_)
402  {
403  neighbourEdgesPtr_.emplace(size());
404  auto& neighbourEdges = *neighbourEdgesPtr_;
405 
406  // Initialize.
407  label nEdgePairs = 0;
408  forAll(*this, patchi)
409  {
410  const polyPatch& pp = operator[](patchi);
411 
412  neighbourEdges[patchi].setSize(pp.nEdges() - pp.nInternalEdges());
413 
414  for (labelPair& edgeInfo : neighbourEdges[patchi])
415  {
416  edgeInfo[0] = -1;
417  edgeInfo[1] = -1;
418  }
419 
420  nEdgePairs += pp.nEdges() - pp.nInternalEdges();
421  }
422 
423  // From mesh edge (expressed as a point pair so as not to construct
424  // point addressing) to patch + relative edge index.
425  EdgeMap<labelPair> pointsToEdge(nEdgePairs);
426 
427  forAll(*this, patchi)
428  {
429  const polyPatch& pp = operator[](patchi);
430 
431  const edgeList& edges = pp.edges();
432 
433  for
434  (
435  label edgei = pp.nInternalEdges();
436  edgei < edges.size();
437  edgei++
438  )
439  {
440  // Edge in patch local points
441  const edge& e = edges[edgei];
442 
443  // Edge in mesh points.
444  edge meshEdge(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
445 
446  auto fnd = pointsToEdge.find(meshEdge);
447 
448  if (!fnd.good())
449  {
450  // First occurrence of mesh edge. Store patch and my
451  // local index.
452  pointsToEdge.insert
453  (
454  meshEdge,
455  labelPair
456  (
457  patchi,
458  edgei - pp.nInternalEdges()
459  )
460  );
461  }
462  else
463  {
464  // Second occurrence. Store.
465  const labelPair& edgeInfo = fnd.val();
466 
467  neighbourEdges[patchi][edgei - pp.nInternalEdges()] =
468  edgeInfo;
469 
470  neighbourEdges[edgeInfo[0]][edgeInfo[1]]
471  = labelPair(patchi, edgei - pp.nInternalEdges());
472 
473  // Found all two occurrences of this edge so remove from
474  // hash to save space. Note that this will give lots of
475  // problems if the polyBoundaryMesh is multiply connected.
476  pointsToEdge.erase(meshEdge);
477  }
478  }
479  }
480 
481  if (pointsToEdge.size())
482  {
484  << "Not all boundary edges of patches match up." << nl
485  << "Is the outside of your mesh multiply connected?"
486  << abort(FatalError);
487  }
488 
489  forAll(*this, patchi)
490  {
491  const polyPatch& pp = operator[](patchi);
492 
493  const labelPairList& nbrEdges = neighbourEdges[patchi];
494 
495  forAll(nbrEdges, i)
496  {
497  const labelPair& edgeInfo = nbrEdges[i];
498 
499  if (edgeInfo[0] == -1 || edgeInfo[1] == -1)
500  {
501  const label edgei = pp.nInternalEdges() + i;
502  const edge& e = pp.edges()[edgei];
503 
505  << "Not all boundary edges of patches match up." << nl
506  << "Edge " << edgei << " on patch " << pp.name()
507  << " end points " << pp.localPoints()[e[0]] << ' '
508  << pp.localPoints()[e[1]] << " is not matched to an"
509  << " edge on any other patch." << nl
510  << "Is the outside of your mesh multiply connected?"
511  << abort(FatalError);
512  }
513  }
514  }
515  }
516 
517  return *neighbourEdgesPtr_;
518 }
519 
520 
522 {
523  if (!patchIDPtr_)
524  {
525  patchIDPtr_.emplace(mesh_.nBoundaryFaces());
526  auto& list = *patchIDPtr_;
527 
528  const polyPatchList& patches = *this;
529 
530  forAll(patches, patchi)
531  {
532  SubList<label>
533  (
534  list,
535  patches[patchi].size(),
536  (patches[patchi].start() - mesh_.nInternalFaces())
537  ) = patchi;
538  }
539  }
540 
541  return *patchIDPtr_;
542 }
543 
544 
545 Foam::label Foam::polyBoundaryMesh::patchID(const label meshFacei) const
546 {
547  const label bndFacei = (meshFacei - mesh_.nInternalFaces());
548 
549  return
550  (
551  (bndFacei >= 0 && bndFacei < mesh_.nBoundaryFaces())
552  ? this->patchID()[bndFacei]
553  : -1
554  );
555 }
556 
557 
559 Foam::polyBoundaryMesh::patchID(const labelUList& meshFaceIndices) const
560 {
561  labelList output(meshFaceIndices.size());
562  forAll(meshFaceIndices, i)
563  {
564  output[i] = patchID(meshFaceIndices[i]);
565  }
566  return output;
567 }
568 
569 
572 {
573  if (!groupIDsPtr_)
574  {
575  calcGroupIDs();
576  }
577 
578  return *groupIDsPtr_;
579 }
580 
581 
583 (
584  const word& groupName,
585  const labelUList& patchIDs
586 )
587 {
588  groupIDsPtr_.reset(nullptr);
589 
590  polyPatchList& patches = *this;
591 
592  boolList pending(patches.size(), true);
593 
594  // Add to specified patches
595  for (const label patchi : patchIDs)
596  {
597  if (pending.test(patchi))
598  {
599  pending.unset(patchi);
600  patches[patchi].addGroup(groupName);
601  }
602  }
603 
604  // Remove from other patches
605  forAll(patches, patchi)
606  {
607  if (pending.test(patchi))
608  {
609  patches[patchi].removeGroup(groupName);
610  }
611  }
612 }
613 
614 
615 Foam::label Foam::polyBoundaryMesh::nNonProcessor() const
616 {
617  const polyPatchList& patches = *this;
618 
619  label count = 0;
620 
621  for (const polyPatch& p : patches)
622  {
623  if (isA<processorPolyPatch>(p))
624  {
625  break;
626  }
627 
628  ++count;
629  }
630 
631  return count;
632 }
633 
634 
636 {
637  const polyPatchList& patches = *this;
638 
639  label count = 0;
640 
641  for (const polyPatch& p : patches)
642  {
643  if (isA<processorPolyPatch>(p))
644  {
645  ++count;
646  }
647  }
648 
649  return count;
650 }
651 
654 {
655  return PtrListOps::get<word>(*this, nameOp<polyPatch>());
656 }
657 
660 {
661  return PtrListOps::get<word>(*this, typeOp<polyPatch>());
662 }
663 
664 
666 {
667  return
668  PtrListOps::get<word>
669  (
670  *this,
671  [](const polyPatch& p) { return p.physicalType(); }
672  );
673 }
674 
675 
677 {
678  return
679  PtrListOps::get<label>
680  (
681  *this,
682  [](const polyPatch& p) { return p.start(); }
683  );
684 }
685 
686 
688 {
689  return
690  PtrListOps::get<label>
691  (
692  *this,
693  [](const polyPatch& p) { return p.size(); }
694  );
695 }
696 
697 
699 {
700  return
701  PtrListOps::get<labelRange>
702  (
703  *this,
704  [](const polyPatch& p) { return p.range(); }
705  );
706 }
707 
710 {
711  return this->groupPatchIDs().sortedToc();
712 }
713 
715 Foam::label Foam::polyBoundaryMesh::start() const noexcept
716 {
717  return mesh_.nInternalFaces();
718 }
719 
721 Foam::label Foam::polyBoundaryMesh::nFaces() const noexcept
722 {
723  return mesh_.nBoundaryFaces();
724 }
725 
728 {
729  return labelRange(mesh_.nInternalFaces(), mesh_.nBoundaryFaces());
730 }
731 
732 
733 Foam::labelRange Foam::polyBoundaryMesh::range(const label patchi) const
734 {
735  if (patchi < 0)
736  {
737  return labelRange(mesh_.nInternalFaces(), 0);
738  }
740  // Will fail if patchi >= size()
741  return (*this)[patchi].range();
742 }
743 
744 
746 (
747  const wordRe& matcher,
748  const bool useGroups
749 ) const
750 {
751  if (matcher.empty())
752  {
753  return labelList();
754  }
755 
756  // Only check groups if requested and they exist
757  const bool checkGroups = (useGroups && this->hasGroupIDs());
758 
759  labelHashSet ids(0);
760 
761  if (matcher.isPattern())
762  {
763  if (checkGroups)
764  {
765  const auto& groupLookup = groupPatchIDs();
766  forAllConstIters(groupLookup, iter)
767  {
768  if (matcher(iter.key()))
769  {
770  // Add patch ids associated with the group
771  ids.insert(iter.val());
772  }
773  }
774  }
775 
776  if (ids.empty())
777  {
778  return PtrListOps::findMatching(*this, matcher);
779  }
780  else
781  {
782  ids.insert(PtrListOps::findMatching(*this, matcher));
783  }
784  }
785  else
786  {
787  // Literal string.
788  // Special version of above for reduced memory footprint.
789 
790  const label patchId = PtrListOps::firstMatching(*this, matcher);
791 
792  if (patchId >= 0)
793  {
794  return labelList(one{}, patchId);
795  }
796  else if (checkGroups)
797  {
798  const auto iter = groupPatchIDs().cfind(matcher);
799 
800  if (iter.good())
801  {
802  // Hash ids associated with the group
803  ids.insert(iter.val());
804  }
805  }
806  }
807 
808  return ids.sortedToc();
809 }
810 
811 
813 (
814  const wordRes& matcher,
815  const bool useGroups
816 ) const
817 {
818  if (matcher.empty())
819  {
820  return labelList();
821  }
822  else if (matcher.size() == 1)
823  {
824  return this->indices(matcher.front(), useGroups);
825  }
826 
827  labelHashSet ids(0);
828 
829  // Only check groups if requested and they exist
830  if (useGroups && this->hasGroupIDs())
831  {
832  ids.reserve(this->size());
833 
834  const auto& groupLookup = groupPatchIDs();
835  forAllConstIters(groupLookup, iter)
836  {
837  if (matcher(iter.key()))
838  {
839  // Add patch ids associated with the group
840  ids.insert(iter.val());
841  }
842  }
843  }
844 
845  if (ids.empty())
846  {
847  return PtrListOps::findMatching(*this, matcher);
848  }
849  else
850  {
851  ids.insert(PtrListOps::findMatching(*this, matcher));
852  }
853 
854  return ids.sortedToc();
855 }
856 
857 
859 (
860  const wordRes& select,
861  const wordRes& ignore,
862  const bool useGroups
863 ) const
864 {
865  if (ignore.empty())
866  {
867  return this->indices(select, useGroups);
868  }
869 
870  const wordRes::filter matcher(select, ignore);
871 
872  labelHashSet ids(0);
873 
874  // Only check groups if requested and they exist
875  if (useGroups && this->hasGroupIDs())
876  {
877  ids.reserve(this->size());
878 
879  const auto& groupLookup = groupPatchIDs();
880  forAllConstIters(groupLookup, iter)
881  {
882  if (matcher(iter.key()))
883  {
884  // Add patch ids associated with the group
885  ids.insert(iter.val());
886  }
887  }
888  }
889 
890  if (ids.empty())
891  {
892  return PtrListOps::findMatching(*this, matcher);
893  }
894  else
895  {
896  ids.insert(PtrListOps::findMatching(*this, matcher));
897  }
898 
899  return ids.sortedToc();
900 }
901 
902 
903 Foam::label Foam::polyBoundaryMesh::findIndex(const wordRe& key) const
904 {
905  if (key.empty())
906  {
907  return -1;
908  }
909  return PtrListOps::firstMatching(*this, key);
910 }
911 
912 
914 (
915  const word& patchName,
916  bool allowNotFound
917 ) const
918 {
919  if (patchName.empty())
920  {
921  return -1;
922  }
923 
924  const label patchId = PtrListOps::firstMatching(*this, patchName);
925 
926  if (patchId >= 0)
927  {
928  return patchId;
929  }
930 
931  if (!allowNotFound)
932  {
934  << "Patch '" << patchName << "' not found. "
935  << "Available patch names";
936 
937  if (polyMesh::defaultRegion != mesh_.name())
938  {
939  FatalError
940  << " in region '" << mesh_.name() << "'";
941  }
942 
943  FatalError
944  << " include: " << names() << endl
945  << exit(FatalError);
946  }
947 
948  // Patch not found
949  if (debug)
950  {
951  Pout<< "label polyBoundaryMesh::findPatchID(const word&) const"
952  << "Patch named " << patchName << " not found. "
953  << "Available patch names: " << names() << endl;
954  }
956  // Not found, return -1
957  return -1;
958 }
959 
960 
962 Foam::polyBoundaryMesh::whichPatchFace(const label meshFacei) const
963 {
964  if (meshFacei < mesh().nInternalFaces())
965  {
966  // Internal face: return (-1, meshFace)
967  return labelPair(-1, meshFacei);
968  }
969  else if (meshFacei >= mesh().nFaces())
970  {
971  // Bounds error: abort
973  << "Face " << meshFacei
974  << " out of bounds. Number of geometric faces " << mesh().nFaces()
975  << abort(FatalError);
976 
977  return labelPair(-1, meshFacei);
978  }
979 
980 
981  // Patches are ordered, use binary search
982  // Find out which patch index and local patch face the specified
983  // mesh face belongs to by comparing label with patch start labels.
984 
985 
986  // TBD: use patchIDPtr_ if it exists?
987 
988  const polyPatchList& patches = *this;
989 
990  const label patchi =
991  findLower
992  (
993  patches,
994  meshFacei,
995  0,
996  // Must include the start in the comparison
997  [](const polyPatch& p, label val) { return (p.start() <= val); }
998  );
999 
1000  if (patchi < 0 || !patches[patchi].range().contains(meshFacei))
1001  {
1002  // If not in any of above, it is trouble!
1004  << "Face " << meshFacei << " not found in any of the patches "
1005  << flatOutput(names()) << nl
1006  << "The patches appear to be inconsistent with the mesh :"
1007  << " internalFaces:" << mesh().nInternalFaces()
1008  << " total number of faces:" << mesh().nFaces()
1009  << abort(FatalError);
1010 
1011  return labelPair(-1, meshFacei);
1012  }
1014  // (patch, local face index)
1015  return labelPair(patchi, meshFacei - patches[patchi].start());
1016 }
1017 
1018 
1020 Foam::polyBoundaryMesh::whichPatchFace(const labelUList& meshFaceIndices) const
1021 {
1022  labelPairList output(meshFaceIndices.size());
1023  forAll(meshFaceIndices, i)
1024  {
1025  output[i] = whichPatchFace(meshFaceIndices[i]);
1026  }
1027  return output;
1028 }
1029 
1030 
1032 (
1033  const UList<wordRe>& select,
1034  const bool warnNotFound,
1035  const bool useGroups
1036 ) const
1037 {
1038  labelHashSet ids(0);
1039  if (select.empty())
1040  {
1041  return ids;
1042  }
1043 
1044  const polyPatchList& patches = *this;
1045 
1046  const label len = patches.size();
1047 
1048  ids.reserve(len);
1049 
1050  // Only check groups if requested and they exist
1051  const bool checkGroups = (useGroups && this->hasGroupIDs());
1052 
1053  for (const wordRe& matcher : select)
1054  {
1055  bool missed = true;
1056 
1057  for (label i = 0; i < len; ++i)
1058  {
1059  if (matcher(patches[i].name()))
1060  {
1061  ids.insert(i);
1062  missed = false;
1063  }
1064  }
1065 
1066  if (missed && checkGroups)
1067  {
1068  // Check group names
1069  if (matcher.isPattern())
1070  {
1071  forAllConstIters(groupPatchIDs(), iter)
1072  {
1073  if (matcher.match(iter.key()))
1074  {
1075  // Hash ids associated with the group
1076  ids.insert(iter.val());
1077  missed = false;
1078  }
1079  }
1080  }
1081  else
1082  {
1083  const auto iter = groupPatchIDs().cfind(matcher);
1084 
1085  if (iter.good())
1086  {
1087  // Hash ids associated with the group
1088  ids.insert(iter.val());
1089  missed = false;
1090  }
1091  }
1092  }
1093 
1094  if (missed && warnNotFound)
1095  {
1096  if (checkGroups)
1097  {
1099  << "Cannot find any patch or group names matching "
1100  << matcher << endl;
1101  }
1102  else
1103  {
1105  << "Cannot find any patch names matching "
1106  << matcher << endl;
1107  }
1108  }
1109  }
1110 
1111  return ids;
1112 }
1113 
1114 
1116 (
1117  const labelUList& patchIDs,
1118  wordList& groups,
1119  labelHashSet& nonGroupPatches
1120 ) const
1121 {
1122  // Current matched groups
1123  DynamicList<word> matchedGroups(1);
1124 
1125  // Current set of unmatched patches
1126  nonGroupPatches = labelHashSet(patchIDs);
1127 
1128  const HashTable<labelList>& groupLookup = this->groupPatchIDs();
1129  forAllConstIters(groupLookup, iter)
1130  {
1131  // Store currently unmatched patches so we can restore
1132  labelHashSet oldNonGroupPatches(nonGroupPatches);
1133 
1134  // Match by deleting patches in group from the current set and seeing
1135  // if all have been deleted.
1136  labelHashSet groupPatchSet(iter.val());
1137 
1138  label nMatch = nonGroupPatches.erase(groupPatchSet);
1139 
1140  if (nMatch == groupPatchSet.size())
1141  {
1142  matchedGroups.push_back(iter.key());
1143  }
1144  else if (nMatch != 0)
1145  {
1146  // No full match. Undo.
1147  nonGroupPatches.transfer(oldNonGroupPatches);
1148  }
1149  }
1150 
1151  groups.transfer(matchedGroups);
1152 }
1153 
1154 
1155 bool Foam::polyBoundaryMesh::checkParallelSync(const bool report) const
1156 {
1157  if (!Pstream::parRun())
1158  {
1159  return false;
1160  }
1161 
1162  const polyBoundaryMesh& bm = *this;
1163 
1164  bool hasError = false;
1165 
1166  // Collect non-proc patches and check proc patches are last.
1167  wordList localNames(bm.size());
1168  wordList localTypes(bm.size());
1169 
1170  label nonProci = 0;
1171 
1172  forAll(bm, patchi)
1173  {
1174  if (!isA<processorPolyPatch>(bm[patchi]))
1175  {
1176  if (nonProci != patchi)
1177  {
1178  // A processor patch in between normal patches!
1179  hasError = true;
1180 
1181  if (debug || report)
1182  {
1183  Pout<< " ***Problem with boundary patch " << patchi
1184  << " name:" << bm[patchi].name()
1185  << " type:" << bm[patchi].type()
1186  << " - seems to be preceeded by processor patches."
1187  << " This is usually a problem." << endl;
1188  }
1189  }
1190  else
1191  {
1192  localNames[nonProci] = bm[patchi].name();
1193  localTypes[nonProci] = bm[patchi].type();
1194  ++nonProci;
1195  }
1196  }
1197  }
1198  localNames.resize(nonProci);
1199  localTypes.resize(nonProci);
1200 
1201  // Check and report error(s) on master
1202  // - don't need indexing on master itself
1203 
1204  const globalIndex procAddr(globalIndex::gatherNonLocal{}, nonProci);
1205 
1206  const wordList allNames(procAddr.gather(localNames));
1207  const wordList allTypes(procAddr.gather(localTypes));
1208 
1209  // Automatically restricted to master
1210  for (const int proci : procAddr.subProcs())
1211  {
1212  const auto procNames(allNames.slice(procAddr.range(proci)));
1213  const auto procTypes(allTypes.slice(procAddr.range(proci)));
1214 
1215  if (procNames != localNames || procTypes != localTypes)
1216  {
1217  hasError = true;
1218 
1219  if (debug || report)
1220  {
1221  Info<< " ***Inconsistent patches across processors, "
1222  "processor0 has patch names:" << localNames
1223  << " patch types:" << localTypes
1224  << " processor" << proci
1225  << " has patch names:" << procNames
1226  << " patch types:" << procTypes
1227  << endl;
1228  }
1229  }
1230  }
1231 
1232  // Reduce (not broadcast) to respect local out-of-order errors (first loop)
1233  return returnReduceOr(hasError);
1234 }
1235 
1236 
1237 bool Foam::polyBoundaryMesh::checkDefinition(const bool report) const
1238 {
1239  label nextPatchStart = mesh().nInternalFaces();
1240  const polyBoundaryMesh& bm = *this;
1241 
1242  bool hasError = false;
1243 
1244  wordHashSet patchNames(2*this->size());
1245 
1246  forAll(bm, patchi)
1247  {
1248  if (bm[patchi].start() != nextPatchStart && !hasError)
1249  {
1250  hasError = true;
1251 
1252  Info<< " ****Problem with boundary patch " << patchi
1253  << " named " << bm[patchi].name()
1254  << " of type " << bm[patchi].type()
1255  << ". The patch should start on face no " << nextPatchStart
1256  << " and the patch specifies " << bm[patchi].start()
1257  << "." << endl
1258  << "Possibly consecutive patches have this same problem."
1259  << " Suppressing future warnings." << endl;
1260  }
1261 
1262  if (!patchNames.insert(bm[patchi].name()) && !hasError)
1263  {
1264  hasError = true;
1265 
1266  Info<< " ****Duplicate boundary patch " << patchi
1267  << " named " << bm[patchi].name()
1268  << " of type " << bm[patchi].type()
1269  << "." << endl
1270  << "Suppressing future warnings." << endl;
1271  }
1272 
1273  nextPatchStart += bm[patchi].size();
1274  }
1275 
1276  Pstream::reduceOr(hasError);
1277 
1278  if (debug || report)
1279  {
1280  if (hasError)
1281  {
1282  Pout<< " ***Boundary definition is in error." << endl;
1283  }
1284  else
1285  {
1286  Info<< " Boundary definition OK." << endl;
1287  }
1288  }
1289 
1290  return hasError;
1291 }
1292 
1293 
1295 {
1296  PstreamBuffers pBufs(Pstream::defaultCommsType);
1297 
1298  if
1299  (
1300  pBufs.commsType() == Pstream::commsTypes::blocking
1301  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
1302  )
1303  {
1304  forAll(*this, patchi)
1305  {
1306  operator[](patchi).initMovePoints(pBufs, p);
1307  }
1308 
1309  pBufs.finishedSends();
1310 
1311  forAll(*this, patchi)
1312  {
1313  operator[](patchi).movePoints(pBufs, p);
1314  }
1315  }
1316  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1317  {
1318  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1319 
1320  // Dummy.
1321  pBufs.finishedSends();
1322 
1323  for (const auto& schedEval : patchSchedule)
1324  {
1325  const label patchi = schedEval.patch;
1326 
1327  if (schedEval.init)
1328  {
1329  operator[](patchi).initMovePoints(pBufs, p);
1330  }
1331  else
1332  {
1333  operator[](patchi).movePoints(pBufs, p);
1334  }
1335  }
1336  }
1337 }
1338 
1339 
1341 {
1342  neighbourEdgesPtr_.reset(nullptr);
1343  patchIDPtr_.reset(nullptr);
1344  groupIDsPtr_.reset(nullptr);
1345 
1346  PstreamBuffers pBufs(Pstream::defaultCommsType);
1347 
1348  if
1349  (
1350  pBufs.commsType() == Pstream::commsTypes::blocking
1351  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
1352  )
1353  {
1354  forAll(*this, patchi)
1355  {
1356  operator[](patchi).initUpdateMesh(pBufs);
1357  }
1358 
1359  pBufs.finishedSends();
1360 
1361  forAll(*this, patchi)
1362  {
1363  operator[](patchi).updateMesh(pBufs);
1364  }
1365  }
1366  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1367  {
1368  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1369 
1370  // Dummy.
1371  pBufs.finishedSends();
1372 
1373  for (const auto& schedEval : patchSchedule)
1374  {
1375  const label patchi = schedEval.patch;
1376 
1377  if (schedEval.init)
1378  {
1379  operator[](patchi).initUpdateMesh(pBufs);
1380  }
1381  else
1382  {
1383  operator[](patchi).updateMesh(pBufs);
1384  }
1385  }
1386  }
1387 }
1388 
1389 
1391 (
1392  const labelUList& oldToNew,
1393  const bool validBoundary
1394 )
1395 {
1396  // Change order of patches
1397  polyPatchList::reorder(oldToNew);
1398 
1399  // Adapt indices
1400  polyPatchList& patches = *this;
1401 
1402  forAll(patches, patchi)
1403  {
1404  patches[patchi].index() = patchi;
1405  }
1406 
1407  if (validBoundary)
1408  {
1409  updateMesh();
1410  }
1411 }
1412 
1413 
1414 void Foam::polyBoundaryMesh::writeEntry(Ostream& os) const
1415 {
1416  const polyPatchList& entries = *this;
1417 
1418  os << entries.size();
1419 
1420  if (entries.empty())
1421  {
1422  // 0-sized : can write with less vertical space
1424  }
1425  else
1426  {
1427  os << nl << token::BEGIN_LIST << incrIndent << nl;
1428 
1429  for (const auto& pp : entries)
1430  {
1431  os.beginBlock(pp.name());
1432  os << pp;
1433  os.endBlock();
1434  }
1436  }
1438 }
1439 
1440 
1442 (
1443  const keyType& keyword,
1444  Ostream& os
1445 ) const
1446 {
1447  const polyPatchList& entries = *this;
1448 
1449  if (!keyword.empty())
1450  {
1451  os.write(keyword);
1452  os << (entries.empty() ? token::SPACE : token::NL);
1453  }
1455  writeEntry(os);
1456 
1457  if (!keyword.empty()) os.endEntry();
1458 }
1459 
1460 
1463  writeEntry(os);
1464  return os.good();
1465 }
1466 
1467 
1469 (
1470  IOstreamOption streamOpt,
1471  const bool writeOnProc
1472 ) const
1473 {
1475  return regIOobject::writeObject(streamOpt, writeOnProc);
1476 }
1477 
1478 
1479 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
1480 
1482 (
1483  const word& patchName
1484 ) const
1485 {
1486  const label patchi = findPatchID(patchName);
1487 
1488  if (patchi < 0)
1489  {
1491  << "Patch named " << patchName << " not found." << nl
1492  << "Available patch names: " << names() << endl
1493  << abort(FatalError);
1494  }
1495 
1496  return operator[](patchi);
1497 }
1498 
1499 
1501 (
1502  const word& patchName
1503 )
1504 {
1505  const label patchi = findPatchID(patchName);
1506 
1507  if (patchi < 0)
1508  {
1510  << "Patch named " << patchName << " not found." << nl
1511  << "Available patch names: " << names() << endl
1512  << abort(FatalError);
1513  }
1515  return operator[](patchi);
1516 }
1517 
1518 
1519 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1520 
1521 Foam::Ostream& Foam::operator<<(Ostream& os, const polyBoundaryMesh& pbm)
1522 {
1523  pbm.writeData(os);
1524  return os;
1525 }
1526 
1527 
1528 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
labelRange range(const label proci) const
Return start/size range of proci data.
Definition: globalIndexI.H:282
label patchId(-1)
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const polyBoundaryMesh & pbm
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:116
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:621
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type unset(const label i)
Unset the bool entry at specified position, always false for out-of-range access. ...
Definition: UList.H:796
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const Field< point_type > & localPoints() const
Return pointField of points in patch.
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
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:598
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
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:195
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
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:52
polyBoundaryMesh(const polyBoundaryMesh &)=delete
No copy construct.
List< bool > select(const label n, const labelUList &locations)
Construct a selection list of bools (all false) with the given pre-size, subsequently add specified l...
Definition: BitOps.C:134
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
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
The writeData member function required by regIOobject.
Newline [isspace].
Definition: token.H:130
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:531
labelPair whichPatchFace(const label meshFacei) const
Lookup mesh face index and return (patchi, patchFacei) tuple or (-1, meshFacei) for internal faces...
const labelList & patchID() const
Per boundary face label the patch index.
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
void reorder(const labelUList &oldToNew, const bool check=false)
Reorder elements. Reordering must be unique (ie, shuffle).
Definition: UPtrList.C:79
wordList groupNames() const
A list of the group names (if any)
label nInternalEdges() const
Number of internal edges.
void clearGeom()
Clear geometry at this level and at patches.
static void reduceOr(bool &value, const label communicator=worldComm)
Logical (or) reduction (MPI_AllReduce)
Begin list [isseparator].
Definition: token.H:161
labelRange range() const noexcept
The face range for all boundary faces.
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:46
List< labelPair > labelPairList
List of labelPair.
Definition: labelPair.H:33
Functions to operate on Pointer Lists.
A simple container for options an IOstream can normally have.
SubList< face > subList
Declare type of subList.
Definition: List.H:144
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write using stream options.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
wordList types() const
Return a list of patch types.
label nFaces() const noexcept
Number of mesh faces.
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:128
void movePoints(const pointField &p)
Correct polyBoundaryMesh after moving points.
scalar range
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:340
label nProcessorPatches() const
The number of processorPolyPatch patches.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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.
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:96
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
A List obtained as a section of another List.
Definition: SubList.H:50
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
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
"scheduled" : (MPI_Send, MPI_Recv)
A class for handling words, derived from Foam::string.
Definition: word.H:63
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:405
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 entries in the list.
Definition: UPtrListI.H:106
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:108
Space [isspace].
Definition: token.H:131
wordList patchNames(nPatches)
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalFaces() const noexcept
Number of internal faces.
End list [isseparator].
Definition: token.H:162
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:108
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:106
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:159
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings...
Definition: wordRe.H:78
const faceList::subList faces() const
Return mesh faces for the entire boundary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
void clear()
Clear the patch list and all demand-driven data.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
label nEdges() const
Number of edges in patch.
void writeEntry(Ostream &os) const
Write as a plain list of entries.
int debug
Static debugging option.
void updateMesh()
Correct polyBoundaryMesh after topology update.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
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:511
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:197
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:472
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
List< word > wordList
List of word.
Definition: fileName.H:59
label nFaces() const noexcept
The number of boundary faces in the underlying mesh.
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:385
#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:281
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:66
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:81
"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:90
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:44
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:74
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
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write using stream options, but always UNCOMPRESSED.
virtual Ostream & endEntry()
Write end entry (&#39;;&#39;) followed by newline.
Definition: Ostream.C:117
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:502
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
List< bool > boolList
A List of bools.
Definition: List.H:60
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
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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
const labelList::subList faceOwner() const
Return face owner for the entire boundary.