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-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 "pointBoundaryMesh.H"
27 #include "polyBoundaryMesh.H"
28 #include "facePointPatch.H"
29 #include "pointMesh.H"
30 #include "PstreamBuffers.H"
31 #include "lduSchedule.H"
32 #include "globalMeshData.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
37 (
38  const pointMesh& m,
39  const polyBoundaryMesh& basicBdry
40 )
41 :
42  pointPatchList(basicBdry.size()),
43  mesh_(m)
44 {
45  reset(basicBdry);
46 }
47 
48 
49 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
50 
52 {
53  return mesh()().boundaryMesh().findPatchID(patchName);
54 }
55 
56 
58 (
59  const wordRe& key,
60  const bool usePatchGroups
61 ) const
62 {
63  return mesh()().boundaryMesh().findIndices(key, usePatchGroups);
64 }
65 
66 
67 void Foam::pointBoundaryMesh::calcGeometry()
68 {
70 
71  if
72  (
75  )
76  {
77  forAll(*this, patchi)
78  {
79  operator[](patchi).initCalcGeometry(pBufs);
80  }
81 
82  pBufs.finishedSends();
83 
84  forAll(*this, patchi)
85  {
86  operator[](patchi).calcGeometry(pBufs);
87  }
88  }
90  {
91  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
92 
93  // Dummy.
94  pBufs.finishedSends();
95 
96  forAll(patchSchedule, patchEvali)
97  {
98  label patchi = patchSchedule[patchEvali].patch;
99 
100  if (patchSchedule[patchEvali].init)
101  {
102  operator[](patchi).initCalcGeometry(pBufs);
103  }
104  else
105  {
106  operator[](patchi).calcGeometry(pBufs);
107  }
108  }
109  }
110 }
111 
112 
114 {
116 
117  if
118  (
121  )
122  {
123  forAll(*this, patchi)
124  {
125  operator[](patchi).initMovePoints(pBufs, p);
126  }
127 
128  pBufs.finishedSends();
129 
130  forAll(*this, patchi)
131  {
132  operator[](patchi).movePoints(pBufs, p);
133  }
134  }
136  {
137  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
138 
139  // Dummy.
140  pBufs.finishedSends();
141 
142  forAll(patchSchedule, patchEvali)
143  {
144  label patchi = patchSchedule[patchEvali].patch;
145 
146  if (patchSchedule[patchEvali].init)
147  {
148  operator[](patchi).initMovePoints(pBufs, p);
149  }
150  else
151  {
152  operator[](patchi).movePoints(pBufs, p);
153  }
154  }
155  }
156 }
157 
158 
160 {
162 
163  if
164  (
167  )
168  {
169  forAll(*this, patchi)
170  {
171  operator[](patchi).initUpdateMesh(pBufs);
172  }
173 
174  pBufs.finishedSends();
175 
176  forAll(*this, patchi)
177  {
178  operator[](patchi).updateMesh(pBufs);
179  }
180  }
182  {
183  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
184 
185  // Dummy.
186  pBufs.finishedSends();
187 
188  forAll(patchSchedule, patchEvali)
189  {
190  label patchi = patchSchedule[patchEvali].patch;
191 
192  if (patchSchedule[patchEvali].init)
193  {
194  operator[](patchi).initUpdateMesh(pBufs);
195  }
196  else
197  {
198  operator[](patchi).updateMesh(pBufs);
199  }
200  }
201  }
202 }
203 
204 
206 {
207  // Set boundary patches
208  pointPatchList& Patches = *this;
209 
210  forAll(Patches, patchi)
211  {
212  Patches.set
213  (
214  patchi,
215  facePointPatch::New(basicBdry[patchi], *this).ptr()
216  );
217  }
218 }
219 
220 
222 (
223  const labelUList& newToOld,
224  const bool validBoundary
225 )
226 {
227  pointPatchList::shuffle(newToOld);
228  if (validBoundary)
229  {
230  updateMesh();
231  }
232 }
233 
234 
235 // ************************************************************************* //
labelList findIndices(const wordRe &, const bool useGroups) const
Find patch indices given a name.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 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.
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
void movePoints(const pointField &)
Correct pointBoundaryMesh after moving points.
label findPatchID(const word &patchName) const
Find patch index given a name.
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
void shuffle(const labelUList &newToOld, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
void updateMesh()
Correct pointBoundaryMesh after topology update.
static autoPtr< facePointPatch > New(const polyPatch &, const pointBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
const globalMeshData & globalData() const
Return parallel info.
Definition: pointMesh.H:109
A class for handling words, derived from string.
Definition: word.H:59
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:74
void reset(const polyBoundaryMesh &)
Create pointBoundaryMesh from polyBoundaryMesh.
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
Foam::polyBoundaryMesh.
const pointMesh & mesh() const
Return the mesh reference.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
pointBoundaryMesh(const pointMesh &, const polyBoundaryMesh &)
Construct from polyBoundaryMesh.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
label patchi
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 shuffle(const labelUList &newToOld)
Reorders elements. Ordering does not have to be done in.
Definition: PtrList.C:248
volScalarField & p
PtrList< pointPatch > pointPatchList
container classes for pointPatch