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 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 Foam::pointBoundaryMesh::pointBoundaryMesh
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 
64 
66 (
67  const keyType& key,
68  const bool useGroups
69 ) const
70 {
71  return mesh()().boundaryMesh().indices(key, useGroups);
72 }
73 
74 
76 {
77  return mesh()().boundaryMesh().findPatchID(patchName);
78 }
79 
80 
81 void Foam::pointBoundaryMesh::calcGeometry()
82 {
84 
85  if
86  (
89  )
90  {
91  forAll(*this, patchi)
92  {
93  operator[](patchi).initGeometry(pBufs);
94  }
95 
96  pBufs.finishedSends();
97 
98  forAll(*this, patchi)
99  {
100  operator[](patchi).calcGeometry(pBufs);
101  }
102  }
104  {
105  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
106 
107  // Dummy.
108  pBufs.finishedSends();
109 
110  forAll(patchSchedule, patchEvali)
111  {
112  label patchi = patchSchedule[patchEvali].patch;
113 
114  if (patchSchedule[patchEvali].init)
115  {
116  operator[](patchi).initGeometry(pBufs);
117  }
118  else
119  {
120  operator[](patchi).calcGeometry(pBufs);
121  }
122  }
123  }
124 }
125 
126 
128 {
130 
131  if
132  (
135  )
136  {
137  forAll(*this, patchi)
138  {
139  operator[](patchi).initMovePoints(pBufs, p);
140  }
141 
142  pBufs.finishedSends();
143 
144  forAll(*this, patchi)
145  {
146  operator[](patchi).movePoints(pBufs, p);
147  }
148  }
150  {
151  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
152 
153  // Dummy.
154  pBufs.finishedSends();
155 
156  forAll(patchSchedule, patchEvali)
157  {
158  label patchi = patchSchedule[patchEvali].patch;
159 
160  if (patchSchedule[patchEvali].init)
161  {
162  operator[](patchi).initMovePoints(pBufs, p);
163  }
164  else
165  {
166  operator[](patchi).movePoints(pBufs, p);
167  }
168  }
169  }
170 }
171 
172 
174 {
176 
177  if
178  (
181  )
182  {
183  forAll(*this, patchi)
184  {
185  operator[](patchi).initUpdateMesh(pBufs);
186  }
187 
188  pBufs.finishedSends();
189 
190  forAll(*this, patchi)
191  {
192  operator[](patchi).updateMesh(pBufs);
193  }
194  }
196  {
197  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
198 
199  // Dummy.
200  pBufs.finishedSends();
201 
202  forAll(patchSchedule, patchEvali)
203  {
204  label patchi = patchSchedule[patchEvali].patch;
205 
206  if (patchSchedule[patchEvali].init)
207  {
208  operator[](patchi).initUpdateMesh(pBufs);
209  }
210  else
211  {
212  operator[](patchi).updateMesh(pBufs);
213  }
214  }
215  }
216 }
217 
218 
219 // ************************************************************************* //
Foam::UPstream::commsTypes::blocking
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::UPstream::commsTypes::nonBlocking
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::lduSchedule
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:76
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
globalMeshData.H
Foam::pointBoundaryMesh::movePoints
void movePoints(const pointField &)
Correct polyBoundaryMesh after moving points.
Definition: pointBoundaryMesh.C:127
Foam::pointPatchList
PtrList< pointPatch > pointPatchList
container classes for pointPatch
Definition: pointPatchList.H:47
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
Foam::pointBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: pointBoundaryMesh.C:75
Foam::UPstream::defaultCommsType
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:273
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
pointBoundaryMesh.H
Foam::pointBoundaryMesh::updateMesh
void updateMesh()
Correct polyBoundaryMesh after topology update.
Definition: pointBoundaryMesh.C:173
Foam::keyType
A class for handling keywords in dictionaries.
Definition: keyType.H:60
Foam::label
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:62
Foam::Field< vector >
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::pointBoundaryMesh::mesh
const pointMesh & mesh() const
Return the mesh reference.
Definition: pointBoundaryMesh.H:96
lduSchedule.H
Foam::PstreamBuffers::finishedSends
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
Definition: PstreamBuffers.C:80
Foam::UPstream::commsTypes::scheduled
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::ProcessorTopology::patchSchedule
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
Definition: ProcessorTopology.H:101
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:50
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: PtrListI.H:143
Foam::pointBoundaryMesh::indices
labelList indices(const keyType &key, const bool useGroups) const
Find patch indices given a name.
Definition: pointBoundaryMesh.C:66
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
facePointPatch.H
Foam::List< label >
PstreamBuffers.H
polyBoundaryMesh.H
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1241
pointMesh.H