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 | Copyright (C) 2011-2013 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  // Set boundary patches
46  pointPatchList& Patches = *this;
47 
48  forAll(Patches, patchI)
49  {
50  Patches.set
51  (
52  patchI,
53  facePointPatch::New(basicBdry[patchI], *this).ptr()
54  );
55  }
56 }
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
62 {
63  return mesh()().boundaryMesh().findPatchID(patchName);
64 }
65 
66 
68 (
69  const keyType& key,
70  const bool usePatchGroups
71 ) const
72 {
73  return mesh()().boundaryMesh().findIndices(key, usePatchGroups);
74 }
75 
76 
78 {
80 
81  if
82  (
85  )
86  {
87  forAll(*this, patchi)
88  {
89  operator[](patchi).initGeometry(pBufs);
90  }
91 
92  pBufs.finishedSends();
93 
94  forAll(*this, patchi)
95  {
96  operator[](patchi).calcGeometry(pBufs);
97  }
98  }
100  {
101  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
102 
103  // Dummy.
104  pBufs.finishedSends();
105 
106  forAll(patchSchedule, patchEvali)
107  {
108  label patchi = patchSchedule[patchEvali].patch;
109 
110  if (patchSchedule[patchEvali].init)
111  {
112  operator[](patchi).initGeometry(pBufs);
113  }
114  else
115  {
116  operator[](patchi).calcGeometry(pBufs);
117  }
118  }
119  }
120 }
121 
122 
124 {
126 
127  if
128  (
131  )
132  {
133  forAll(*this, patchi)
134  {
135  operator[](patchi).initMovePoints(pBufs, p);
136  }
137 
138  pBufs.finishedSends();
139 
140  forAll(*this, patchi)
141  {
142  operator[](patchi).movePoints(pBufs, p);
143  }
144  }
146  {
147  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
148 
149  // Dummy.
150  pBufs.finishedSends();
151 
152  forAll(patchSchedule, patchEvali)
153  {
154  label patchi = patchSchedule[patchEvali].patch;
155 
156  if (patchSchedule[patchEvali].init)
157  {
158  operator[](patchi).initMovePoints(pBufs, p);
159  }
160  else
161  {
162  operator[](patchi).movePoints(pBufs, p);
163  }
164  }
165  }
166 }
167 
168 
170 {
172 
173  if
174  (
177  )
178  {
179  forAll(*this, patchi)
180  {
181  operator[](patchi).initUpdateMesh(pBufs);
182  }
183 
184  pBufs.finishedSends();
185 
186  forAll(*this, patchi)
187  {
188  operator[](patchi).updateMesh(pBufs);
189  }
190  }
192  {
193  const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
194 
195  // Dummy.
196  pBufs.finishedSends();
197 
198  forAll(patchSchedule, patchEvali)
199  {
200  label patchi = patchSchedule[patchEvali].patch;
201 
202  if (patchSchedule[patchEvali].init)
203  {
204  operator[](patchi).initUpdateMesh(pBufs);
205  }
206  else
207  {
208  operator[](patchi).updateMesh(pBufs);
209  }
210  }
211  }
212 }
213 
214 
215 // ************************************************************************* //
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::UPstream::scheduled
@ scheduled
Definition: UPstream.H:67
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
globalMeshData.H
Foam::pointBoundaryMesh::movePoints
void movePoints(const pointField &)
Correct polyBoundaryMesh after moving points.
Definition: pointBoundaryMesh.C:123
Foam::pointPatchList
PtrList< pointPatch > pointPatchList
container classes for pointPatch
Definition: pointPatchList.H:45
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:85
Foam::pointBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: pointBoundaryMesh.C:61
Foam::UPstream::defaultCommsType
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:252
Foam::pointBoundaryMesh::findIndices
labelList findIndices(const keyType &, const bool useGroups) const
Find patch indices given a name.
Definition: pointBoundaryMesh.C:68
pointBoundaryMesh.H
Foam::pointBoundaryMesh::updateMesh
void updateMesh()
Correct polyBoundaryMesh after topology update.
Definition: pointBoundaryMesh.C:169
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
Foam::keyType
A class for handling keywords in dictionaries.
Definition: keyType.H:56
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::UPstream::nonBlocking
@ nonBlocking
Definition: UPstream.H:68
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::pointBoundaryMesh::pointBoundaryMesh
pointBoundaryMesh(const pointBoundaryMesh &)
Disallow default bitwise copy construct.
Foam::pointBoundaryMesh::mesh
const pointMesh & mesh() const
Return the mesh reference.
Definition: pointBoundaryMesh.H:93
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:82
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::ProcessorTopology::patchSchedule
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
Definition: ProcessorTopology.H:99
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::pointBoundaryMesh::calcGeometry
void calcGeometry()
Calculate the geometry for the patches (transformation tensors etc.)
Definition: pointBoundaryMesh.C:77
facePointPatch.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
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:59
Foam::boundaryMesh::findPatchID
label findPatchID(const polyPatchList &, const word &) const
Get index of polypatch by name.
Definition: boundaryMesh.C:254
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
pointMesh.H