polyBoundaryMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "polyBoundaryMesh.H"
27 #include "polyMesh.H"
28 #include "primitiveMesh.H"
29 #include "processorPolyPatch.H"
30 #include "stringListOps.H"
31 #include "PstreamBuffers.H"
32 #include "lduSchedule.H"
33 #include "globalMeshData.H"
34 #include "stringListOps.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const IOobject& io,
49  const polyMesh& mesh
50 )
51 :
52  polyPatchList(),
53  regIOobject(io),
54  mesh_(mesh)
55 {
56  if
57  (
60  )
61  {
63  {
65  << "Specified IOobject::MUST_READ_IF_MODIFIED but class"
66  << " does not support automatic rereading."
67  << endl;
68  }
69 
70 
71  polyPatchList& patches = *this;
72 
73  // Read polyPatchList
74  Istream& is = readStream(typeName);
75 
76  PtrList<entry> patchEntries(is);
77  patches.setSize(patchEntries.size());
78 
80  {
81  patches.set
82  (
83  patchi,
85  (
86  patchEntries[patchi].keyword(),
87  patchEntries[patchi].dict(),
88  patchi,
89  *this
90  )
91  );
92  }
93 
94  // Check state of IOstream
95  is.check
96  (
97  "polyBoundaryMesh::polyBoundaryMesh"
98  "(const IOobject&, const polyMesh&)"
99  );
100 
101  close();
102  }
103 }
104 
105 
107 (
108  const IOobject& io,
109  const polyMesh& pm,
110  const label size
111 )
112 :
113  polyPatchList(size),
114  regIOobject(io),
115  mesh_(pm)
116 {}
117 
118 
120 (
121  const IOobject& io,
122  const polyMesh& pm,
123  const polyPatchList& ppl
124 )
125 :
126  polyPatchList(),
127  regIOobject(io),
128  mesh_(pm)
129 {
130  if
131  (
132  (this->readOpt() == IOobject::READ_IF_PRESENT && this->headerOk())
133  || this->readOpt() == IOobject::MUST_READ
135  )
136  {
137 
139  {
141  << "Specified IOobject::MUST_READ_IF_MODIFIED but class"
142  << " does not support automatic rereading."
143  << endl;
144  }
145 
146  polyPatchList& patches = *this;
147 
148  // Read polyPatchList
149  Istream& is = readStream(typeName);
150 
151  PtrList<entry> patchEntries(is);
152  patches.setSize(patchEntries.size());
153 
155  {
156  patches.set
157  (
158  patchi,
160  (
161  patchEntries[patchi].keyword(),
162  patchEntries[patchi].dict(),
163  patchi,
164  *this
165  )
166  );
167  }
168 
169  // Check state of IOstream
170  is.check
171  (
172  "polyBoundaryMesh::polyBoundaryMesh"
173  "(const IOobject&, const polyMesh&, const polyPatchList&)"
174  );
175 
176  close();
177  }
178  else
179  {
180  polyPatchList& patches = *this;
181  patches.setSize(ppl.size());
183  {
184  patches.set(patchi, ppl[patchi].clone(*this).ptr());
185  }
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
191 
193 {}
194 
195 
197 {
198  forAll(*this, patchi)
199  {
200  if (this->set(patchi))
201  {
202  operator[](patchi).clearGeom();
203  }
204  }
205 }
206 
207 
209 {
210  nbrEdgesPtr_.clear();
211  patchIndicesPtr_.clear();
212  patchFaceIndicesPtr_.clear();
213  groupPatchIndicesPtr_.clear();
214 
215  forAll(*this, patchi)
216  {
217  if (this->set(patchi))
218  {
219  operator[](patchi).clearAddressing();
220  }
221  }
222 }
223 
224 
225 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
226 
227 void Foam::polyBoundaryMesh::calcGeometry()
228 {
230 
231  if
232  (
235  )
236  {
237  forAll(*this, patchi)
238  {
239  operator[](patchi).initCalcGeometry(pBufs);
240  }
241 
242  pBufs.finishedSends();
243 
244  forAll(*this, patchi)
245  {
246  operator[](patchi).calcGeometry(pBufs);
247  }
248  }
250  {
251  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
252 
253  // Dummy.
254  pBufs.finishedSends();
255 
256  forAll(patchSchedule, patchEvali)
257  {
258  const label patchi = patchSchedule[patchEvali].patch;
259 
260  if (patchSchedule[patchEvali].init)
261  {
262  operator[](patchi).initCalcGeometry(pBufs);
263  }
264  else
265  {
266  operator[](patchi).calcGeometry(pBufs);
267  }
268  }
269  }
270 }
271 
272 
275 {
276  if (Pstream::parRun())
277  {
279  << "Neighbour edge addressing not correct across parallel"
280  << " boundaries." << endl;
281  }
282 
283  if (!nbrEdgesPtr_.valid())
284  {
285  nbrEdgesPtr_.reset(new List<labelPairList>(size()));
286  List<labelPairList>& nbrEdges = nbrEdgesPtr_();
287 
288  // Initialise.
289  label nEdgePairs = 0;
290  forAll(*this, patchi)
291  {
292  const polyPatch& pp = operator[](patchi);
293 
294  nbrEdges[patchi].setSize(pp.nEdges() - pp.nInternalEdges());
295 
296  forAll(nbrEdges[patchi], i)
297  {
298  labelPair& edgeInfo = nbrEdges[patchi][i];
299 
300  edgeInfo[0] = -1;
301  edgeInfo[1] = -1;
302  }
303 
304  nEdgePairs += pp.nEdges() - pp.nInternalEdges();
305  }
306 
307  // From mesh edge (expressed as a point pair so as not to construct
308  // point addressing) to patch + relative edge index.
309  HashTable<labelPair, edge, Hash<edge>> pointsToEdge(nEdgePairs);
310 
311  forAll(*this, patchi)
312  {
313  const polyPatch& pp = operator[](patchi);
314 
315  const edgeList& edges = pp.edges();
316 
317  for
318  (
319  label edgei = pp.nInternalEdges();
320  edgei < edges.size();
321  edgei++
322  )
323  {
324  // Edge in patch local points
325  const edge& e = edges[edgei];
326 
327  // Edge in mesh points.
328  edge meshEdge(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
329 
331  pointsToEdge.find(meshEdge);
332 
333  if (fnd == pointsToEdge.end())
334  {
335  // First occurrence of mesh edge. Store patch and my
336  // local index.
337  pointsToEdge.insert
338  (
339  meshEdge,
340  labelPair
341  (
342  patchi,
343  edgei - pp.nInternalEdges()
344  )
345  );
346  }
347  else
348  {
349  // Second occurrence. Store.
350  const labelPair& edgeInfo = fnd();
351 
352  nbrEdges[patchi][edgei - pp.nInternalEdges()] =
353  edgeInfo;
354 
355  nbrEdges[edgeInfo[0]][edgeInfo[1]]
356  = labelPair(patchi, edgei - pp.nInternalEdges());
357 
358  // Found all two occurrences of this edge so remove from
359  // hash to save space. Note that this will give lots of
360  // problems if the polyBoundaryMesh is multiply connected.
361  pointsToEdge.erase(meshEdge);
362  }
363  }
364  }
365 
366  if (pointsToEdge.size())
367  {
369  << "Not all boundary edges of patches match up." << nl
370  << "Is the outside of your mesh multiply connected?"
371  << abort(FatalError);
372  }
373 
374  forAll(*this, patchi)
375  {
376  const polyPatch& pp = operator[](patchi);
377 
378  const labelPairList& nbrEdgesp = nbrEdges[patchi];
379 
380  forAll(nbrEdgesp, i)
381  {
382  const labelPair& edgeInfo = nbrEdgesp[i];
383 
384  if (edgeInfo[0] == -1 || edgeInfo[1] == -1)
385  {
386  label edgeI = pp.nInternalEdges() + i;
387  const edge& e = pp.edges()[edgeI];
388 
390  << "Not all boundary edges of patches match up." << nl
391  << "Edge " << edgeI << " on patch " << pp.name()
392  << " end points " << pp.localPoints()[e[0]] << ' '
393  << pp.localPoints()[e[1]] << " is not matched to an"
394  << " edge on any other patch." << nl
395  << "Is the outside of your mesh multiply connected?"
396  << abort(FatalError);
397  }
398  }
399  }
400  }
401 
402  return nbrEdgesPtr_();
403 }
404 
405 
407 {
408  if (!patchIndicesPtr_.valid())
409  {
410  patchIndicesPtr_.reset
411  (
412  new labelList
413  (
414  mesh_.nFaces()
415  - mesh_.nInternalFaces()
416  )
417  );
418  labelList& patchIndices = patchIndicesPtr_();
419 
420  const polyBoundaryMesh& bm = *this;
421 
422  forAll(bm, patchi)
423  {
424  label bFacei = bm[patchi].start() - mesh_.nInternalFaces();
425  forAll(bm[patchi], i)
426  {
427  patchIndices[bFacei++] = patchi;
428  }
429  }
430  }
431  return patchIndicesPtr_();
432 }
433 
434 
436 {
437  if (!patchFaceIndicesPtr_.valid())
438  {
439  patchFaceIndicesPtr_.reset
440  (
441  new labelList
442  (
443  mesh_.nFaces()
444  - mesh_.nInternalFaces()
445  )
446  );
447  labelList& patchFaceID = patchFaceIndicesPtr_();
448 
449  const polyBoundaryMesh& bm = *this;
450 
451  forAll(bm, patchi)
452  {
453  label bFacei = bm[patchi].start() - mesh_.nInternalFaces();
454  forAll(bm[patchi], i)
455  {
456  patchFaceID[bFacei++] = i;
457  }
458  }
459  }
460  return patchFaceIndicesPtr_();
461 }
462 
463 
466 {
467  if (!groupPatchIndicesPtr_.valid())
468  {
469  groupPatchIndicesPtr_.reset(new HashTable<labelList, word>(10));
470  HashTable<labelList, word>& groupPatchIndices = groupPatchIndicesPtr_();
471 
472  const polyBoundaryMesh& bm = *this;
473 
474  forAll(bm, patchi)
475  {
476  const wordList& groups = bm[patchi].inGroups();
477 
478  forAll(groups, i)
479  {
480  const word& name = groups[i];
481 
483  (
484  groupPatchIndices.find
485  (
486  name
487  )
488  );
489 
490  if (iter != groupPatchIndices.end())
491  {
492  iter().append(patchi);
493  }
494  else
495  {
496  groupPatchIndices.insert(name, labelList(1, patchi));
497  }
498  }
499  }
500 
501  // Remove patch names from patchGroups
502  forAll(bm, patchi)
503  {
504  if (groupPatchIndices.erase(bm[patchi].name()))
505  {
507  << "Removing patchGroup '" << bm[patchi].name()
508  << "' which clashes with patch " << patchi
509  << " of the same name."
510  << endl;
511  }
512  }
513  }
514 
515  return groupPatchIndicesPtr_();
516 }
517 
518 
520 (
521  const word& groupName,
522  const labelList& patchIndices
523 )
524 {
525  groupPatchIndicesPtr_.clear();
526 
527  polyPatchList& patches = *this;
528 
529  boolList donePatch(patches.size(), false);
530 
531  // Add to specified patches
532  forAll(patchIndices, i)
533  {
534  label patchi = patchIndices[i];
535  polyPatch& pp = patches[patchi];
536 
537  if (!pp.inGroup(groupName))
538  {
539  pp.inGroups().append(groupName);
540  }
541  donePatch[patchi] = true;
542  }
543 
544  // Remove from other patches
546  {
547  if (!donePatch[patchi])
548  {
549  polyPatch& pp = patches[patchi];
550 
551  label newI = 0;
552  if (pp.inGroup(groupName))
553  {
554  wordList& groups = pp.inGroups();
555 
556  forAll(groups, i)
557  {
558  if (groups[i] != groupName)
559  {
560  groups[newI++] = groups[i];
561  }
562  }
563  groups.setSize(newI);
564  }
565  }
566  }
567 }
568 
569 
571 {
572  const polyPatchList& patches = *this;
573 
574  wordList t(patches.size());
575 
577  {
578  t[patchi] = patches[patchi].name();
579  }
580 
581  return t;
582 }
583 
584 
586 {
587  return toc();
588 }
589 
590 
592 {
593  const polyPatchList& patches = *this;
594 
595  wordList t(patches.size());
596 
598  {
599  t[patchi] = patches[patchi].type();
600  }
601 
602  return t;
603 }
604 
605 
607 {
608  const polyPatchList& patches = *this;
609 
610  wordList t(patches.size());
611 
613  {
614  t[patchi] = patches[patchi].physicalType();
615  }
616 
617  return t;
618 }
619 
620 
622 (
623  const wordRe& key,
624  const bool usePatchGroups
625 ) const
626 {
627  DynamicList<label> indices;
628 
629  if (!key.empty())
630  {
631  if (key.isPattern())
632  {
633  indices = findStrings(key, this->names());
634 
635  if (usePatchGroups && groupPatchIndices().size())
636  {
637  labelHashSet indexSet(indices);
638 
639  const wordList allGroupNames = groupPatchIndices().toc();
640  labelList groupIndices = findStrings(key, allGroupNames);
641  forAll(groupIndices, i)
642  {
643  const word& grpName = allGroupNames[groupIndices[i]];
644  const labelList& patchIndices =
645  groupPatchIndices()[grpName];
646  forAll(patchIndices, j)
647  {
648  if (indexSet.insert(patchIndices[j]))
649  {
650  indices.append(patchIndices[j]);
651  }
652  }
653  }
654  }
655  }
656  else
657  {
658  // Literal string. Special version of above to avoid
659  // unnecessary memory allocations
660 
661  indices.setCapacity(1);
662  forAll(*this, i)
663  {
664  if (key == operator[](i).name())
665  {
666  indices.append(i);
667  break;
668  }
669  }
670 
671  if (usePatchGroups && groupPatchIndices().size())
672  {
674  groupPatchIndices().find(key);
675 
676  if (iter != groupPatchIndices().end())
677  {
678  labelHashSet indexSet(indices);
679 
680  const labelList& patchIndices = iter();
681  forAll(patchIndices, j)
682  {
683  if (indexSet.insert(patchIndices[j]))
684  {
685  indices.append(patchIndices[j]);
686  }
687  }
688  }
689  }
690  }
691  }
692 
693  return indices;
694 }
695 
696 
698 {
699  const polyPatchList& patches = *this;
700 
702  {
703  if (patches[patchi].name() == patchName)
704  {
705  return patchi;
706  }
707  }
708 
709  // Patch not found
710  if (debug)
711  {
712  Pout<< "label polyBoundaryMesh::findIndex(const word&) const"
713  << "Patch named " << patchName << " not found. "
714  << "List of available patch names: " << names() << endl;
715  }
716 
717  // Not found, return -1
718  return -1;
719 }
720 
721 
723 {
724  // Find out which patch the current face belongs to by comparing label
725  // with patch start labels.
726  // If the face is internal, return -1;
727  // if it is off the end of the list, abort
728  if (faceIndex < mesh().nInternalFaces())
729  {
730  return -1;
731  }
732  else if (faceIndex >= mesh().nFaces())
733  {
735  << "given label " << faceIndex
736  << " greater than the number of geometric faces " << mesh().nFaces()
737  << abort(FatalError);
738  }
739 
740 
741  forAll(*this, patchi)
742  {
743  const polyPatch& bp = operator[](patchi);
744 
745  if
746  (
747  faceIndex >= bp.start()
748  && faceIndex < bp.start() + bp.size()
749  )
750  {
751  return patchi;
752  }
753  }
754 
755  // If not in any of above, it is trouble!
757  << "Cannot find face " << faceIndex << " in any of the patches "
758  << names() << nl
759  << "It seems your patches are not consistent with the mesh :"
760  << " internalFaces:" << mesh().nInternalFaces()
761  << " total number of faces:" << mesh().nFaces()
762  << abort(FatalError);
763 
764  return -1;
765 }
766 
767 
769 (
770  const UList<wordRe>& patchNames,
771  const bool warnNotFound,
772  const bool usePatchGroups
773 ) const
774 {
775  const wordList allPatchNames(this->names());
776  labelHashSet ids(size());
777 
778  forAll(patchNames, i)
779  {
780  const wordRe& patchName = patchNames[i];
781 
782  // Treat the given patch names as wild-cards and search the set
783  // of all patch names for matches
784  labelList patchIndices = findStrings(patchName, allPatchNames);
785 
786  forAll(patchIndices, j)
787  {
788  ids.insert(patchIndices[j]);
789  }
790 
791  if (patchIndices.empty())
792  {
793  if (usePatchGroups)
794  {
795  const wordList allGroupNames = groupPatchIndices().toc();
796 
797  // Regard as group name
798  labelList groupIndices = findStrings(patchName, allGroupNames);
799 
800  forAll(groupIndices, i)
801  {
802  const word& name = allGroupNames[groupIndices[i]];
803  const labelList& extraPatchIndices =
804  groupPatchIndices()[name];
805 
806  forAll(extraPatchIndices, extraI)
807  {
808  ids.insert(extraPatchIndices[extraI]);
809  }
810  }
811 
812  if (groupIndices.empty() && warnNotFound)
813  {
815  << "Cannot find any patch or group names matching "
816  << patchName
817  << endl;
818  }
819  }
820  else if (warnNotFound)
821  {
823  << "Cannot find any patch names matching " << patchName
824  << endl;
825  }
826  }
827  }
828 
829  return ids;
830 }
831 
832 
834 (
835  const dictionary& dict,
836  const bool optional
837 ) const
838 {
840 
841  if (dict.found("patch"))
842  {
843  patchNames = List<wordRe>(1, dict.lookup("patch"));
844  }
845  else if (dict.found("patches"))
846  {
847  patchNames = dict.lookup<wordReList>("patches");
848  }
849  else
850  {
851  if (!optional)
852  {
854  << "Neither 'patch' nor 'patches' specified"
855  << exit(FatalIOError);
856  }
857  }
858 
859  labelHashSet patchIDs(patchSet(patchNames));
860 
861  if (!optional && !patchIDs.size())
862  {
864  << "Cannot find any patch named " << patchNames << endl
865  << "Valid patch names are " << names() << endl;
866  }
867 
868  return patchIDs;
869 }
870 
871 
873 (
874  const labelUList& patchIndices,
875  wordList& groups,
876  labelHashSet& nonGroupPatches
877 ) const
878 {
879  // Current matched groups
880  DynamicList<word> matchedGroups(1);
881 
882  // Current set of unmatched patches
883  nonGroupPatches = labelHashSet(patchIndices);
884 
885  const HashTable<labelList, word>& groupPatchIndices =
886  this->groupPatchIndices();
887 
888  for
889  (
891  groupPatchIndices.begin();
892  iter != groupPatchIndices.end();
893  ++iter
894  )
895  {
896  // Store currently unmatched patches so we can restore
897  labelHashSet oldNonGroupPatches(nonGroupPatches);
898 
899  // Match by deleting patches in group from the current set and seeing
900  // if all have been deleted.
901  labelHashSet groupPatchSet(iter());
902 
903  label nMatch = nonGroupPatches.erase(groupPatchSet);
904 
905  if (nMatch == groupPatchSet.size())
906  {
907  matchedGroups.append(iter.key());
908  }
909  else if (nMatch != 0)
910  {
911  // No full match. Undo.
912  nonGroupPatches.transfer(oldNonGroupPatches);
913  }
914  }
915 
916  groups.transfer(matchedGroups);
917 }
918 
919 
920 bool Foam::polyBoundaryMesh::checkParallelSync(const bool report) const
921 {
922  if (!Pstream::parRun())
923  {
924  return false;
925  }
926 
927 
928  const polyBoundaryMesh& bm = *this;
929 
930  bool hasError = false;
931 
932  // Collect non-proc patches and check proc patches are last.
933  wordList names(bm.size());
934  wordList types(bm.size());
935 
936  label nonProci = 0;
937 
938  forAll(bm, patchi)
939  {
940  if (!isA<processorPolyPatch>(bm[patchi]))
941  {
942  if (nonProci != patchi)
943  {
944  // There is processor patch in between normal patches.
945  hasError = true;
946 
947  if (debug || report)
948  {
949  Pout<< " ***Problem with boundary patch " << patchi
950  << " named " << bm[patchi].name()
951  << " of type " << bm[patchi].type()
952  << ". The patch seems to be preceded by processor"
953  << " patches. This is can give problems."
954  << endl;
955  }
956  }
957  else
958  {
959  names[nonProci] = bm[patchi].name();
960  types[nonProci] = bm[patchi].type();
961  nonProci++;
962  }
963  }
964  }
965  names.setSize(nonProci);
966  types.setSize(nonProci);
967 
968  List<wordList> allNames(Pstream::nProcs());
969  allNames[Pstream::myProcNo()] = names;
970  Pstream::gatherList(allNames);
971  Pstream::scatterList(allNames);
972 
973  List<wordList> allTypes(Pstream::nProcs());
974  allTypes[Pstream::myProcNo()] = types;
975  Pstream::gatherList(allTypes);
976  Pstream::scatterList(allTypes);
977 
978  // Have every processor check but only master print error.
979 
980  for (label proci = 1; proci < allNames.size(); ++proci)
981  {
982  if
983  (
984  (allNames[proci] != allNames[0])
985  || (allTypes[proci] != allTypes[0])
986  )
987  {
988  hasError = true;
989 
990  if (debug || (report && Pstream::master()))
991  {
992  Info<< " ***Inconsistent patches across processors, "
993  "processor 0 has patch names:" << allNames[0]
994  << " patch types:" << allTypes[0]
995  << " processor " << proci << " has patch names:"
996  << allNames[proci]
997  << " patch types:" << allTypes[proci]
998  << endl;
999  }
1000  }
1001  }
1002 
1003  return hasError;
1004 }
1005 
1006 
1007 bool Foam::polyBoundaryMesh::checkDefinition(const bool report) const
1008 {
1009  const polyBoundaryMesh& bm = *this;
1010 
1011  bool hasError = false;
1012 
1013  HashSet<word> patchNames(2*size());
1014 
1015  forAll(bm, patchi)
1016  {
1017  if (!patchNames.insert(bm[patchi].name()) && !hasError)
1018  {
1019  hasError = true;
1020 
1021  Info<< " ****Duplicate boundary patch " << patchi
1022  << " named " << bm[patchi].name()
1023  << " of type " << bm[patchi].type()
1024  << "." << endl
1025  << "Suppressing future warnings." << endl;
1026  }
1027  }
1028 
1029  reduce(hasError, orOp<bool>());
1030 
1031  if (debug || report)
1032  {
1033  if (hasError)
1034  {
1035  Pout<< " ***Boundary definition is in error." << endl;
1036  }
1037  else
1038  {
1039  Info<< " Boundary definition OK." << endl;
1040  }
1041  }
1042 
1043  return hasError;
1044 }
1045 
1046 
1048 {
1050 
1051  if
1052  (
1055  )
1056  {
1057  forAll(*this, patchi)
1058  {
1059  operator[](patchi).initMovePoints(pBufs, p);
1060  }
1061 
1062  pBufs.finishedSends();
1063 
1064  forAll(*this, patchi)
1065  {
1066  operator[](patchi).movePoints(pBufs, p);
1067  }
1068  }
1070  {
1071  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1072 
1073  // Dummy.
1074  pBufs.finishedSends();
1075 
1076  forAll(patchSchedule, patchEvali)
1077  {
1078  const label patchi = patchSchedule[patchEvali].patch;
1079 
1080  if (patchSchedule[patchEvali].init)
1081  {
1082  operator[](patchi).initMovePoints(pBufs, p);
1083  }
1084  else
1085  {
1086  operator[](patchi).movePoints(pBufs, p);
1087  }
1088  }
1089  }
1090 }
1091 
1092 
1094 {
1095  nbrEdgesPtr_.clear();
1096  patchIndicesPtr_.clear();
1097  patchFaceIndicesPtr_.clear();
1098  groupPatchIndicesPtr_.clear();
1099 
1101 
1102  if
1103  (
1106  )
1107  {
1108  forAll(*this, patchi)
1109  {
1110  operator[](patchi).initTopoChange(pBufs);
1111  }
1112 
1113  pBufs.finishedSends();
1114 
1115  forAll(*this, patchi)
1116  {
1117  operator[](patchi).topoChange(pBufs);
1118  }
1119  }
1121  {
1122  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1123 
1124  // Dummy.
1125  pBufs.finishedSends();
1126 
1127  forAll(patchSchedule, patchEvali)
1128  {
1129  const label patchi = patchSchedule[patchEvali].patch;
1130 
1131  if (patchSchedule[patchEvali].init)
1132  {
1133  operator[](patchi).initTopoChange(pBufs);
1134  }
1135  else
1136  {
1137  operator[](patchi).topoChange(pBufs);
1138  }
1139  }
1140  }
1141 }
1142 
1143 
1145 (
1146  const wordUList& newNames,
1147  const bool validBoundary
1148 )
1149 {
1150  polyPatchList& patches = *this;
1152  {
1153  if (patches.set(patchi))
1154  {
1155  patches[patchi].rename(newNames);
1156  }
1157  }
1158 
1159  if (validBoundary)
1160  {
1161  topoChange();
1162  }
1163 }
1164 
1165 
1167 (
1168  const labelUList& newToOld,
1169  const bool validBoundary
1170 )
1171 {
1172  // Change order of patches
1173  polyPatchList::shuffle(newToOld);
1174 
1175  // Adapt indices
1176  polyPatchList& patches = *this;
1178  {
1179  if (patches.set(patchi))
1180  {
1181  patches[patchi].reorder(newToOld);
1182  }
1183  }
1184 
1185  if (validBoundary)
1186  {
1187  topoChange();
1188  }
1189 }
1190 
1191 
1193 {
1194  const polyPatchList& patches = *this;
1195 
1196  os << patches.size() << nl << token::BEGIN_LIST << incrIndent << nl;
1197 
1199  {
1200  os << indent << patches[patchi].name() << nl
1201  << indent << token::BEGIN_BLOCK << nl
1203  << indent << token::END_BLOCK << endl;
1204  }
1205 
1206  os << decrIndent << token::END_LIST;
1207 
1208  // Check state of IOstream
1209  os.check("polyBoundaryMesh::writeData(Ostream& os) const");
1210 
1211  return os.good();
1212 }
1213 
1214 
1216 (
1220  const bool write
1221 ) const
1222 {
1224 }
1225 
1226 
1227 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
1228 
1230 (
1231  const word& patchName
1232 ) const
1233 {
1234  const label patchi = findIndex(patchName);
1235 
1236  if (patchi < 0)
1237  {
1239  << "Patch named " << patchName << " not found." << nl
1240  << "Available patch names: " << names() << endl
1241  << abort(FatalError);
1242  }
1243 
1244  return operator[](patchi);
1245 }
1246 
1247 
1249 (
1250  const word& patchName
1251 )
1252 {
1253  const label patchi = findIndex(patchName);
1254 
1255  if (patchi < 0)
1256  {
1258  << "Patch named " << patchName << " not found." << nl
1259  << "Available patch names: " << names() << endl
1260  << abort(FatalError);
1261  }
1262 
1263  return operator[](patchi);
1264 }
1265 
1266 
1267 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1268 
1270 {
1271  pbm.writeData(os);
1272  return os;
1273 }
1274 
1275 
1276 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
void setCapacity(const label)
Alter the size of the underlying storage.
Definition: DynamicListI.H:130
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
An STL-conforming const_iterator.
Definition: HashTable.H:498
An STL-conforming iterator.
Definition: HashTable.H:443
An STL-conforming hash table.
Definition: HashTable.H:127
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:445
void transfer(HashTable< T, Key, Hash > &)
Transfer the contents of the argument table into this table.
Definition: HashTable.C:587
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:167
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
readOption readOpt() const
Definition: IOobject.H:357
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:280
const word & name() const
Return name.
Definition: IOobject.H:307
Version number type.
Definition: IOstream.H:97
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:87
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
label nEdges() const
Return number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
void reorder(const labelUList &oldToNew)
Reorders elements. Ordering does not have to be done in.
Definition: PtrList.C:221
void shuffle(const labelUList &newToOld)
Reorders elements. Ordering does not have to be done in.
Definition: PtrList.C:272
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
An STL iterator.
Definition: UPtrList.H:206
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
friend Ostream & operator(Ostream &, const UPtrList< T > &)
Write UPtrList to Ostream.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
bool inGroup(const word &) const
Test if in group.
const wordList & inGroups() const
Return the optional groups patch belongs to.
const word & name() const
Return name.
Foam::polyBoundaryMesh.
bool checkDefinition(const bool report=false) const
Check boundary definition. Return true if in error.
const labelList & patchIndices() const
Boundary face patch indices.
wordList physicalTypes() const
Return a list of physical types.
void topoChange()
Correct polyBoundaryMesh after topology update.
bool writeData(Ostream &) const
writeData member function required by regIOobject
label findIndex(const word &patchName) const
Find patch index given a name.
const HashTable< labelList, word > & groupPatchIndices() const
Per patch group the patch indices.
wordList types() const
Return a list of patch types.
void clearGeom()
Clear geometry at this level and at patches.
void setGroup(const word &groupName, const labelList &patchIndices)
Set/add group with patches.
const labelList & patchFaceIndices() const
Boundary face patch face indices.
~polyBoundaryMesh()
Destructor.
bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write using given format, version and form uncompression.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
void clearAddressing()
Clear addressing at this level and at patches.
void matchGroups(const labelUList &patchIndices, wordList &groups, labelHashSet &nonGroupPatches) const
Match the patches to groups. Returns all the (fully matched) groups.
void movePoints(const pointField &)
Correct polyBoundaryMesh after moving points.
void renamePatches(const wordUList &newNames, const bool validBoundary)
Rename the patches. If validBoundary is set this calls topoChange()
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the patch set corresponding to the given names.
labelHashSet findIndices() const
Find patch indices for a given polyPatch type.
void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorders the patches. Ordering does not have to be done in.
polyBoundaryMesh(const IOobject &, const polyMesh &)
Read constructor given IOobject and a polyMesh reference.
wordList names() const
Return the list of patch names.
wordList toc() const
Return the list of patch names.
const List< labelPairList > & nbrEdges() const
Per patch the edges on the neighbouring patch. Is for every external.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
label nInternalFaces() const
label nFaces() const
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
void close()
Close Istream.
virtual bool writeObject(IOstream::streamFormat, IOstream::versionNumber, IOstream::compressionType, const bool write) const
Write using given format, version and compression.
bool headerOk()
Read and check header info.
Definition: regIOobject.C:453
Istream & readStream(const word &, const bool read=true)
Return Istream and check object type against that given.
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
@ BEGIN_LIST
Definition: token.H:109
@ END_LIST
Definition: token.H:110
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:77
static bool isPattern(const string &)
Test string for regular expression meta characters.
Definition: wordReI.H:34
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const doubleScalar e
Definition: doubleScalar.H:106
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:74
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
edge meshEdge(const PrimitivePatch< FaceList, PointField > &p, const label edgei)
defineTypeNameAndDebug(combustionModel, 0)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
wordList patchNames(nPatches)
dictionary dict
volScalarField & p
Operations on lists of strings.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112