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-2020 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 {
40  defineTypeNameAndDebug(polyBoundaryMesh, 0);
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 Foam::polyBoundaryMesh::polyBoundaryMesh
47 (
48  const IOobject& io,
49  const polyMesh& mesh
50 )
51 :
52  polyPatchList(),
53  regIOobject(io),
54  mesh_(mesh)
55 {
56  if
57  (
58  readOpt() == IOobject::MUST_READ
59  || readOpt() == IOobject::MUST_READ_IF_MODIFIED
60  )
61  {
62  if (readOpt() == IOobject::MUST_READ_IF_MODIFIED)
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 
79  forAll(patches, patchi)
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 
106 Foam::polyBoundaryMesh::polyBoundaryMesh
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 
119 Foam::polyBoundaryMesh::polyBoundaryMesh
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
134  || this->readOpt() == IOobject::MUST_READ_IF_MODIFIED
135  )
136  {
137 
138  if (readOpt() == IOobject::MUST_READ_IF_MODIFIED)
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 
154  forAll(patches, patchi)
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());
182  forAll(patches, patchi)
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  patchIDPtr_.clear();
212  groupPatchIDsPtr_.clear();
213 
214  forAll(*this, patchi)
215  {
216  if (this->set(patchi))
217  {
218  operator[](patchi).clearAddressing();
219  }
220  }
221 }
222 
223 
224 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225 
226 void Foam::polyBoundaryMesh::calcGeometry()
227 {
229 
230  if
231  (
234  )
235  {
236  forAll(*this, patchi)
237  {
238  operator[](patchi).initCalcGeometry(pBufs);
239  }
240 
241  pBufs.finishedSends();
242 
243  forAll(*this, patchi)
244  {
245  operator[](patchi).calcGeometry(pBufs);
246  }
247  }
249  {
250  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
251 
252  // Dummy.
253  pBufs.finishedSends();
254 
255  forAll(patchSchedule, patchEvali)
256  {
257  const label patchi = patchSchedule[patchEvali].patch;
258 
259  if (patchSchedule[patchEvali].init)
260  {
261  operator[](patchi).initCalcGeometry(pBufs);
262  }
263  else
264  {
265  operator[](patchi).calcGeometry(pBufs);
266  }
267  }
268  }
269 }
270 
271 
274 {
275  if (Pstream::parRun())
276  {
278  << "Neighbour edge addressing not correct across parallel"
279  << " boundaries." << endl;
280  }
281 
282  if (!nbrEdgesPtr_.valid())
283  {
284  nbrEdgesPtr_.reset(new List<labelPairList>(size()));
285  List<labelPairList>& nbrEdges = nbrEdgesPtr_();
286 
287  // Initialize.
288  label nEdgePairs = 0;
289  forAll(*this, patchi)
290  {
291  const polyPatch& pp = operator[](patchi);
292 
293  nbrEdges[patchi].setSize(pp.nEdges() - pp.nInternalEdges());
294 
295  forAll(nbrEdges[patchi], i)
296  {
297  labelPair& edgeInfo = nbrEdges[patchi][i];
298 
299  edgeInfo[0] = -1;
300  edgeInfo[1] = -1;
301  }
302 
303  nEdgePairs += pp.nEdges() - pp.nInternalEdges();
304  }
305 
306  // From mesh edge (expressed as a point pair so as not to construct
307  // point addressing) to patch + relative edge index.
308  HashTable<labelPair, edge, Hash<edge>> pointsToEdge(nEdgePairs);
309 
310  forAll(*this, patchi)
311  {
312  const polyPatch& pp = operator[](patchi);
313 
314  const edgeList& edges = pp.edges();
315 
316  for
317  (
318  label edgei = pp.nInternalEdges();
319  edgei < edges.size();
320  edgei++
321  )
322  {
323  // Edge in patch local points
324  const edge& e = edges[edgei];
325 
326  // Edge in mesh points.
327  edge meshEdge(pp.meshPoints()[e[0]], pp.meshPoints()[e[1]]);
328 
330  pointsToEdge.find(meshEdge);
331 
332  if (fnd == pointsToEdge.end())
333  {
334  // First occurrence of mesh edge. Store patch and my
335  // local index.
336  pointsToEdge.insert
337  (
338  meshEdge,
339  labelPair
340  (
341  patchi,
342  edgei - pp.nInternalEdges()
343  )
344  );
345  }
346  else
347  {
348  // Second occurrence. Store.
349  const labelPair& edgeInfo = fnd();
350 
351  nbrEdges[patchi][edgei - pp.nInternalEdges()] =
352  edgeInfo;
353 
354  nbrEdges[edgeInfo[0]][edgeInfo[1]]
355  = labelPair(patchi, edgei - pp.nInternalEdges());
356 
357  // Found all two occurrences of this edge so remove from
358  // hash to save space. Note that this will give lots of
359  // problems if the polyBoundaryMesh is multiply connected.
360  pointsToEdge.erase(meshEdge);
361  }
362  }
363  }
364 
365  if (pointsToEdge.size())
366  {
368  << "Not all boundary edges of patches match up." << nl
369  << "Is the outside of your mesh multiply connected?"
370  << abort(FatalError);
371  }
372 
373  forAll(*this, patchi)
374  {
375  const polyPatch& pp = operator[](patchi);
376 
377  const labelPairList& nbrEdgesp = nbrEdges[patchi];
378 
379  forAll(nbrEdgesp, i)
380  {
381  const labelPair& edgeInfo = nbrEdgesp[i];
382 
383  if (edgeInfo[0] == -1 || edgeInfo[1] == -1)
384  {
385  label edgeI = pp.nInternalEdges() + i;
386  const edge& e = pp.edges()[edgeI];
387 
389  << "Not all boundary edges of patches match up." << nl
390  << "Edge " << edgeI << " on patch " << pp.name()
391  << " end points " << pp.localPoints()[e[0]] << ' '
392  << pp.localPoints()[e[1]] << " is not matched to an"
393  << " edge on any other patch." << nl
394  << "Is the outside of your mesh multiply connected?"
395  << abort(FatalError);
396  }
397  }
398  }
399  }
400 
401  return nbrEdgesPtr_();
402 }
403 
404 
406 {
407  if (!patchIDPtr_.valid())
408  {
409  patchIDPtr_.reset
410  (
411  new labelList
412  (
413  mesh_.nFaces()
414  - mesh_.nInternalFaces()
415  )
416  );
417  labelList& patchID = patchIDPtr_();
418 
419  const polyBoundaryMesh& bm = *this;
420 
421  forAll(bm, patchi)
422  {
423  label bFacei = bm[patchi].start() - mesh_.nInternalFaces();
424  forAll(bm[patchi], i)
425  {
426  patchID[bFacei++] = patchi;
427  }
428  }
429  }
430  return patchIDPtr_();
431 }
432 
433 
436 {
437  if (!groupPatchIDsPtr_.valid())
438  {
439  groupPatchIDsPtr_.reset(new HashTable<labelList, word>(10));
440  HashTable<labelList, word>& groupPatchIDs = groupPatchIDsPtr_();
441 
442  const polyBoundaryMesh& bm = *this;
443 
444  forAll(bm, patchi)
445  {
446  const wordList& groups = bm[patchi].inGroups();
447 
448  forAll(groups, i)
449  {
450  const word& name = groups[i];
451 
452  HashTable<labelList, word>::iterator iter = groupPatchIDs.find
453  (
454  name
455  );
456 
457  if (iter != groupPatchIDs.end())
458  {
459  iter().append(patchi);
460  }
461  else
462  {
463  groupPatchIDs.insert(name, labelList(1, patchi));
464  }
465  }
466  }
467 
468  // Remove patch names from patchGroups
469  forAll(bm, patchi)
470  {
471  if (groupPatchIDs.erase(bm[patchi].name()))
472  {
474  << "Removing patchGroup '" << bm[patchi].name()
475  << "' which clashes with patch " << patchi
476  << " of the same name."
477  << endl;
478  }
479  }
480  }
481 
482  return groupPatchIDsPtr_();
483 }
484 
485 
487 (
488  const word& groupName,
489  const labelList& patchIDs
490 )
491 {
492  groupPatchIDsPtr_.clear();
493 
494  polyPatchList& patches = *this;
495 
496  boolList donePatch(patches.size(), false);
497 
498  // Add to specified patches
499  forAll(patchIDs, i)
500  {
501  label patchi = patchIDs[i];
502  polyPatch& pp = patches[patchi];
503 
504  if (!pp.inGroup(groupName))
505  {
506  pp.inGroups().append(groupName);
507  }
508  donePatch[patchi] = true;
509  }
510 
511  // Remove from other patches
512  forAll(patches, patchi)
513  {
514  if (!donePatch[patchi])
515  {
516  polyPatch& pp = patches[patchi];
517 
518  label newI = 0;
519  if (pp.inGroup(groupName))
520  {
521  wordList& groups = pp.inGroups();
522 
523  forAll(groups, i)
524  {
525  if (groups[i] != groupName)
526  {
527  groups[newI++] = groups[i];
528  }
529  }
530  groups.setSize(newI);
531  }
532  }
533  }
534 }
535 
536 
538 {
539  const polyPatchList& patches = *this;
540 
541  wordList t(patches.size());
542 
543  forAll(patches, patchi)
544  {
545  t[patchi] = patches[patchi].name();
546  }
547 
548  return t;
549 }
550 
551 
553 {
554  const polyPatchList& patches = *this;
555 
556  wordList t(patches.size());
557 
558  forAll(patches, patchi)
559  {
560  t[patchi] = patches[patchi].type();
561  }
562 
563  return t;
564 }
565 
566 
568 {
569  const polyPatchList& patches = *this;
570 
571  wordList t(patches.size());
572 
573  forAll(patches, patchi)
574  {
575  t[patchi] = patches[patchi].physicalType();
576  }
577 
578  return t;
579 }
580 
581 
583 (
584  const keyType& key,
585  const bool usePatchGroups
586 ) const
587 {
588  DynamicList<label> indices;
589 
590  if (!key.empty())
591  {
592  if (key.isPattern())
593  {
594  indices = findStrings(key, this->names());
595 
596  if (usePatchGroups && groupPatchIDs().size())
597  {
598  labelHashSet indexSet(indices);
599 
600  const wordList allGroupNames = groupPatchIDs().toc();
601  labelList groupIDs = findStrings(key, allGroupNames);
602  forAll(groupIDs, i)
603  {
604  const word& grpName = allGroupNames[groupIDs[i]];
605  const labelList& patchIDs = groupPatchIDs()[grpName];
606  forAll(patchIDs, j)
607  {
608  if (indexSet.insert(patchIDs[j]))
609  {
610  indices.append(patchIDs[j]);
611  }
612  }
613  }
614  }
615  }
616  else
617  {
618  // Literal string. Special version of above to avoid
619  // unnecessary memory allocations
620 
621  indices.setCapacity(1);
622  forAll(*this, i)
623  {
624  if (key == operator[](i).name())
625  {
626  indices.append(i);
627  break;
628  }
629  }
630 
631  if (usePatchGroups && groupPatchIDs().size())
632  {
634  groupPatchIDs().find(key);
635 
636  if (iter != groupPatchIDs().end())
637  {
638  labelHashSet indexSet(indices);
639 
640  const labelList& patchIDs = iter();
641  forAll(patchIDs, j)
642  {
643  if (indexSet.insert(patchIDs[j]))
644  {
645  indices.append(patchIDs[j]);
646  }
647  }
648  }
649  }
650  }
651  }
652 
653  return move(indices);
654 }
655 
656 
658 {
659  if (!key.empty())
660  {
661  if (key.isPattern())
662  {
663  labelList indices = this->findIndices(key);
664 
665  // return first element
666  if (!indices.empty())
667  {
668  return indices[0];
669  }
670  }
671  else
672  {
673  forAll(*this, i)
674  {
675  if (key == operator[](i).name())
676  {
677  return i;
678  }
679  }
680  }
681  }
682 
683  // not found
684  return -1;
685 }
686 
687 
689 {
690  const polyPatchList& patches = *this;
691 
692  forAll(patches, patchi)
693  {
694  if (patches[patchi].name() == patchName)
695  {
696  return patchi;
697  }
698  }
699 
700  // Patch not found
701  if (debug)
702  {
703  Pout<< "label polyBoundaryMesh::findPatchID(const word&) const"
704  << "Patch named " << patchName << " not found. "
705  << "List of available patch names: " << names() << endl;
706  }
707 
708  // Not found, return -1
709  return -1;
710 }
711 
712 
714 {
715  // Find out which patch the current face belongs to by comparing label
716  // with patch start labels.
717  // If the face is internal, return -1;
718  // if it is off the end of the list, abort
719  if (faceIndex < mesh().nInternalFaces())
720  {
721  return -1;
722  }
723  else if (faceIndex >= mesh().nFaces())
724  {
726  << "given label " << faceIndex
727  << " greater than the number of geometric faces " << mesh().nFaces()
728  << abort(FatalError);
729  }
730 
731 
732  forAll(*this, patchi)
733  {
734  const polyPatch& bp = operator[](patchi);
735 
736  if
737  (
738  faceIndex >= bp.start()
739  && faceIndex < bp.start() + bp.size()
740  )
741  {
742  return patchi;
743  }
744  }
745 
746  // If not in any of above, it is trouble!
748  << "Cannot find face " << faceIndex << " in any of the patches "
749  << names() << nl
750  << "It seems your patches are not consistent with the mesh :"
751  << " internalFaces:" << mesh().nInternalFaces()
752  << " total number of faces:" << mesh().nFaces()
753  << abort(FatalError);
754 
755  return -1;
756 }
757 
758 
760 (
761  const UList<wordRe>& patchNames,
762  const bool warnNotFound,
763  const bool usePatchGroups
764 ) const
765 {
766  const wordList allPatchNames(this->names());
767  labelHashSet ids(size());
768 
769  forAll(patchNames, i)
770  {
771  const wordRe& patchName = patchNames[i];
772 
773  // Treat the given patch names as wild-cards and search the set
774  // of all patch names for matches
775  labelList patchIDs = findStrings(patchName, allPatchNames);
776 
777  forAll(patchIDs, j)
778  {
779  ids.insert(patchIDs[j]);
780  }
781 
782  if (patchIDs.empty())
783  {
784  if (usePatchGroups)
785  {
786  const wordList allGroupNames = groupPatchIDs().toc();
787 
788  // Regard as group name
789  labelList groupIDs = findStrings(patchName, allGroupNames);
790 
791  forAll(groupIDs, i)
792  {
793  const word& name = allGroupNames[groupIDs[i]];
794  const labelList& extraPatchIDs = groupPatchIDs()[name];
795 
796  forAll(extraPatchIDs, extraI)
797  {
798  ids.insert(extraPatchIDs[extraI]);
799  }
800  }
801 
802  if (groupIDs.empty() && warnNotFound)
803  {
805  << "Cannot find any patch or group names matching "
806  << patchName
807  << endl;
808  }
809  }
810  else if (warnNotFound)
811  {
813  << "Cannot find any patch names matching " << patchName
814  << endl;
815  }
816  }
817  }
818 
819  return ids;
820 }
821 
822 
824 (
825  const labelUList& patchIDs,
826  wordList& groups,
827  labelHashSet& nonGroupPatches
828 ) const
829 {
830  // Current matched groups
831  DynamicList<word> matchedGroups(1);
832 
833  // Current set of unmatched patches
834  nonGroupPatches = labelHashSet(patchIDs);
835 
836  const HashTable<labelList, word>& groupPatchIDs = this->groupPatchIDs();
837  for
838  (
840  groupPatchIDs.begin();
841  iter != groupPatchIDs.end();
842  ++iter
843  )
844  {
845  // Store currently unmatched patches so we can restore
846  labelHashSet oldNonGroupPatches(nonGroupPatches);
847 
848  // Match by deleting patches in group from the current set and seeing
849  // if all have been deleted.
850  labelHashSet groupPatchSet(iter());
851 
852  label nMatch = nonGroupPatches.erase(groupPatchSet);
853 
854  if (nMatch == groupPatchSet.size())
855  {
856  matchedGroups.append(iter.key());
857  }
858  else if (nMatch != 0)
859  {
860  // No full match. Undo.
861  nonGroupPatches.transfer(oldNonGroupPatches);
862  }
863  }
864 
865  groups.transfer(matchedGroups);
866 }
867 
868 
869 bool Foam::polyBoundaryMesh::checkParallelSync(const bool report) const
870 {
871  if (!Pstream::parRun())
872  {
873  return false;
874  }
875 
876 
877  const polyBoundaryMesh& bm = *this;
878 
879  bool hasError = false;
880 
881  // Collect non-proc patches and check proc patches are last.
882  wordList names(bm.size());
883  wordList types(bm.size());
884 
885  label nonProci = 0;
886 
887  forAll(bm, patchi)
888  {
889  if (!isA<processorPolyPatch>(bm[patchi]))
890  {
891  if (nonProci != patchi)
892  {
893  // There is processor patch in between normal patches.
894  hasError = true;
895 
896  if (debug || report)
897  {
898  Pout<< " ***Problem with boundary patch " << patchi
899  << " named " << bm[patchi].name()
900  << " of type " << bm[patchi].type()
901  << ". The patch seems to be preceded by processor"
902  << " patches. This is can give problems."
903  << endl;
904  }
905  }
906  else
907  {
908  names[nonProci] = bm[patchi].name();
909  types[nonProci] = bm[patchi].type();
910  nonProci++;
911  }
912  }
913  }
914  names.setSize(nonProci);
915  types.setSize(nonProci);
916 
917  List<wordList> allNames(Pstream::nProcs());
918  allNames[Pstream::myProcNo()] = names;
919  Pstream::gatherList(allNames);
920  Pstream::scatterList(allNames);
921 
922  List<wordList> allTypes(Pstream::nProcs());
923  allTypes[Pstream::myProcNo()] = types;
924  Pstream::gatherList(allTypes);
925  Pstream::scatterList(allTypes);
926 
927  // Have every processor check but only master print error.
928 
929  for (label proci = 1; proci < allNames.size(); ++proci)
930  {
931  if
932  (
933  (allNames[proci] != allNames[0])
934  || (allTypes[proci] != allTypes[0])
935  )
936  {
937  hasError = true;
938 
939  if (debug || (report && Pstream::master()))
940  {
941  Info<< " ***Inconsistent patches across processors, "
942  "processor 0 has patch names:" << allNames[0]
943  << " patch types:" << allTypes[0]
944  << " processor " << proci << " has patch names:"
945  << allNames[proci]
946  << " patch types:" << allTypes[proci]
947  << endl;
948  }
949  }
950  }
951 
952  return hasError;
953 }
954 
955 
956 bool Foam::polyBoundaryMesh::checkDefinition(const bool report) const
957 {
958  label nextPatchStart = mesh().nInternalFaces();
959  const polyBoundaryMesh& bm = *this;
960 
961  bool hasError = false;
962 
963  HashSet<word> patchNames(2*size());
964 
965  forAll(bm, patchi)
966  {
967  if (bm[patchi].start() != nextPatchStart && !hasError)
968  {
969  hasError = true;
970 
971  Info<< " ****Problem with boundary patch " << patchi
972  << " named " << bm[patchi].name()
973  << " of type " << bm[patchi].type()
974  << ". The patch should start on face no " << nextPatchStart
975  << " and the patch specifies " << bm[patchi].start()
976  << "." << endl
977  << "Possibly consecutive patches have this same problem."
978  << " Suppressing future warnings." << endl;
979  }
980 
981  if (!patchNames.insert(bm[patchi].name()) && !hasError)
982  {
983  hasError = true;
984 
985  Info<< " ****Duplicate boundary patch " << patchi
986  << " named " << bm[patchi].name()
987  << " of type " << bm[patchi].type()
988  << "." << endl
989  << "Suppressing future warnings." << endl;
990  }
991 
992  nextPatchStart += bm[patchi].size();
993  }
994 
995  reduce(hasError, orOp<bool>());
996 
997  if (debug || report)
998  {
999  if (hasError)
1000  {
1001  Pout<< " ***Boundary definition is in error." << endl;
1002  }
1003  else
1004  {
1005  Info<< " Boundary definition OK." << endl;
1006  }
1007  }
1008 
1009  return hasError;
1010 }
1011 
1012 
1014 {
1016 
1017  if
1018  (
1021  )
1022  {
1023  forAll(*this, patchi)
1024  {
1025  operator[](patchi).initMovePoints(pBufs, p);
1026  }
1027 
1028  pBufs.finishedSends();
1029 
1030  forAll(*this, patchi)
1031  {
1032  operator[](patchi).movePoints(pBufs, p);
1033  }
1034  }
1036  {
1037  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1038 
1039  // Dummy.
1040  pBufs.finishedSends();
1041 
1042  forAll(patchSchedule, patchEvali)
1043  {
1044  const label patchi = patchSchedule[patchEvali].patch;
1045 
1046  if (patchSchedule[patchEvali].init)
1047  {
1048  operator[](patchi).initMovePoints(pBufs, p);
1049  }
1050  else
1051  {
1052  operator[](patchi).movePoints(pBufs, p);
1053  }
1054  }
1055  }
1056 }
1057 
1058 
1060 {
1061  nbrEdgesPtr_.clear();
1062  patchIDPtr_.clear();
1063  groupPatchIDsPtr_.clear();
1064 
1066 
1067  if
1068  (
1071  )
1072  {
1073  forAll(*this, patchi)
1074  {
1075  operator[](patchi).initUpdateMesh(pBufs);
1076  }
1077 
1078  pBufs.finishedSends();
1079 
1080  forAll(*this, patchi)
1081  {
1082  operator[](patchi).updateMesh(pBufs);
1083  }
1084  }
1086  {
1087  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
1088 
1089  // Dummy.
1090  pBufs.finishedSends();
1091 
1092  forAll(patchSchedule, patchEvali)
1093  {
1094  const label patchi = patchSchedule[patchEvali].patch;
1095 
1096  if (patchSchedule[patchEvali].init)
1097  {
1098  operator[](patchi).initUpdateMesh(pBufs);
1099  }
1100  else
1101  {
1102  operator[](patchi).updateMesh(pBufs);
1103  }
1104  }
1105  }
1106 }
1107 
1108 
1111  const labelUList& newToOld,
1112  const bool validBoundary
1113 )
1114 {
1115  // Change order of patches
1116  polyPatchList::shuffle(newToOld);
1117 
1118  // Adapt indices
1119  polyPatchList& patches = *this;
1120 
1121  forAll(patches, patchi)
1122  {
1123  if (patches.set(patchi))
1124  {
1125  patches[patchi].index() = patchi;
1126  }
1127  }
1128 
1129  if (validBoundary)
1130  {
1131  updateMesh();
1132  }
1133 }
1134 
1135 
1137 {
1138  const polyPatchList& patches = *this;
1139 
1140  os << patches.size() << nl << token::BEGIN_LIST << incrIndent << nl;
1141 
1142  forAll(patches, patchi)
1143  {
1144  os << indent << patches[patchi].name() << nl
1145  << indent << token::BEGIN_BLOCK << nl
1146  << incrIndent << patches[patchi] << decrIndent
1147  << indent << token::END_BLOCK << endl;
1148  }
1149 
1150  os << decrIndent << token::END_LIST;
1151 
1152  // Check state of IOstream
1153  os.check("polyBoundaryMesh::writeData(Ostream& os) const");
1154 
1155  return os.good();
1156 }
1157 
1158 
1164  const bool write
1165 ) const
1166 {
1167  return regIOobject::writeObject(fmt, ver, IOstream::UNCOMPRESSED, write);
1168 }
1169 
1170 
1171 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
1172 
1175  const word& patchName
1176 ) const
1177 {
1178  const label patchi = findPatchID(patchName);
1179 
1180  if (patchi < 0)
1181  {
1183  << "Patch named " << patchName << " not found." << nl
1184  << "Available patch names: " << names() << endl
1185  << abort(FatalError);
1186  }
1187 
1188  return operator[](patchi);
1189 }
1190 
1191 
1194  const word& patchName
1195 )
1196 {
1197  const label patchi = findPatchID(patchName);
1198 
1199  if (patchi < 0)
1200  {
1202  << "Patch named " << patchName << " not found." << nl
1203  << "Available patch names: " << names() << endl
1204  << abort(FatalError);
1205  }
1206 
1207  return operator[](patchi);
1208 }
1209 
1210 
1211 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1212 
1214 {
1215  pbm.writeData(os);
1216  return os;
1217 }
1218 
1219 
1220 // ************************************************************************* //
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
A class for handling keywords in dictionaries.
Definition: keyType.H:66
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
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 writeData(Ostream &) const
writeData member function required by regIOobject
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:303
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
An STL-conforming const_iterator.
Definition: HashTable.H:481
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
label nInternalFaces() const
bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write using given format, version and form uncompression.
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
label nFaces() const
bool inGroup(const word &) const
Test if in group.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void shuffle(const labelUList &newToOld, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const labelList & patchID() const
Per boundary face label the patch index.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
label findPatchID(const word &patchName) const
Find patch index given a name.
label nInternalEdges() const
Number of internal edges.
void clearGeom()
Clear geometry at this level and at patches.
const wordList & inGroups() const
Return the optional groups patch belongs to.
Operations on lists of strings.
void transfer(HashTable< T, Key, Hash > &)
Transfer the contents of the argument table into this table.
Definition: HashTable.C:513
patches[0]
const Field< PointType > & localPoints() const
Return pointField of points in patch.
wordList types() const
Return a list of patch types.
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:371
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
An STL iterator.
Definition: UPtrList.H:185
friend Ostream & operator(Ostream &, const UPtrList< T > &)
Write UPtrList to Ostream.
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:142
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
An STL-conforming iterator.
Definition: HashTable.H:426
void matchGroups(const labelUList &patchIDs, wordList &groups, labelHashSet &nonGroupPatches) const
Match the patches to groups. Returns all the (fully matched) groups.
T clone(const T &t)
Definition: List.H:54
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
dynamicFvMesh & mesh
~polyBoundaryMesh()
Destructor.
void movePoints(const pointField &)
Correct polyBoundaryMesh after moving points.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A class for handling words, derived from string.
Definition: word.H:59
bool checkDefinition(const bool report=false) const
Check boundary definition. Return true if in error.
wordList names() const
Return a list of patch names.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void setCapacity(const label)
Alter the size of the underlying storage.
Definition: DynamicListI.H:130
wordList patchNames(nPatches)
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
List< label > labelList
A List of labels.
Definition: labelList.H:56
wordList physicalTypes() const
Return a list of physical types.
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
An STL-conforming hash table.
Definition: HashTable.H:61
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const List< labelPairList > & nbrEdges() const
Per patch the edges on the neighbouring patch. Is for every external.
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:74
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
label nEdges() const
Return number of edges in patch.
static const char nl
Definition: Ostream.H:260
void updateMesh()
Correct polyBoundaryMesh after topology update.
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
void setSize(const label)
Reset size of List.
Definition: List.C:281
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void clearAddressing()
Clear addressing at this level and at patches.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
Ostream & operator<<(Ostream &, const ensightPart &)
labelList findIndices(const keyType &, const bool usePatchGroups=true) const
Return patch indices for all matches. Optionally matches patchGroups.
Version number type.
Definition: IOstream.H:96
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
messageStream Info
virtual bool writeObject(IOstream::streamFormat, IOstream::versionNumber, IOstream::compressionType, const bool write) const
Write using given format, version and compression.
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
void shuffle(const labelUList &newToOld)
Reorders elements. Ordering does not have to be done in.
Definition: PtrList.C:248
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
label findIndex(const keyType &) const
Return patch index for the first match, return -1 if not found.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void setGroup(const word &groupName, const labelList &patchIDs)
Set/add group with patches.
Namespace for OpenFOAM.
bool isPattern() const
Should be treated as a match rather than a literal string.
Definition: keyTypeI.H:97
const HashTable< labelList, word > & groupPatchIDs() const
Per patch group the patch indices.