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-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "polyBoundaryMesh.H"
30 #include "polyMesh.H"
31 #include "primitiveMesh.H"
32 #include "processorPolyPatch.H"
33 #include "PstreamBuffers.H"
34 #include "lduSchedule.H"
35 #include "globalMeshData.H"
36 #include "stringListOps.H"
37 #include "edgeHashes.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(polyBoundaryMesh, 0);
44 }
45 
46 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50  // Templated implementation for types(), names(), etc - file-scope
51  template<class ListType, class GetOp>
52  static ListType getMethodImpl
53  (
54  const polyPatchList& list,
55  const GetOp& getop
56  )
57  {
58  const label len = list.size();
59 
60  ListType output(len);
61 
62  for (label i = 0; i < len; ++i)
63  {
64  output[i] = getop(list[i]);
65  }
66 
67  return output;
68  }
69 
70 
71  // Templated implementation for indices() - file-scope
72  template<class UnaryMatchPredicate>
73  static labelList indicesImpl
74  (
75  const polyPatchList& list,
76  const UnaryMatchPredicate& matcher
77  )
78  {
79  const label len = list.size();
80 
81  labelList output(len);
82 
83  label count = 0;
84  for (label i = 0; i < len; ++i)
85  {
86  if (matcher(list[i].name()))
87  {
88  output[count++] = i;
89  }
90  }
91 
92  output.resize(count);
93 
94  return output;
95  }
96 
97 
98  // Templated implementation for findIndex() - file-scope
99  template<class UnaryMatchPredicate>
100  label findIndexImpl
101  (
102  const polyPatchList& list,
103  const UnaryMatchPredicate& matcher
104  )
105  {
106  const label len = list.size();
107 
108  for (label i = 0; i < len; ++i)
109  {
110  if (matcher(list[i].name()))
111  {
112  return i;
113  }
114  }
115 
116  return -1;
117  }
118 
119 } // End namespace Foam
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
124 Foam::polyBoundaryMesh::polyBoundaryMesh
125 (
126  const IOobject& io,
127  const polyMesh& mesh
128 )
129 :
130  polyPatchList(),
131  regIOobject(io),
132  mesh_(mesh)
133 {
134  if
135  (
136  readOpt() == IOobject::MUST_READ
137  || readOpt() == IOobject::MUST_READ_IF_MODIFIED
138  )
139  {
140  // Warn for MUST_READ_IF_MODIFIED
141  warnNoRereading<polyBoundaryMesh>();
142 
143  polyPatchList& patches = *this;
144 
145  // Read polyPatchList
146  Istream& is = readStream(typeName);
147 
148  PtrList<entry> patchEntries(is);
149  patches.setSize(patchEntries.size());
150 
151  forAll(patches, patchi)
152  {
153  patches.set
154  (
155  patchi,
157  (
158  patchEntries[patchi].keyword(),
159  patchEntries[patchi].dict(),
160  patchi,
161  *this
162  )
163  );
164  }
165 
166  is.check(FUNCTION_NAME);
167 
168  close();
169  }
170 }
171 
172 
173 Foam::polyBoundaryMesh::polyBoundaryMesh
174 (
175  const IOobject& io,
176  const polyMesh& pm,
177  const label size
178 )
179 :
180  polyPatchList(size),
181  regIOobject(io),
182  mesh_(pm)
183 {}
184 
185 
186 Foam::polyBoundaryMesh::polyBoundaryMesh
187 (
188  const IOobject& io,
189  const polyMesh& pm,
190  const polyPatchList& ppl
191 )
192 :
193  polyPatchList(),
194  regIOobject(io),
195  mesh_(pm)
196 {
197  if
198  (
199  (this->readOpt() == IOobject::READ_IF_PRESENT && this->headerOk())
200  || this->readOpt() == IOobject::MUST_READ
201  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
202  )
203  {
204  // Warn for MUST_READ_IF_MODIFIED
205  warnNoRereading<polyBoundaryMesh>();
206 
207  polyPatchList& patches = *this;
208 
209  // Read polyPatchList
210  Istream& is = readStream(typeName);
211 
212  PtrList<entry> patchEntries(is);
213  patches.setSize(patchEntries.size());
214 
215  forAll(patches, patchi)
216  {
217  patches.set
218  (
219  patchi,
221  (
222  patchEntries[patchi].keyword(),
223  patchEntries[patchi].dict(),
224  patchi,
225  *this
226  )
227  );
228  }
229 
230  is.check(FUNCTION_NAME);
231 
232  close();
233  }
234  else
235  {
236  polyPatchList& patches = *this;
237  patches.setSize(ppl.size());
238  forAll(patches, patchi)
239  {
240  patches.set(patchi, ppl[patchi].clone(*this).ptr());
241  }
242  }
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
247 
249 {
250  polyPatchList& patches = *this;
251 
252  for (polyPatch& p : patches)
253  {
254  p.clearGeom();
255  }
256 }
257 
258 
260 {
261  neighbourEdgesPtr_.clear();
262  patchIDPtr_.clear();
263  groupPatchIDsPtr_.clear();
264 
265  polyPatchList& patches = *this;
266 
267  for (polyPatch& p : patches)
268  {
269  p.clearAddressing();
270  }
271 }
272 
273 
274 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
275 
276 void Foam::polyBoundaryMesh::calcGeometry()
277 {
279 
280  if
281  (
284  )
285  {
286  forAll(*this, patchi)
287  {
288  operator[](patchi).initGeometry(pBufs);
289  }
290 
291  pBufs.finishedSends();
292 
293  forAll(*this, patchi)
294  {
295  operator[](patchi).calcGeometry(pBufs);
296  }
297  }
299  {
300  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
301 
302  // Dummy.
303  pBufs.finishedSends();
304 
305  for (const auto& patchEval : patchSchedule)
306  {
307  const label patchi = patchEval.patch;
308 
309  if (patchEval.init)
310  {
311  operator[](patchi).initGeometry(pBufs);
312  }
313  else
314  {
315  operator[](patchi).calcGeometry(pBufs);
316  }
317  }
318  }
319 }
320 
321 
324 {
325  if (Pstream::parRun())
326  {
328  << "Neighbour edge addressing not correct across parallel"
329  << " boundaries." << endl;
330  }
331 
332  if (!neighbourEdgesPtr_)
333  {
334  neighbourEdgesPtr_.reset(new List<labelPairList>(size()));
335  List<labelPairList>& neighbourEdges = neighbourEdgesPtr_();
336 
337  // Initialize.
338  label nEdgePairs = 0;
339  forAll(*this, patchi)
340  {
341  const polyPatch& pp = operator[](patchi);
342 
343  neighbourEdges[patchi].setSize(pp.nEdges() - pp.nInternalEdges());
344 
345  for (labelPair& edgeInfo : neighbourEdges[patchi])
346  {
347  edgeInfo[0] = -1;
348  edgeInfo[1] = -1;
349  }
350 
351  nEdgePairs += pp.nEdges() - pp.nInternalEdges();
352  }
353 
354  // From mesh edge (expressed as a point pair so as not to construct
355  // point addressing) to patch + relative edge index.
356  EdgeMap<labelPair> pointsToEdge(nEdgePairs);
357 
358  forAll(*this, patchi)
359  {
360  const polyPatch& pp = operator[](patchi);
361 
362  const edgeList& edges = pp.edges();
363 
364  for
365  (
366  label edgei = pp.nInternalEdges();
367  edgei < edges.size();
368  edgei++
369  )
370  {
371  // Edge in patch local points
372  const edge& e = edges[edgei];
373 
374  // Edge in mesh points.
375  edge meshEdge(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
376 
377  auto fnd = pointsToEdge.find(meshEdge);
378 
379  if (!fnd.found())
380  {
381  // First occurrence of mesh edge. Store patch and my
382  // local index.
383  pointsToEdge.insert
384  (
385  meshEdge,
386  labelPair
387  (
388  patchi,
389  edgei - pp.nInternalEdges()
390  )
391  );
392  }
393  else
394  {
395  // Second occurrence. Store.
396  const labelPair& edgeInfo = fnd.val();
397 
398  neighbourEdges[patchi][edgei - pp.nInternalEdges()] =
399  edgeInfo;
400 
401  neighbourEdges[edgeInfo[0]][edgeInfo[1]]
402  = labelPair(patchi, edgei - pp.nInternalEdges());
403 
404  // Found all two occurrences of this edge so remove from
405  // hash to save space. Note that this will give lots of
406  // problems if the polyBoundaryMesh is multiply connected.
407  pointsToEdge.erase(meshEdge);
408  }
409  }
410  }
411 
412  if (pointsToEdge.size())
413  {
415  << "Not all boundary edges of patches match up." << nl
416  << "Is the outside of your mesh multiply connected?"
417  << abort(FatalError);
418  }
419 
420  forAll(*this, patchi)
421  {
422  const polyPatch& pp = operator[](patchi);
423 
424  const labelPairList& nbrEdges = neighbourEdges[patchi];
425 
426  forAll(nbrEdges, i)
427  {
428  const labelPair& edgeInfo = nbrEdges[i];
429 
430  if (edgeInfo[0] == -1 || edgeInfo[1] == -1)
431  {
432  const label edgei = pp.nInternalEdges() + i;
433  const edge& e = pp.edges()[edgei];
434 
436  << "Not all boundary edges of patches match up." << nl
437  << "Edge " << edgei << " on patch " << pp.name()
438  << " end points " << pp.localPoints()[e[0]] << ' '
439  << pp.localPoints()[e[1]] << " is not matched to an"
440  << " edge on any other patch." << nl
441  << "Is the outside of your mesh multiply connected?"
442  << abort(FatalError);
443  }
444  }
445  }
446  }
447 
448  return *neighbourEdgesPtr_;
449 }
450 
451 
453 {
454  if (!patchIDPtr_)
455  {
456  patchIDPtr_.reset(new labelList(mesh_.nBoundaryFaces()));
457  labelList& list = *patchIDPtr_;
458 
459  const polyPatchList& patches = *this;
460 
461  forAll(patches, patchi)
462  {
464  (
465  list,
466  patches[patchi].size(),
467  (patches[patchi].start() - mesh_.nInternalFaces())
468  ) = patchi;
469  }
470  }
471 
472  return *patchIDPtr_;
473 }
474 
475 
478 {
479  if (!groupPatchIDsPtr_)
480  {
481  groupPatchIDsPtr_.reset(new HashTable<labelList>(16));
482  auto& groupPatchIDs = *groupPatchIDsPtr_;
483 
484  const polyBoundaryMesh& patches = *this;
485 
486  forAll(patches, patchi)
487  {
488  const wordList& groups = patches[patchi].inGroups();
489 
490  for (const word& groupName : groups)
491  {
492  auto iter = groupPatchIDs.find(groupName);
493 
494  if (iter.found())
495  {
496  (*iter).append(patchi);
497  }
498  else
499  {
500  groupPatchIDs.insert(groupName, labelList(one{}, patchi));
501  }
502  }
503  }
504 
505  // Remove patch names from patchGroups
506  forAll(patches, patchi)
507  {
508  if (groupPatchIDs.erase(patches[patchi].name()))
509  {
511  << "Removed patchGroup '" << patches[patchi].name()
512  << "' which clashes with patch " << patchi
513  << " of the same name."
514  << endl;
515  }
516  }
517  }
518 
519  return *groupPatchIDsPtr_;
520 }
521 
522 
524 (
525  const word& groupName,
526  const labelUList& patchIDs
527 )
528 {
529  groupPatchIDsPtr_.clear();
530 
531  polyPatchList& patches = *this;
532 
533  boolList donePatch(patches.size(), false);
534 
535  // Add to specified patches
536  for (const label patchi : patchIDs)
537  {
538  polyPatch& pp = patches[patchi];
539 
540  if (!pp.inGroup(groupName))
541  {
542  pp.inGroups().append(groupName);
543  }
544  donePatch[patchi] = true;
545  }
546 
547  // Remove from other patches
548  forAll(patches, patchi)
549  {
550  if (!donePatch[patchi])
551  {
552  polyPatch& pp = patches[patchi];
553 
554  label newI = 0;
555  if (pp.inGroup(groupName))
556  {
557  wordList& groups = pp.inGroups();
558 
559  forAll(groups, i)
560  {
561  if (groups[i] != groupName)
562  {
563  groups[newI++] = groups[i];
564  }
565  }
566  groups.setSize(newI);
567  }
568  }
569  }
570 }
571 
572 
574 {
575  const polyPatchList& patches = *this;
576 
577  label nonProc = 0;
578 
579  for (const polyPatch& p : patches)
580  {
581  if (isA<processorPolyPatch>(p))
582  {
583  break;
584  }
585 
586  ++nonProc;
587  }
588 
589  return nonProc;
590 }
591 
592 
594 {
595  return getMethodImpl<wordList>(*this, getNameOp<polyPatch>());
596 }
597 
598 
600 {
601  return getMethodImpl<wordList>(*this, getTypeOp<polyPatch>());
602 }
603 
604 
606 {
607  return
608  getMethodImpl<wordList>
609  (
610  *this,
611  [](const polyPatch& p) { return p.physicalType(); }
612  );
613 }
614 
615 
617 {
618  return
619  getMethodImpl<labelList>
620  (
621  *this,
622  [](const polyPatch& p) { return p.start(); }
623  );
624 }
625 
626 
628 {
629  return
630  getMethodImpl<labelList>
631  (
632  *this,
633  [](const polyPatch& p) { return p.size(); }
634  );
635 }
636 
637 
638 Foam::label Foam::polyBoundaryMesh::start() const
639 {
640  return mesh_.nInternalFaces();
641 }
642 
643 
645 {
646  return mesh_.nBoundaryFaces();
647 }
648 
649 
651 {
652  return labelRange(mesh_.nInternalFaces(), mesh_.nBoundaryFaces());
653 }
654 
655 
657 {
658  if (patchi < 0)
659  {
660  return labelRange(mesh_.nInternalFaces(), 0);
661  }
662 
663  // Will fail if patchi >= size()
664  return (*this)[patchi].range();
665 }
666 
667 
669 (
670  const keyType& key,
671  const bool useGroups
672 ) const
673 {
674  if (key.empty())
675  {
676  return labelList();
677  }
678 
679  DynamicList<label> patchIndices;
680 
681  if (key.isPattern())
682  {
683  const regExp matcher(key);
684 
685  patchIndices = indicesImpl(*this, matcher);
686 
687  // Only examine patch groups if requested and when they exist.
688  if (useGroups && !groupPatchIDs().empty())
689  {
690  const wordList groupNames
691  (
692  groupPatchIDs().tocKeys(matcher)
693  );
694 
695  if (groupNames.size())
696  {
697  labelHashSet groupIndices;
698 
699  for (const word& grpName : groupNames)
700  {
701  // Hash the patch ids for the group
702  groupIndices.insert( groupPatchIDs()[grpName] );
703  }
704 
705  groupIndices.erase(patchIndices); // Skip existing
706  patchIndices.append(groupIndices.sortedToc());
707  }
708  }
709  }
710  else
711  {
712  // Literal string.
713  // Special version of above for reduced memory footprint
714 
715  const word& matcher = key;
716 
717  const label patchId = findIndexImpl(*this, matcher);
718 
719  if (patchId >= 0)
720  {
721  patchIndices.append(patchId);
722  }
723 
724  // Only examine patch groups if requested and when they exist.
725  if (useGroups && !groupPatchIDs().empty())
726  {
727  const auto iter = groupPatchIDs().cfind(key);
728 
729  if (iter.found())
730  {
731  // Hash the patch ids for the group
732  labelHashSet groupIndices(iter.val());
733 
734  groupIndices.erase(patchIndices); // Skip existing
735  patchIndices.append(groupIndices.sortedToc());
736  }
737  }
738  }
739 
740  return patchIndices;
741 }
742 
743 
744 Foam::label Foam::polyBoundaryMesh::findIndex(const keyType& key) const
745 {
746  if (key.empty())
747  {
748  return -1;
749  }
750  else if (key.isPattern())
751  {
752  // Find as regex
753  regExp matcher(key);
754  return findIndexImpl(*this, matcher);
755  }
756  else
757  {
758  // Find as literal string
759  const word& matcher = key;
760  return findIndexImpl(*this, matcher);
761  }
762 }
763 
764 
766 (
767  const word& patchName,
768  bool allowNotFound
769 ) const
770 {
771  if (patchName.empty())
772  {
773  return -1;
774  }
775 
776  const label patchId = findIndexImpl(*this, patchName);
777 
778  if (patchId >= 0)
779  {
780  return patchId;
781  }
782 
783  if (!allowNotFound)
784  {
785  string regionStr;
786  if (mesh_.name() != polyMesh::defaultRegion)
787  {
788  regionStr = "in region '" + mesh_.name() + "' ";
789  }
790 
792  << "Patch '" << patchName << "' not found. "
793  << "Available patch names " << regionStr << "include: " << names()
794  << exit(FatalError);
795  }
796 
797  // Patch not found
798  if (debug)
799  {
800  Pout<< "label polyBoundaryMesh::findPatchID(const word&) const"
801  << "Patch named " << patchName << " not found. "
802  << "List of available patch names: " << names() << endl;
803  }
804 
805  // Not found, return -1
806  return -1;
807 }
808 
809 
810 Foam::label Foam::polyBoundaryMesh::whichPatch(const label faceIndex) const
811 {
812  // Find out which patch the current face belongs to by comparing label
813  // with patch start labels.
814  // If the face is internal, return -1;
815  // if it is off the end of the list, abort
816  if (faceIndex < mesh().nInternalFaces())
817  {
818  return -1;
819  }
820  else if (faceIndex >= mesh().nFaces())
821  {
823  << "Face " << faceIndex
824  << " out of bounds. Number of geometric faces " << mesh().nFaces()
825  << abort(FatalError);
826  }
827 
828 
829  // Patches are ordered, use binary search
830 
831  const polyPatchList& patches = *this;
832 
833  const label patchi =
834  findLower
835  (
836  patches,
837  faceIndex,
838  0,
839  // Must include the start in the comparison
840  [](const polyPatch& p, label val) { return (p.start() <= val); }
841  );
842 
843  if (patchi < 0 || !patches[patchi].range().found(faceIndex))
844  {
845  // If not in any of above, it is trouble!
847  << "Face " << faceIndex << " not found in any of the patches "
848  << flatOutput(names()) << nl
849  << "The patches appear to be inconsistent with the mesh :"
850  << " internalFaces:" << mesh().nInternalFaces()
851  << " total number of faces:" << mesh().nFaces()
852  << abort(FatalError);
853 
854  return -1;
855  }
856 
857  return patchi;
858 }
859 
860 
862 (
863  const UList<wordRe>& patchNames,
864  const bool warnNotFound,
865  const bool useGroups
866 ) const
867 {
868  const wordList allPatchNames(this->names());
869  labelHashSet ids(size());
870 
871  for (const wordRe& patchName : patchNames)
872  {
873  // Treat the given patch names as wild-cards and search the set
874  // of all patch names for matches
875  labelList patchIndices = findStrings(patchName, allPatchNames);
876  ids.insert(patchIndices);
877 
878  if (patchIndices.empty())
879  {
880  if (useGroups)
881  {
882  // Treat as group name or regex for group name
883 
884  const wordList groupNames
885  (
886  groupPatchIDs().tocKeys(patchName)
887  );
888 
889  for (const word& grpName : groupNames)
890  {
891  ids.insert( groupPatchIDs()[grpName] );
892  }
893 
894  if (groupNames.empty() && warnNotFound)
895  {
897  << "Cannot find any patch or group names matching "
898  << patchName
899  << endl;
900  }
901  }
902  else if (warnNotFound)
903  {
905  << "Cannot find any patch names matching " << patchName
906  << endl;
907  }
908  }
909  }
910 
911  return ids;
912 }
913 
914 
916 (
917  const labelUList& patchIDs,
918  wordList& groups,
919  labelHashSet& nonGroupPatches
920 ) const
921 {
922  // Current matched groups
923  DynamicList<word> matchedGroups(1);
924 
925  // Current set of unmatched patches
926  nonGroupPatches = labelHashSet(patchIDs);
927 
928  const HashTable<labelList>& groupPatchIDs = this->groupPatchIDs();
929  forAllConstIters(groupPatchIDs, iter)
930  {
931  // Store currently unmatched patches so we can restore
932  labelHashSet oldNonGroupPatches(nonGroupPatches);
933 
934  // Match by deleting patches in group from the current set and seeing
935  // if all have been deleted.
936  labelHashSet groupPatchSet(iter());
937 
938  label nMatch = nonGroupPatches.erase(groupPatchSet);
939 
940  if (nMatch == groupPatchSet.size())
941  {
942  matchedGroups.append(iter.key());
943  }
944  else if (nMatch != 0)
945  {
946  // No full match. Undo.
947  nonGroupPatches.transfer(oldNonGroupPatches);
948  }
949  }
950 
951  groups.transfer(matchedGroups);
952 }
953 
954 
955 bool Foam::polyBoundaryMesh::checkParallelSync(const bool report) const
956 {
957  if (!Pstream::parRun())
958  {
959  return false;
960  }
961 
962  const polyBoundaryMesh& bm = *this;
963 
964  bool hasError = false;
965 
966  // Collect non-proc patches and check proc patches are last.
967  wordList names(bm.size());
968  wordList types(bm.size());
969 
970  label nonProci = 0;
971 
972  forAll(bm, patchi)
973  {
974  if (!isA<processorPolyPatch>(bm[patchi]))
975  {
976  if (nonProci != patchi)
977  {
978  // There is processor patch in between normal patches.
979  hasError = true;
980 
981  if (debug || report)
982  {
983  Pout<< " ***Problem with boundary patch " << patchi
984  << " named " << bm[patchi].name()
985  << " of type " << bm[patchi].type()
986  << ". The patch seems to be preceeded by processor"
987  << " patches. This is can give problems."
988  << endl;
989  }
990  }
991  else
992  {
993  names[nonProci] = bm[patchi].name();
994  types[nonProci] = bm[patchi].type();
995  nonProci++;
996  }
997  }
998  }
999  names.setSize(nonProci);
1000  types.setSize(nonProci);
1001 
1002  List<wordList> allNames(Pstream::nProcs());
1003  allNames[Pstream::myProcNo()] = names;
1004  Pstream::gatherList(allNames);
1005  Pstream::scatterList(allNames);
1006 
1007  List<wordList> allTypes(Pstream::nProcs());
1008  allTypes[Pstream::myProcNo()] = types;
1009  Pstream::gatherList(allTypes);
1010  Pstream::scatterList(allTypes);
1011 
1012  // Have every processor check but print error on master
1013  // (in processor sequence).
1014 
1015  for (const int proci : Pstream::subProcs())
1016  {
1017  if
1018  (
1019  (allNames[proci] != allNames.first())
1020  || (allTypes[proci] != allTypes.first())
1021  )
1022  {
1023  hasError = true;
1024 
1025  if (debug || (report && Pstream::master()))
1026  {
1027  Info<< " ***Inconsistent patches across processors, "
1028  "processor 0 has patch names:"
1029  << allNames.first()
1030  << " patch types:" << allTypes.first()
1031  << " processor " << proci
1032  << " has patch names:" << allNames[proci]
1033  << " patch types:" << allTypes[proci]
1034  << endl;
1035  }
1036  }
1037  }
1038 
1039  return hasError;
1040 }
1041 
1042 
1043 bool Foam::polyBoundaryMesh::checkDefinition(const bool report) const
1044 {
1045  label nextPatchStart = mesh().nInternalFaces();
1046  const polyBoundaryMesh& bm = *this;
1047 
1048  bool hasError = false;
1049 
1050  wordHashSet patchNames(2*size());
1051 
1052  forAll(bm, patchi)
1053  {
1054  if (bm[patchi].start() != nextPatchStart && !hasError)
1055  {
1056  hasError = true;
1057 
1058  Info<< " ****Problem with boundary patch " << patchi
1059  << " named " << bm[patchi].name()
1060  << " of type " << bm[patchi].type()
1061  << ". The patch should start on face no " << nextPatchStart
1062  << " and the patch specifies " << bm[patchi].start()
1063  << "." << endl
1064  << "Possibly consecutive patches have this same problem."
1065  << " Suppressing future warnings." << endl;
1066  }
1067 
1068  if (!patchNames.insert(bm[patchi].name()) && !hasError)
1069  {
1070  hasError = true;
1071 
1072  Info<< " ****Duplicate boundary patch " << patchi
1073  << " named " << bm[patchi].name()
1074  << " of type " << bm[patchi].type()
1075  << "." << endl
1076  << "Suppressing future warnings." << endl;
1077  }
1078 
1079  nextPatchStart += bm[patchi].size();
1080  }
1081 
1082  reduce(hasError, orOp<bool>());
1083 
1084  if (debug || report)
1085  {
1086  if (hasError)
1087  {
1088  Pout<< " ***Boundary definition is in error." << endl;
1089  }
1090  else
1091  {
1092  Info<< " Boundary definition OK." << endl;
1093  }
1094  }
1095 
1096  return hasError;
1097 }
1098 
1099 
1101 {
1103 
1104  if
1105  (
1108  )
1109  {
1110  forAll(*this, patchi)
1111  {
1112  operator[](patchi).initMovePoints(pBufs, p);
1113  }
1114 
1115  pBufs.finishedSends();
1116 
1117  forAll(*this, patchi)
1118  {
1119  operator[](patchi).movePoints(pBufs, p);
1120  }
1121  }
1123  {
1124  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1125 
1126  // Dummy.
1127  pBufs.finishedSends();
1128 
1129  for (const auto& patchEval : patchSchedule)
1130  {
1131  const label patchi = patchEval.patch;
1132 
1133  if (patchEval.init)
1134  {
1135  operator[](patchi).initMovePoints(pBufs, p);
1136  }
1137  else
1138  {
1139  operator[](patchi).movePoints(pBufs, p);
1140  }
1141  }
1142  }
1143 }
1144 
1145 
1147 {
1148  neighbourEdgesPtr_.clear();
1149  patchIDPtr_.clear();
1150  groupPatchIDsPtr_.clear();
1151 
1153 
1154  if
1155  (
1158  )
1159  {
1160  forAll(*this, patchi)
1161  {
1162  operator[](patchi).initUpdateMesh(pBufs);
1163  }
1164 
1165  pBufs.finishedSends();
1166 
1167  forAll(*this, patchi)
1168  {
1169  operator[](patchi).updateMesh(pBufs);
1170  }
1171  }
1173  {
1174  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1175 
1176  // Dummy.
1177  pBufs.finishedSends();
1178 
1179  for (const auto& patchEval : patchSchedule)
1180  {
1181  const label patchi = patchEval.patch;
1182 
1183  if (patchEval.init)
1184  {
1185  operator[](patchi).initUpdateMesh(pBufs);
1186  }
1187  else
1188  {
1189  operator[](patchi).updateMesh(pBufs);
1190  }
1191  }
1192  }
1193 }
1194 
1195 
1198  const labelUList& oldToNew,
1199  const bool validBoundary
1200 )
1201 {
1202  // Change order of patches
1203  polyPatchList::reorder(oldToNew);
1204 
1205  // Adapt indices
1206  polyPatchList& patches = *this;
1207 
1208  forAll(patches, patchi)
1209  {
1210  patches[patchi].index() = patchi;
1211  }
1212 
1213  if (validBoundary)
1214  {
1215  updateMesh();
1216  }
1217 }
1218 
1219 
1221 {
1222  const polyPatchList& patches = *this;
1223 
1224  os << patches.size() << nl << token::BEGIN_LIST << incrIndent << nl;
1225 
1226  for (const polyPatch& pp : patches)
1227  {
1228  os.beginBlock(pp.name());
1229  os << pp;
1230  os.endBlock();
1231  }
1232 
1233  os << decrIndent << token::END_LIST;
1234 
1235  os.check(FUNCTION_NAME);
1236  return os.good();
1237 }
1238 
1239 
1242  IOstreamOption streamOpt,
1243  const bool valid
1244 ) const
1245 {
1247  return regIOobject::writeObject(streamOpt, valid);
1248 }
1249 
1250 
1251 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
1252 
1253 const Foam::polyPatch& Foam::polyBoundaryMesh::operator[]
1255  const word& patchName
1256 ) const
1257 {
1258  const label patchi = findPatchID(patchName);
1259 
1260  if (patchi < 0)
1261  {
1263  << "Patch named " << patchName << " not found." << nl
1264  << "Available patch names: " << names() << endl
1265  << abort(FatalError);
1266  }
1267 
1268  return operator[](patchi);
1269 }
1270 
1271 
1272 Foam::polyPatch& Foam::polyBoundaryMesh::operator[]
1274  const word& patchName
1275 )
1276 {
1277  const label patchi = findPatchID(patchName);
1278 
1279  if (patchi < 0)
1280  {
1282  << "Patch named " << patchName << " not found." << nl
1283  << "Available patch names: " << names() << endl
1284  << abort(FatalError);
1285  }
1286 
1287  return operator[](patchi);
1288 }
1289 
1290 
1291 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1292 
1294 {
1295  pbm.writeData(os);
1296  return os;
1297 }
1298 
1299 
1300 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOstreamOption::UNCOMPRESSED
compression = false
Definition: IOstreamOption.H:79
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::getTypeOp
General get operation to extract the 'type' from an object as a word.
Definition: ops.H:291
Foam::HashTable< T, edge, Hash< edge > >::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::UPstream::commsTypes::blocking
Foam::polyBoundaryMesh::patchSizes
labelList patchSizes() const
Return a list of patch sizes.
Definition: polyBoundaryMesh.C:627
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::UPstream::commsTypes::nonBlocking
Foam::polyBoundaryMesh::nFaces
label nFaces() const
The number of boundary faces in the underlying mesh.
Definition: polyBoundaryMesh.C:644
Foam::output
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:70
Foam::polyBoundaryMesh::matchGroups
void matchGroups(const labelUList &patchIDs, wordList &groups, labelHashSet &nonGroupPatches) const
Match the patches to groups.
Definition: polyBoundaryMesh.C:916
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::lduSchedule
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:76
Foam::polyBoundaryMesh::physicalTypes
wordList physicalTypes() const
Return a list of physical types.
Definition: polyBoundaryMesh.C:605
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:190
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:318
Foam::DynamicList< label >
Foam::keyType::isPattern
bool isPattern() const
The keyType is treated as a pattern, not as literal string.
Definition: keyTypeI.H:128
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
globalMeshData.H
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::polyBoundaryMesh::findIndex
label findIndex(const keyType &key) const
Return patch index for the first match, return -1 if not found.
Definition: polyBoundaryMesh.C:744
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatch.H:322
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:87
Foam::polyBoundaryMesh::groupPatchIDs
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
Definition: polyBoundaryMesh.C:477
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
Foam::polyBoundaryMesh::reorder
void reorder(const labelUList &oldToNew, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
Definition: polyBoundaryMesh.C:1197
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::HashTable< T, edge, Hash< edge > >::insert
bool insert(const edge &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:61
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::UPstream::defaultCommsType
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:283
Foam::polyBoundaryMesh::clearGeom
void clearGeom()
Clear geometry at this level and at patches.
Definition: polyBoundaryMesh.C:248
Foam::Ostream::beginBlock
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:91
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:327
Foam::polyBoundaryMesh::clearAddressing
void clearAddressing()
Clear addressing at this level and at patches.
Definition: polyBoundaryMesh.C:259
Foam::wordRe
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:72
Foam::polyBoundaryMesh::start
label start() const
The start label of the boundary faces in the polyMesh face list.
Definition: polyBoundaryMesh.C:638
Foam::polyBoundaryMesh::types
wordList types() const
Return a list of patch types.
Definition: polyBoundaryMesh.C:599
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:593
Foam::polyBoundaryMesh::range
labelRange range() const
The face range for all boundary faces.
Definition: polyBoundaryMesh.C:650
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::findIndexImpl
label findIndexImpl(const faPatchList &list, const UnaryMatchPredicate &matcher)
Definition: faBoundaryMesh.C:95
Foam::indicesImpl
static labelList indicesImpl(const faPatchList &list, const UnaryMatchPredicate &matcher)
Definition: faBoundaryMesh.C:68
Foam::labelPair
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:54
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::keyType
A class for handling keywords in dictionaries.
Definition: keyType.H:60
Foam::polyBoundaryMesh::patchID
const labelList & patchID() const
Per boundary face label the patch index.
Definition: polyBoundaryMesh.C:452
Foam::polyBoundaryMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
Definition: polyBoundaryMesh.C:955
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:516
Foam::Field< vector >
Foam::findLower
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Foam::patchIdentifier::inGroups
const wordList & inGroups() const
The (optional) groups that the patch belongs to.
Definition: patchIdentifier.H:170
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::HashTable< T, edge, Hash< edge > >::erase
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:392
Foam::patchIdentifier::inGroup
bool inGroup(const word &name) const
True if the patch is in named group.
Definition: patchIdentifier.H:182
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::polyBoundaryMesh::checkDefinition
bool checkDefinition(const bool report=false) const
Check boundary definition. Return true if in error.
Definition: polyBoundaryMesh.C:1043
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::polyPatchList
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:47
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:474
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:459
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::PtrList::setSize
void setSize(const label newLen)
Same as resize()
Definition: PtrListI.H:105
Foam::polyBoundaryMesh::patchStarts
labelList patchStarts() const
Return a list of patch start face indices.
Definition: polyBoundaryMesh.C:616
Foam::labelRange
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:55
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:810
Foam::reorder
ListType reorder(const labelUList &oldToNew, const ListType &input, const bool prune=false)
Reorder the elements of a list.
Definition: ListOpsTemplates.C:80
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
patchNames
wordList patchNames(nPatches)
Foam::regExpCxx
Wrapper around C++11 regular expressions.
Definition: regExpCxx.H:72
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:766
Foam::PtrList::append
void append(T *ptr)
Append an element to the end of the list.
Definition: PtrListI.H:120
lduSchedule.H
Foam::PstreamBuffers::finishedSends
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
Definition: PstreamBuffers.C:80
Foam::UPstream::commsTypes::scheduled
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Ostream::endBlock
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:109
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
Foam::FatalError
error FatalError
processorPolyPatch.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::regIOobject::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
Definition: regIOobjectWrite.C:37
Foam::ProcessorTopology::patchSchedule
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
Definition: ProcessorTopology.H:101
Foam::getNameOp
General get operation to extract the 'name' from an object as a word.
Definition: ops.H:278
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:217
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:334
Foam::findStrings
labelList findStrings(const regExp &matcher, const UList< StringType > &input, const bool invert=false)
Return list indices for strings matching the regular expression.
Definition: stringListOps.H:76
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::HashTable::sortedToc
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:136
Foam::HashTable< T, edge, Hash< edge > >::find
iterator find(const edge &key)
Find and return an iterator set at the hashed entry.
Definition: HashTableI.H:114
Foam::HashTable::transfer
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
Definition: HashTable.C:671
Foam::polyBoundaryMesh::movePoints
void movePoints(const pointField &p)
Correct polyBoundaryMesh after moving points.
Definition: polyBoundaryMesh.C:1100
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:339
Foam::polyBoundaryMesh::indices
labelList indices(const keyType &key, const bool useGroups=true) const
Return patch indices for all matches.
Definition: polyBoundaryMesh.C:669
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::polyBoundaryMesh::updateMesh
void updateMesh()
Correct polyBoundaryMesh after topology update.
Definition: polyBoundaryMesh.C:1146
Foam::polyBoundaryMesh::nNonProcessor
label nNonProcessor() const
The number of patches before the first processor patch.
Definition: polyBoundaryMesh.C:573
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatch.C:203
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::polyBoundaryMesh::writeData
virtual bool writeData(Ostream &os) const
writeData member function required by regIOobject
Definition: polyBoundaryMesh.C:1220
Foam::EdgeMap
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:51
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:71
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::polyBoundaryMesh::neighbourEdges
const List< labelPairList > & neighbourEdges() const
Per patch the edges on the neighbouring patch.
Definition: polyBoundaryMesh.C:323
edgeHashes.H
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:464
range
scalar range
Definition: LISASMDCalcMethod1.H:12
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Pair< label >
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::polyBoundaryMesh::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
Definition: polyBoundaryMesh.C:1241
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::polyBoundaryMesh::setGroup
void setGroup(const word &groupName, const labelUList &patchIDs)
Set/add group with patches.
Definition: polyBoundaryMesh.C:524
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:862
PstreamBuffers.H
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:181
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
polyBoundaryMesh.H
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:270
Foam::token::END_LIST
End list [isseparator].
Definition: token.H:123
stringListOps.H
Operations on lists of strings.
patchId
label patchId(-1)
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1295
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:224
Foam::IOstreamOption::compression
compressionType compression() const noexcept
Get the stream compression.
Definition: IOstreamOption.H:315
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:409
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:310
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::token::BEGIN_LIST
Begin list [isseparator].
Definition: token.H:122
Foam::orOp
Definition: ops.H:234
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
The patch name.
Definition: patchIdentifier.H:134
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:446
Foam::getMethodImpl
static ListType getMethodImpl(const faPatchList &list, const GetOp &getop)
Definition: faBoundaryMesh.C:47