pointBoundaryMesh.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-2023 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 "pointBoundaryMesh.H"
30 #include "polyBoundaryMesh.H"
31 #include "facePointPatch.H"
32 #include "pointMesh.H"
33 #include "PstreamBuffers.H"
34 #include "lduSchedule.H"
35 #include "globalMeshData.H"
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 void Foam::pointBoundaryMesh::addPatches(const polyBoundaryMesh& pbm)
40 {
41  // Set boundary patches
42  pointPatchList& patches = *this;
43 
44  patches.resize_null(pbm.size());
45 
46  forAll(patches, patchi)
47  {
48  // NB: needs ptr() to get *pointPatch instead of *facePointPatch
49  patches.set(patchi, facePointPatch::New(pbm[patchi], *this).ptr());
50  }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 Foam::pointBoundaryMesh::pointBoundaryMesh
57 (
58  const pointMesh& m,
59  const polyBoundaryMesh& pbm
60 )
61 :
63  mesh_(m)
64 {
65  addPatches(pbm);
66 }
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
72 (
73  const wordRe& matcher,
74  const bool useGroups
75 ) const
76 {
77  return mesh().boundaryMesh().indices(matcher, useGroups);
78 }
79 
80 
82 (
83  const wordRes& matcher,
84  const bool useGroups
85 ) const
86 {
87  return mesh().boundaryMesh().indices(matcher, useGroups);
88 }
89 
90 
92 (
93  const wordRes& select,
94  const wordRes& ignore,
95  const bool useGroups
96 ) const
97 {
98  return mesh().boundaryMesh().indices(select, ignore, useGroups);
99 }
100 
101 
102 Foam::label Foam::pointBoundaryMesh::findPatchID(const word& patchName) const
103 {
104  return mesh().boundaryMesh().findPatchID(patchName);
105 }
106 
107 
108 void Foam::pointBoundaryMesh::calcGeometry()
109 {
110  PstreamBuffers pBufs(Pstream::defaultCommsType);
111 
112  if
113  (
114  pBufs.commsType() == Pstream::commsTypes::blocking
115  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
116  )
117  {
118  forAll(*this, patchi)
119  {
120  operator[](patchi).initGeometry(pBufs);
121  }
122 
123  pBufs.finishedSends();
124 
125  forAll(*this, patchi)
126  {
127  operator[](patchi).calcGeometry(pBufs);
128  }
129  }
130  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
131  {
132  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
133 
134  // Dummy.
135  pBufs.finishedSends();
136 
137  for (const auto& schedEval : patchSchedule)
138  {
139  const label patchi = schedEval.patch;
140 
141  if (schedEval.init)
142  {
143  operator[](patchi).initGeometry(pBufs);
144  }
145  else
146  {
147  operator[](patchi).calcGeometry(pBufs);
148  }
149  }
150  }
151 }
152 
153 
155 {
156  PstreamBuffers pBufs(Pstream::defaultCommsType);
157 
158  if
159  (
160  pBufs.commsType() == Pstream::commsTypes::blocking
161  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
162  )
163  {
164  forAll(*this, patchi)
165  {
166  operator[](patchi).initMovePoints(pBufs, p);
167  }
168 
169  pBufs.finishedSends();
170 
171  forAll(*this, patchi)
172  {
173  operator[](patchi).movePoints(pBufs, p);
174  }
175  }
176  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
177  {
178  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
179 
180  // Dummy.
181  pBufs.finishedSends();
182 
183  for (const auto& schedEval : patchSchedule)
184  {
185  const label patchi = schedEval.patch;
186 
187  if (schedEval.init)
188  {
189  operator[](patchi).initMovePoints(pBufs, p);
190  }
191  else
192  {
193  operator[](patchi).movePoints(pBufs, p);
194  }
195  }
196  }
197 }
198 
199 
201 {
202  PstreamBuffers pBufs(Pstream::defaultCommsType);
203 
204  if
205  (
206  pBufs.commsType() == Pstream::commsTypes::blocking
207  || pBufs.commsType() == Pstream::commsTypes::nonBlocking
208  )
209  {
210  forAll(*this, patchi)
211  {
212  operator[](patchi).initUpdateMesh(pBufs);
213  }
214 
215  pBufs.finishedSends();
216 
217  forAll(*this, patchi)
218  {
219  operator[](patchi).updateMesh(pBufs);
220  }
221  }
222  else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
223  {
224  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
225 
226  // Dummy.
227  pBufs.finishedSends();
228 
229  for (const auto& schedEval : patchSchedule)
230  {
231  const label patchi = schedEval.patch;
232 
233  if (schedEval.init)
234  {
235  operator[](patchi).initUpdateMesh(pBufs);
236  }
237  else
238  {
239  operator[](patchi).updateMesh(pBufs);
240  }
241  }
242  }
243 }
244 
245 
246 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const polyBoundaryMesh & pbm
PtrList< pointPatch > pointPatchList
Store lists of pointPatch as a PtrList.
Definition: pointPatch.H:50
"blocking" : (MPI_Bsend, MPI_Recv)
void movePoints(const pointField &)
Correct polyBoundaryMesh after moving points.
List< bool > select(const label n, const labelUList &locations)
Construct a selection list of bools (all false) with the given pre-size, subsequently add specified l...
Definition: BitOps.C:134
label findPatchID(const word &patchName) const
Find patch index given a name.
labelList indices(const wordRe &matcher, const bool useGroups) const
Return (sorted) patch indices for all matches.
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:46
void updateMesh()
Correct polyBoundaryMesh after topology update.
static autoPtr< facePointPatch > New(const polyPatch &, const pointBoundaryMesh &)
Return a pointer to a new patch created on freestore from polyPatch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) patch indices for all matches.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
"scheduled" : (MPI_Send, MPI_Recv)
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings...
Definition: wordRe.H:78
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:385
const lduSchedule & patchSchedule() const noexcept
Order in which the patches should be initialised/evaluated corresponding to the schedule.
const polyBoundaryMesh & patches
"nonBlocking" : (MPI_Isend, MPI_Irecv)
volScalarField & p