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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
40(
41 const pointMesh& m,
42 const polyBoundaryMesh& basicBdry
43)
44:
45 pointPatchList(basicBdry.size()),
46 mesh_(m)
47{
48 // Set boundary patches
49 pointPatchList& Patches = *this;
50
51 forAll(Patches, patchi)
52 {
53 Patches.set
54 (
55 patchi,
56 facePointPatch::New(basicBdry[patchi], *this).ptr()
57 );
58 }
59}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
65(
66 const wordRe& matcher,
67 const bool useGroups
68) const
69{
70 return mesh()().boundaryMesh().indices(matcher, useGroups);
71}
72
73
75(
76 const wordRes& matcher,
77 const bool useGroups
78) const
79{
80 return mesh()().boundaryMesh().indices(matcher, useGroups);
81}
82
83
84Foam::label Foam::pointBoundaryMesh::findPatchID(const word& patchName) const
85{
86 return mesh()().boundaryMesh().findPatchID(patchName);
87}
88
89
90void Foam::pointBoundaryMesh::calcGeometry()
91{
93
94 if
95 (
96 pBufs.commsType() == Pstream::commsTypes::blocking
97 || pBufs.commsType() == Pstream::commsTypes::nonBlocking
98 )
99 {
100 forAll(*this, patchi)
101 {
102 operator[](patchi).initGeometry(pBufs);
103 }
104
105 pBufs.finishedSends();
106
107 forAll(*this, patchi)
108 {
109 operator[](patchi).calcGeometry(pBufs);
110 }
111 }
112 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
113 {
114 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
115
116 // Dummy.
117 pBufs.finishedSends();
118
119 for (const auto& schedEval : patchSchedule)
120 {
121 const label patchi = schedEval.patch;
122
123 if (schedEval.init)
124 {
125 operator[](patchi).initGeometry(pBufs);
126 }
127 else
128 {
129 operator[](patchi).calcGeometry(pBufs);
130 }
131 }
132 }
133}
134
135
137{
139
140 if
141 (
144 )
145 {
146 forAll(*this, patchi)
147 {
148 operator[](patchi).initMovePoints(pBufs, p);
149 }
150
151 pBufs.finishedSends();
152
153 forAll(*this, patchi)
154 {
155 operator[](patchi).movePoints(pBufs, p);
156 }
157 }
158 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
159 {
160 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
161
162 // Dummy.
163 pBufs.finishedSends();
164
165 for (const auto& schedEval : patchSchedule)
166 {
167 const label patchi = schedEval.patch;
168
169 if (schedEval.init)
170 {
171 operator[](patchi).initMovePoints(pBufs, p);
172 }
173 else
174 {
175 operator[](patchi).movePoints(pBufs, p);
176 }
177 }
178 }
179}
180
181
183{
185
186 if
187 (
190 )
191 {
192 forAll(*this, patchi)
193 {
194 operator[](patchi).initUpdateMesh(pBufs);
195 }
196
197 pBufs.finishedSends();
198
199 forAll(*this, patchi)
200 {
201 operator[](patchi).updateMesh(pBufs);
202 }
203 }
204 else if (pBufs.commsType() == Pstream::commsTypes::scheduled)
205 {
206 const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
207
208 // Dummy.
209 pBufs.finishedSends();
210
211 for (const auto& schedEval : patchSchedule)
212 {
213 const label patchi = schedEval.patch;
214
215 if (schedEval.init)
216 {
217 operator[](patchi).initUpdateMesh(pBufs);
218 }
219 else
220 {
221 operator[](patchi).updateMesh(pBufs);
222 }
223 }
224 }
225}
226
227
228// ************************************************************************* //
Buffers for inter-processor communications streams (UOPstream, UIPstream).
UPstream::commsTypes commsType() const noexcept
The communications type of the stream.
void finishedSends(const bool wait=true)
Mark sends as done.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
const labelUList & indices() const noexcept
Return the list of sorted indices (updated every sort).
Definition: SortListI.H:54
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
@ nonBlocking
"nonBlocking"
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
const lduSchedule & patchSchedule() const noexcept
void movePoints()
Update for new mesh geometry.
Foam::pointBoundaryMesh.
label findPatchID(const word &patchName) const
Find patch index given a name.
void updateMesh()
Correct polyBoundaryMesh after topology update.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:83
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333