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-2024 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::readIOcontents(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  readIOcontents(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 (!readIOcontents(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 (!readIOcontents(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  // Make sure messages don't interact by having unique tag
304  const int oldTag = UPstream::incrMsgType();
305 
306  PstreamBuffers pBufs(Pstream::defaultCommsType);
307 
308  if
309  (
310  pBufs.commsType() == Pstream::commsTypes::buffered
311  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
312  )
313  {
314  forAll(*this, patchi)
315  {
316  operator[](patchi).initGeometry(pBufs);
317  }
318 
319  pBufs.finishedSends();
320 
321  forAll(*this, patchi)
322  {
323  operator[](patchi).calcGeometry(pBufs);
324  }
325  }
326  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
327  {
328  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
329 
330  // Dummy.
331  pBufs.finishedSends();
332 
333  for (const auto& schedEval : patchSchedule)
334  {
335  const label patchi = schedEval.patch;
336 
337  if (schedEval.init)
338  {
339  operator[](patchi).initGeometry(pBufs);
340  }
341  else
342  {
343  operator[](patchi).calcGeometry(pBufs);
344  }
345  }
346  }
347 
348  // Reset tag
349  UPstream::msgType(oldTag);
350 }
351 
352 
354 {
355  return faceList::subList
356  (
357  mesh_.faces(),
358  mesh_.nBoundaryFaces(),
359  mesh_.nInternalFaces()
360  );
361 }
362 
363 
365 {
366  return labelList::subList
367  (
368  mesh_.faceOwner(),
369  mesh_.nBoundaryFaces(),
370  mesh_.nInternalFaces()
371  );
372 }
373 
374 // Potentially useful to simplify logic elsewhere?
375 // const Foam::labelList::subList Foam::polyBoundaryMesh::faceNeighbour() const
376 // {
377 // return labelList::subList();
378 // }
379 
380 
383 {
384  const polyPatchList& patches = *this;
385 
387 
388  forAll(patches, patchi)
389  {
390  list.set(patchi, &patches[patchi].faceCells());
391  }
392 
393  return list;
394 }
395 
396 
399 {
400  if (Pstream::parRun())
401  {
403  << "Neighbour edge addressing not correct across parallel"
404  << " boundaries." << endl;
405  }
406 
407  if (!neighbourEdgesPtr_)
408  {
409  neighbourEdgesPtr_.emplace(size());
410  auto& neighbourEdges = *neighbourEdgesPtr_;
411 
412  // Initialize.
413  label nEdgePairs = 0;
414  forAll(*this, patchi)
415  {
416  const polyPatch& pp = operator[](patchi);
417 
418  neighbourEdges[patchi].setSize(pp.nEdges() - pp.nInternalEdges());
419 
420  for (labelPair& edgeInfo : neighbourEdges[patchi])
421  {
422  edgeInfo[0] = -1;
423  edgeInfo[1] = -1;
424  }
425 
426  nEdgePairs += pp.nEdges() - pp.nInternalEdges();
427  }
428 
429  // From mesh edge (expressed as a point pair so as not to construct
430  // point addressing) to patch + relative edge index.
431  EdgeMap<labelPair> pointsToEdge(nEdgePairs);
432 
433  forAll(*this, patchi)
434  {
435  const polyPatch& pp = operator[](patchi);
436 
437  const edgeList& edges = pp.edges();
438 
439  for
440  (
441  label edgei = pp.nInternalEdges();
442  edgei < edges.size();
443  edgei++
444  )
445  {
446  // Edge in patch local points
447  const edge& e = edges[edgei];
448 
449  // Edge in mesh points.
450  edge meshEdge(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
451 
452  auto fnd = pointsToEdge.find(meshEdge);
453 
454  if (!fnd.good())
455  {
456  // First occurrence of mesh edge. Store patch and my
457  // local index.
458  pointsToEdge.insert
459  (
460  meshEdge,
461  labelPair
462  (
463  patchi,
464  edgei - pp.nInternalEdges()
465  )
466  );
467  }
468  else
469  {
470  // Second occurrence. Store.
471  const labelPair& edgeInfo = fnd.val();
472 
473  neighbourEdges[patchi][edgei - pp.nInternalEdges()] =
474  edgeInfo;
475 
476  neighbourEdges[edgeInfo[0]][edgeInfo[1]]
477  = labelPair(patchi, edgei - pp.nInternalEdges());
478 
479  // Found all two occurrences of this edge so remove from
480  // hash to save space. Note that this will give lots of
481  // problems if the polyBoundaryMesh is multiply connected.
482  pointsToEdge.erase(meshEdge);
483  }
484  }
485  }
486 
487  if (pointsToEdge.size())
488  {
490  << "Not all boundary edges of patches match up." << nl
491  << "Is the outside of your mesh multiply connected?"
492  << abort(FatalError);
493  }
494 
495  forAll(*this, patchi)
496  {
497  const polyPatch& pp = operator[](patchi);
498 
499  const labelPairList& nbrEdges = neighbourEdges[patchi];
500 
501  forAll(nbrEdges, i)
502  {
503  const labelPair& edgeInfo = nbrEdges[i];
504 
505  if (edgeInfo[0] == -1 || edgeInfo[1] == -1)
506  {
507  const label edgei = pp.nInternalEdges() + i;
508  const edge& e = pp.edges()[edgei];
509 
511  << "Not all boundary edges of patches match up." << nl
512  << "Edge " << edgei << " on patch " << pp.name()
513  << " end points " << pp.localPoints()[e[0]] << ' '
514  << pp.localPoints()[e[1]] << " is not matched to an"
515  << " edge on any other patch." << nl
516  << "Is the outside of your mesh multiply connected?"
517  << abort(FatalError);
518  }
519  }
520  }
521  }
522 
523  return *neighbourEdgesPtr_;
524 }
525 
526 
528 {
529  if (!patchIDPtr_)
530  {
531  patchIDPtr_.emplace(mesh_.nBoundaryFaces());
532  auto& list = *patchIDPtr_;
533 
534  const polyPatchList& patches = *this;
535 
536  forAll(patches, patchi)
537  {
538  SubList<label>
539  (
540  list,
541  patches[patchi].size(),
542  (patches[patchi].start() - mesh_.nInternalFaces())
543  ) = patchi;
544  }
545  }
546 
547  return *patchIDPtr_;
548 }
549 
550 
551 Foam::label Foam::polyBoundaryMesh::patchID(const label meshFacei) const
552 {
553  const label bndFacei = (meshFacei - mesh_.nInternalFaces());
554 
555  return
556  (
557  (bndFacei >= 0 && bndFacei < mesh_.nBoundaryFaces())
558  ? this->patchID()[bndFacei]
559  : -1
560  );
561 }
562 
563 
565 Foam::polyBoundaryMesh::patchID(const labelUList& meshFaceIndices) const
566 {
567  labelList output(meshFaceIndices.size());
568  forAll(meshFaceIndices, i)
569  {
570  output[i] = patchID(meshFaceIndices[i]);
571  }
572  return output;
573 }
574 
575 
578 {
579  if (!groupIDsPtr_)
580  {
581  calcGroupIDs();
582  }
583 
584  return *groupIDsPtr_;
585 }
586 
587 
589 (
590  const word& groupName,
591  const labelUList& patchIDs
592 )
593 {
594  groupIDsPtr_.reset(nullptr);
595 
596  polyPatchList& patches = *this;
597 
598  boolList pending(patches.size(), true);
599 
600  // Add to specified patches
601  for (const label patchi : patchIDs)
602  {
603  if (pending.test(patchi))
604  {
605  pending.unset(patchi);
606  patches[patchi].addGroup(groupName);
607  }
608  }
609 
610  // Remove from other patches
611  forAll(patches, patchi)
612  {
613  if (pending.test(patchi))
614  {
615  patches[patchi].removeGroup(groupName);
616  }
617  }
618 }
619 
620 
621 Foam::label Foam::polyBoundaryMesh::nNonProcessor() const
622 {
623  const polyPatchList& patches = *this;
624 
625  label count = 0;
626 
627  for (const polyPatch& p : patches)
628  {
629  if (isA<processorPolyPatch>(p))
630  {
631  break;
632  }
633 
634  ++count;
635  }
636 
637  return count;
638 }
639 
640 
642 {
643  const polyPatchList& patches = *this;
644 
645  label count = 0;
646 
647  for (const polyPatch& p : patches)
648  {
649  if (isA<processorPolyPatch>(p))
650  {
651  ++count;
652  }
653  }
654 
655  return count;
656 }
657 
660 {
661  return PtrListOps::get<word>(*this, nameOp<polyPatch>());
662 }
663 
666 {
667  return PtrListOps::get<word>(*this, typeOp<polyPatch>());
668 }
669 
670 
672 {
673  return
674  PtrListOps::get<word>
675  (
676  *this,
677  [](const polyPatch& p) { return p.physicalType(); }
678  );
679 }
680 
681 
683 {
684  return
685  PtrListOps::get<label>
686  (
687  *this,
688  [](const polyPatch& p) { return p.start(); }
689  );
690 }
691 
692 
694 {
695  return
696  PtrListOps::get<label>
697  (
698  *this,
699  [](const polyPatch& p) { return p.size(); }
700  );
701 }
702 
703 
705 {
706  return
707  PtrListOps::get<labelRange>
708  (
709  *this,
710  [](const polyPatch& p) { return p.range(); }
711  );
712 }
713 
716 {
717  return this->groupPatchIDs().sortedToc();
718 }
719 
721 Foam::label Foam::polyBoundaryMesh::start() const noexcept
722 {
723  return mesh_.nInternalFaces();
724 }
725 
727 Foam::label Foam::polyBoundaryMesh::nFaces() const noexcept
728 {
729  return mesh_.nBoundaryFaces();
730 }
731 
734 {
735  return labelRange(mesh_.nInternalFaces(), mesh_.nBoundaryFaces());
736 }
737 
738 
739 Foam::labelRange Foam::polyBoundaryMesh::range(const label patchi) const
740 {
741  if (patchi < 0)
742  {
743  return labelRange(mesh_.nInternalFaces(), 0);
744  }
746  // Will fail if patchi >= size()
747  return (*this)[patchi].range();
748 }
749 
750 
752 (
753  const wordRe& matcher,
754  const bool useGroups
755 ) const
756 {
757  if (matcher.empty())
758  {
759  return labelList();
760  }
761 
762  // Only check groups if requested and they exist
763  const bool checkGroups = (useGroups && this->hasGroupIDs());
764 
765  labelHashSet ids(0);
766 
767  if (matcher.isPattern())
768  {
769  if (checkGroups)
770  {
771  const auto& groupLookup = groupPatchIDs();
772  forAllConstIters(groupLookup, iter)
773  {
774  if (matcher(iter.key()))
775  {
776  // Add patch ids associated with the group
777  ids.insert(iter.val());
778  }
779  }
780  }
781 
782  if (ids.empty())
783  {
784  return PtrListOps::findMatching(*this, matcher);
785  }
786  else
787  {
788  ids.insert(PtrListOps::findMatching(*this, matcher));
789  }
790  }
791  else
792  {
793  // Literal string.
794  // Special version of above for reduced memory footprint.
795 
796  const label patchId = PtrListOps::firstMatching(*this, matcher);
797 
798  if (patchId >= 0)
799  {
800  return labelList(one{}, patchId);
801  }
802  else if (checkGroups)
803  {
804  const auto iter = groupPatchIDs().cfind(matcher);
805 
806  if (iter.good())
807  {
808  // Hash ids associated with the group
809  ids.insert(iter.val());
810  }
811  }
812  }
813 
814  return ids.sortedToc();
815 }
816 
817 
819 (
820  const wordRes& matcher,
821  const bool useGroups
822 ) const
823 {
824  if (matcher.empty())
825  {
826  return labelList();
827  }
828  else if (matcher.size() == 1)
829  {
830  return this->indices(matcher.front(), useGroups);
831  }
832 
833  labelHashSet ids(0);
834 
835  // Only check groups if requested and they exist
836  if (useGroups && this->hasGroupIDs())
837  {
838  ids.reserve(this->size());
839 
840  const auto& groupLookup = groupPatchIDs();
841  forAllConstIters(groupLookup, iter)
842  {
843  if (matcher(iter.key()))
844  {
845  // Add patch ids associated with the group
846  ids.insert(iter.val());
847  }
848  }
849  }
850 
851  if (ids.empty())
852  {
853  return PtrListOps::findMatching(*this, matcher);
854  }
855  else
856  {
857  ids.insert(PtrListOps::findMatching(*this, matcher));
858  }
859 
860  return ids.sortedToc();
861 }
862 
863 
865 (
866  const wordRes& select,
867  const wordRes& ignore,
868  const bool useGroups
869 ) const
870 {
871  if (ignore.empty())
872  {
873  return this->indices(select, useGroups);
874  }
875 
876  const wordRes::filter matcher(select, ignore);
877 
878  labelHashSet ids(0);
879 
880  // Only check groups if requested and they exist
881  if (useGroups && this->hasGroupIDs())
882  {
883  ids.reserve(this->size());
884 
885  const auto& groupLookup = groupPatchIDs();
886  forAllConstIters(groupLookup, iter)
887  {
888  if (matcher(iter.key()))
889  {
890  // Add patch ids associated with the group
891  ids.insert(iter.val());
892  }
893  }
894  }
895 
896  if (ids.empty())
897  {
898  return PtrListOps::findMatching(*this, matcher);
899  }
900  else
901  {
902  ids.insert(PtrListOps::findMatching(*this, matcher));
903  }
904 
905  return ids.sortedToc();
906 }
907 
908 
909 Foam::label Foam::polyBoundaryMesh::findIndex(const wordRe& key) const
910 {
911  if (key.empty())
912  {
913  return -1;
914  }
915  return PtrListOps::firstMatching(*this, key);
916 }
917 
918 
920 (
921  const word& patchName,
922  bool allowNotFound
923 ) const
924 {
925  if (patchName.empty())
926  {
927  return -1;
928  }
929 
930  const label patchId = PtrListOps::firstMatching(*this, patchName);
931 
932  if (patchId >= 0)
933  {
934  return patchId;
935  }
936 
937  if (!allowNotFound)
938  {
940  << "Patch '" << patchName << "' not found. "
941  << "Available patch names";
942 
943  if (polyMesh::defaultRegion != mesh_.name())
944  {
945  FatalError
946  << " in region '" << mesh_.name() << "'";
947  }
948 
949  FatalError
950  << " include: " << names() << endl
951  << exit(FatalError);
952  }
953 
954  // Patch not found
955  if (debug)
956  {
957  Pout<< "label polyBoundaryMesh::findPatchID(const word&) const"
958  << "Patch named " << patchName << " not found. "
959  << "Available patch names: " << names() << endl;
960  }
962  // Not found, return -1
963  return -1;
964 }
965 
966 
968 Foam::polyBoundaryMesh::whichPatchFace(const label meshFacei) const
969 {
970  if (meshFacei < mesh().nInternalFaces())
971  {
972  // Internal face: return (-1, meshFace)
973  return labelPair(-1, meshFacei);
974  }
975  else if (meshFacei >= mesh().nFaces())
976  {
977  // Bounds error: abort
979  << "Face " << meshFacei
980  << " out of bounds. Number of geometric faces " << mesh().nFaces()
981  << abort(FatalError);
982 
983  return labelPair(-1, meshFacei);
984  }
985 
986 
987  const polyPatchList& patches = *this;
988 
989  // Do we have cached patch info?
990  if (patchIDPtr_)
991  {
992  const label patchi =
993  this->patchID()[meshFacei - mesh().nInternalFaces()];
994 
995  // (patch, local face index)
996  return labelPair(patchi, meshFacei - patches[patchi].start());
997  }
998 
999 
1000  // Patches are ordered, use binary search
1001  // Find out which patch index and local patch face the specified
1002  // mesh face belongs to by comparing label with patch start labels.
1003 
1004  const label patchi =
1005  findLower
1006  (
1007  patches,
1008  meshFacei,
1009  0,
1010  // Must include the start in the comparison
1011  [](const polyPatch& p, label val) { return (p.start() <= val); }
1012  );
1013 
1014  if (patchi < 0 || !patches[patchi].range().contains(meshFacei))
1015  {
1016  // If not in any of above, it is trouble!
1018  << "Face " << meshFacei << " not found in any of the patches "
1019  << flatOutput(names()) << nl
1020  << "The patches appear to be inconsistent with the mesh :"
1021  << " internalFaces:" << mesh().nInternalFaces()
1022  << " total number of faces:" << mesh().nFaces()
1023  << abort(FatalError);
1024 
1025  return labelPair(-1, meshFacei);
1026  }
1028  // (patch, local face index)
1029  return labelPair(patchi, meshFacei - patches[patchi].start());
1030 }
1031 
1032 
1034 Foam::polyBoundaryMesh::whichPatchFace(const labelUList& meshFaceIndices) const
1035 {
1036  labelPairList output(meshFaceIndices.size());
1037  forAll(meshFaceIndices, i)
1038  {
1039  output[i] = whichPatchFace(meshFaceIndices[i]);
1040  }
1041  return output;
1042 }
1043 
1044 
1046 (
1047  const UList<wordRe>& select,
1048  const bool warnNotFound,
1049  const bool useGroups
1050 ) const
1051 {
1052  labelHashSet ids(0);
1053  if (select.empty())
1054  {
1055  return ids;
1056  }
1057 
1058  const polyPatchList& patches = *this;
1059 
1060  const label len = patches.size();
1061 
1062  ids.reserve(len);
1063 
1064  // Only check groups if requested and they exist
1065  const bool checkGroups = (useGroups && this->hasGroupIDs());
1066 
1067  for (const wordRe& matcher : select)
1068  {
1069  bool missed = true;
1070 
1071  for (label i = 0; i < len; ++i)
1072  {
1073  if (matcher(patches[i].name()))
1074  {
1075  ids.insert(i);
1076  missed = false;
1077  }
1078  }
1079 
1080  if (missed && checkGroups)
1081  {
1082  // Check group names
1083  if (matcher.isPattern())
1084  {
1085  forAllConstIters(groupPatchIDs(), iter)
1086  {
1087  if (matcher.match(iter.key()))
1088  {
1089  // Hash ids associated with the group
1090  ids.insert(iter.val());
1091  missed = false;
1092  }
1093  }
1094  }
1095  else
1096  {
1097  const auto iter = groupPatchIDs().cfind(matcher);
1098 
1099  if (iter.good())
1100  {
1101  // Hash ids associated with the group
1102  ids.insert(iter.val());
1103  missed = false;
1104  }
1105  }
1106  }
1107 
1108  if (missed && warnNotFound)
1109  {
1110  if (checkGroups)
1111  {
1113  << "Cannot find any patch or group names matching "
1114  << matcher << endl;
1115  }
1116  else
1117  {
1119  << "Cannot find any patch names matching "
1120  << matcher << endl;
1121  }
1122  }
1123  }
1124 
1125  return ids;
1126 }
1127 
1128 
1130 (
1131  const labelUList& patchIDs,
1132  wordList& groups,
1133  labelHashSet& nonGroupPatches
1134 ) const
1135 {
1136  // Current matched groups
1137  DynamicList<word> matchedGroups(1);
1138 
1139  // Current set of unmatched patches
1140  nonGroupPatches = labelHashSet(patchIDs);
1141 
1142  const HashTable<labelList>& groupLookup = this->groupPatchIDs();
1143  forAllConstIters(groupLookup, iter)
1144  {
1145  // Store currently unmatched patches so we can restore
1146  labelHashSet oldNonGroupPatches(nonGroupPatches);
1147 
1148  // Match by deleting patches in group from the current set and seeing
1149  // if all have been deleted.
1150  labelHashSet groupPatchSet(iter.val());
1151 
1152  label nMatch = nonGroupPatches.erase(groupPatchSet);
1153 
1154  if (nMatch == groupPatchSet.size())
1155  {
1156  matchedGroups.push_back(iter.key());
1157  }
1158  else if (nMatch != 0)
1159  {
1160  // No full match. Undo.
1161  nonGroupPatches.transfer(oldNonGroupPatches);
1162  }
1163  }
1164 
1165  groups.transfer(matchedGroups);
1166 }
1167 
1168 
1169 bool Foam::polyBoundaryMesh::checkParallelSync(const bool report) const
1170 {
1171  if (!Pstream::parRun())
1172  {
1173  return false;
1174  }
1175 
1176  const polyBoundaryMesh& bm = *this;
1177 
1178  bool hasError = false;
1179 
1180  // Collect non-proc patches and check proc patches are last.
1181  wordList localNames(bm.size());
1182  wordList localTypes(bm.size());
1183 
1184  label nonProci = 0;
1185 
1186  forAll(bm, patchi)
1187  {
1188  if (!isA<processorPolyPatch>(bm[patchi]))
1189  {
1190  if (nonProci != patchi)
1191  {
1192  // A processor patch in between normal patches!
1193  hasError = true;
1194 
1195  if (debug || report)
1196  {
1197  Pout<< " ***Problem with boundary patch " << patchi
1198  << " name:" << bm[patchi].name()
1199  << " type:" << bm[patchi].type()
1200  << " - seems to be preceeded by processor patches."
1201  << " This is usually a problem." << endl;
1202  }
1203  }
1204  else
1205  {
1206  localNames[nonProci] = bm[patchi].name();
1207  localTypes[nonProci] = bm[patchi].type();
1208  ++nonProci;
1209  }
1210  }
1211  }
1212  localNames.resize(nonProci);
1213  localTypes.resize(nonProci);
1214 
1215  // Check and report error(s) on master
1216  // - don't need indexing on master itself
1217 
1218  const globalIndex procAddr(globalIndex::gatherNonLocal{}, nonProci);
1219 
1220  const wordList allNames(procAddr.gather(localNames));
1221  const wordList allTypes(procAddr.gather(localTypes));
1222 
1223  // Automatically restricted to master
1224  for (const int proci : procAddr.subProcs())
1225  {
1226  const auto procNames(allNames.slice(procAddr.range(proci)));
1227  const auto procTypes(allTypes.slice(procAddr.range(proci)));
1228 
1229  if (procNames != localNames || procTypes != localTypes)
1230  {
1231  hasError = true;
1232 
1233  if (debug || report)
1234  {
1235  Info<< " ***Inconsistent patches across processors, "
1236  "processor0 has patch names:" << localNames
1237  << " patch types:" << localTypes
1238  << " processor" << proci
1239  << " has patch names:" << procNames
1240  << " patch types:" << procTypes
1241  << endl;
1242  }
1243  }
1244  }
1245 
1246  // Reduce (not broadcast) to respect local out-of-order errors (first loop)
1247  return returnReduceOr(hasError);
1248 }
1249 
1250 
1251 bool Foam::polyBoundaryMesh::checkDefinition(const bool report) const
1252 {
1253  label nextPatchStart = mesh().nInternalFaces();
1254  const polyBoundaryMesh& bm = *this;
1255 
1256  bool hasError = false;
1257 
1258  wordHashSet patchNames(2*this->size());
1259 
1260  forAll(bm, patchi)
1261  {
1262  if (bm[patchi].start() != nextPatchStart && !hasError)
1263  {
1264  hasError = true;
1265 
1266  Info<< " ****Problem with boundary patch " << patchi
1267  << " named " << bm[patchi].name()
1268  << " of type " << bm[patchi].type()
1269  << ". The patch should start on face no " << nextPatchStart
1270  << " and the patch specifies " << bm[patchi].start()
1271  << "." << endl
1272  << "Possibly consecutive patches have this same problem."
1273  << " Suppressing future warnings." << endl;
1274  }
1275 
1276  if (!patchNames.insert(bm[patchi].name()) && !hasError)
1277  {
1278  hasError = true;
1279 
1280  Info<< " ****Duplicate boundary patch " << patchi
1281  << " named " << bm[patchi].name()
1282  << " of type " << bm[patchi].type()
1283  << "." << endl
1284  << "Suppressing future warnings." << endl;
1285  }
1286 
1287  nextPatchStart += bm[patchi].size();
1288  }
1289 
1290  Pstream::reduceOr(hasError);
1291 
1292  if (debug || report)
1293  {
1294  if (hasError)
1295  {
1296  Pout<< " ***Boundary definition is in error." << endl;
1297  }
1298  else
1299  {
1300  Info<< " Boundary definition OK." << endl;
1301  }
1302  }
1303 
1304  return hasError;
1305 }
1306 
1307 
1309 {
1310  PstreamBuffers pBufs(Pstream::defaultCommsType);
1311 
1312  if
1313  (
1314  pBufs.commsType() == Pstream::commsTypes::buffered
1315  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
1316  )
1317  {
1318  forAll(*this, patchi)
1319  {
1320  operator[](patchi).initMovePoints(pBufs, p);
1321  }
1322 
1323  pBufs.finishedSends();
1324 
1325  forAll(*this, patchi)
1326  {
1327  operator[](patchi).movePoints(pBufs, p);
1328  }
1329  }
1330  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1331  {
1332  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1333 
1334  // Dummy.
1335  pBufs.finishedSends();
1336 
1337  for (const auto& schedEval : patchSchedule)
1338  {
1339  const label patchi = schedEval.patch;
1340 
1341  if (schedEval.init)
1342  {
1343  operator[](patchi).initMovePoints(pBufs, p);
1344  }
1345  else
1346  {
1347  operator[](patchi).movePoints(pBufs, p);
1348  }
1349  }
1350  }
1351 }
1352 
1353 
1355 {
1356  neighbourEdgesPtr_.reset(nullptr);
1357  patchIDPtr_.reset(nullptr);
1358  groupIDsPtr_.reset(nullptr);
1359 
1360  PstreamBuffers pBufs(Pstream::defaultCommsType);
1361 
1362  if
1363  (
1364  pBufs.commsType() == Pstream::commsTypes::buffered
1365  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
1366  )
1367  {
1368  forAll(*this, patchi)
1369  {
1370  operator[](patchi).initUpdateMesh(pBufs);
1371  }
1372 
1373  pBufs.finishedSends();
1374 
1375  forAll(*this, patchi)
1376  {
1377  operator[](patchi).updateMesh(pBufs);
1378  }
1379  }
1380  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
1381  {
1382  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1383 
1384  // Dummy.
1385  pBufs.finishedSends();
1386 
1387  for (const auto& schedEval : patchSchedule)
1388  {
1389  const label patchi = schedEval.patch;
1390 
1391  if (schedEval.init)
1392  {
1393  operator[](patchi).initUpdateMesh(pBufs);
1394  }
1395  else
1396  {
1397  operator[](patchi).updateMesh(pBufs);
1398  }
1399  }
1400  }
1401 }
1402 
1403 
1405 (
1406  const labelUList& oldToNew,
1407  const bool validBoundary
1408 )
1409 {
1410  // Change order of patches
1411  polyPatchList::reorder(oldToNew);
1412 
1413  // Adapt indices
1414  polyPatchList& patches = *this;
1415 
1416  forAll(patches, patchi)
1417  {
1418  patches[patchi].index() = patchi;
1419  }
1420 
1421  // Clear group-to-patch addressing. Note: could re-calculate
1422  groupIDsPtr_.reset(nullptr);
1423 
1424  if (validBoundary)
1425  {
1426  updateMesh();
1427  }
1428 }
1429 
1430 
1431 void Foam::polyBoundaryMesh::writeEntry(Ostream& os) const
1432 {
1433  const polyPatchList& entries = *this;
1434 
1435  os << entries.size();
1436 
1437  if (entries.empty())
1438  {
1439  // 0-sized : can write with less vertical space
1441  }
1442  else
1443  {
1444  os << nl << token::BEGIN_LIST << incrIndent << nl;
1445 
1446  for (const auto& pp : entries)
1447  {
1448  os.beginBlock(pp.name());
1449  os << pp;
1450  os.endBlock();
1451  }
1453  }
1455 }
1456 
1457 
1459 (
1460  const word& keyword,
1461  Ostream& os
1462 ) const
1463 {
1464  const polyPatchList& entries = *this;
1465 
1466  if (!keyword.empty())
1467  {
1468  os.write(keyword);
1469  os << (entries.empty() ? token::SPACE : token::NL);
1470  }
1472  writeEntry(os);
1473 
1474  if (!keyword.empty()) os.endEntry();
1475 }
1476 
1477 
1480  writeEntry(os);
1481  return os.good();
1482 }
1483 
1484 
1486 (
1487  IOstreamOption streamOpt,
1488  const bool writeOnProc
1489 ) const
1490 {
1492  return regIOobject::writeObject(streamOpt, writeOnProc);
1493 }
1494 
1495 
1496 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
1497 
1499 (
1500  const word& patchName
1501 ) const
1502 {
1503  const label patchi = findPatchID(patchName);
1504 
1505  if (patchi < 0)
1506  {
1508  << "Patch named " << patchName << " not found." << nl
1509  << "Available patch names: " << names() << endl
1510  << abort(FatalError);
1511  }
1512 
1513  return operator[](patchi);
1514 }
1515 
1516 
1518 (
1519  const word& patchName
1520 )
1521 {
1522  const label patchi = findPatchID(patchName);
1523 
1524  if (patchi < 0)
1525  {
1527  << "Patch named " << patchName << " not found." << nl
1528  << "Available patch names: " << names() << endl
1529  << abort(FatalError);
1530  }
1532  return operator[](patchi);
1533 }
1534 
1535 
1536 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1537 
1538 Foam::Ostream& Foam::operator<<(Ostream& os, const polyBoundaryMesh& pbm)
1539 {
1540  pbm.writeData(os);
1541  return os;
1542 }
1543 
1544 
1545 // ************************************************************************* //
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:629
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:805
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...
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
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1274
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:608
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:675
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:1061
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
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1252
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:134
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 standard) : (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:406
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:496
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:397
#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:68
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:81
"nonBlocking" (immediate) : (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:75
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
"buffered" : (MPI_Bsend, MPI_Recv)
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:180
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.