All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
pointBoundaryMesh.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-2023 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 "pointBoundaryMesh.H"
27 #include "polyBoundaryMesh.H"
28 #include "processorPolyPatch.H"
29 #include "facePointPatch.H"
30 #include "pointMesh.H"
31 #include "PstreamBuffers.H"
32 #include "lduSchedule.H"
33 #include "globalMeshData.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const pointMesh& m,
40  const polyBoundaryMesh& pbm
41 )
42 :
43  pointPatchList(pbm.size()),
44  mesh_(m)
45 {
46  // Set boundary patches
47  pointPatchList& Patches = *this;
48 
49  forAll(Patches, patchi)
50  {
51  Patches.set
52  (
53  patchi,
54  facePointPatch::New(pbm[patchi], *this).ptr()
55  );
56  }
57  // reset(pbm);
58 }
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
64 {
65  return mesh()().boundaryMesh().findIndex(patchName);
66 }
67 
68 
70 (
71  const wordRe& key,
72  const bool usePatchGroups
73 ) const
74 {
75  return mesh()().boundaryMesh().findIndices(key, usePatchGroups);
76 }
77 
78 
79 void Foam::pointBoundaryMesh::calcGeometry()
80 {
82 
83  if
84  (
87  )
88  {
89  forAll(*this, patchi)
90  {
91  operator[](patchi).initCalcGeometry(pBufs);
92  }
93 
94  pBufs.finishedSends();
95 
96  forAll(*this, patchi)
97  {
98  operator[](patchi).calcGeometry(pBufs);
99  }
100  }
102  {
103  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
104 
105  // Dummy.
106  pBufs.finishedSends();
107 
108  forAll(patchSchedule, patchEvali)
109  {
110  label patchi = patchSchedule[patchEvali].patch;
111 
112  if (patchSchedule[patchEvali].init)
113  {
114  operator[](patchi).initCalcGeometry(pBufs);
115  }
116  else
117  {
118  operator[](patchi).calcGeometry(pBufs);
119  }
120  }
121  }
122 }
123 
124 
126 {
128 
129  if
130  (
133  )
134  {
135  forAll(*this, patchi)
136  {
137  operator[](patchi).initMovePoints(pBufs, p);
138  }
139 
140  pBufs.finishedSends();
141 
142  forAll(*this, patchi)
143  {
144  operator[](patchi).movePoints(pBufs, p);
145  }
146  }
148  {
149  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
150 
151  // Dummy.
152  pBufs.finishedSends();
153 
154  forAll(patchSchedule, patchEvali)
155  {
156  label patchi = patchSchedule[patchEvali].patch;
157 
158  if (patchSchedule[patchEvali].init)
159  {
160  operator[](patchi).initMovePoints(pBufs, p);
161  }
162  else
163  {
164  operator[](patchi).movePoints(pBufs, p);
165  }
166  }
167  }
168 }
169 
170 
172 {
174 
175  if
176  (
179  )
180  {
181  forAll(*this, patchi)
182  {
183  operator[](patchi).initTopoChange(pBufs);
184  }
185 
186  pBufs.finishedSends();
187 
188  forAll(*this, patchi)
189  {
190  operator[](patchi).topoChange(pBufs);
191  }
192  }
194  {
195  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
196 
197  // Dummy.
198  pBufs.finishedSends();
199 
200  forAll(patchSchedule, patchEvali)
201  {
202  label patchi = patchSchedule[patchEvali].patch;
203 
204  if (patchSchedule[patchEvali].init)
205  {
206  operator[](patchi).initTopoChange(pBufs);
207  }
208  else
209  {
210  operator[](patchi).topoChange(pBufs);
211  }
212  }
213  }
214 }
215 
216 
218 {
219  const polyBoundaryMesh& boundaryMesh = mesh()().boundaryMesh();
220  pointPatchList& Patches = *this;
221 
222  // Reset the number of patches in case the decomposition changed
223  Patches.setSize(boundaryMesh.size());
224 
225  forAll(boundaryMesh, patchi)
226  {
227  // Construct new processor patches in case the decomposition changed
228  if (isA<processorPolyPatch>(boundaryMesh[patchi]))
229  {
230  Patches.set
231  (
232  patchi,
233  facePointPatch::New(boundaryMesh[patchi], *this).ptr()
234  );
235  }
236  }
237 }
238 
239 
241 (
242  const labelUList& newToOld,
243  const bool validBoundary
244 )
245 {
246  pointPatchList::shuffle(newToOld);
247  if (validBoundary)
248  {
249  topoChange();
250  }
251 }
252 
253 
254 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
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 shuffle(const labelUList &newToOld)
Reorders elements. Ordering does not have to be done in.
Definition: PtrList.C:272
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static autoPtr< facePointPatch > New(const polyPatch &, const pointBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
void topoChange()
Correct pointBoundaryMesh after topology update.
label findIndex(const word &patchName) const
Find patch index given a name.
labelList findIndices(const wordRe &, const bool useGroups) const
Find patch indices given a name.
void movePoints(const pointField &)
Correct pointBoundaryMesh after moving points.
pointBoundaryMesh(const pointMesh &, const polyBoundaryMesh &)
Construct from polyBoundaryMesh.
void shuffle(const labelUList &newToOld, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
void reset()
Reset pointBoundaryMesh with respect to the updated polyBoundaryMesh.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
Foam::polyBoundaryMesh.
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:77
A class for handling words, derived from string.
Definition: word.H:62
label patchi
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
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:74
volScalarField & p