ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / OpenFOAM / meshes / pointMesh / pointBoundaryMesh / pointBoundaryMesh.C
blob113a9719a025f3c1eefe479063a2b95e3863384e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
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"
34 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
36 Foam::pointBoundaryMesh::pointBoundaryMesh
38     const pointMesh& m,
39     const polyBoundaryMesh& basicBdry
42     pointPatchList(basicBdry.size()),
43     mesh_(m)
45     // Set boundary patches
46     pointPatchList& Patches = *this;
48     forAll(Patches, patchI)
49     {
50         Patches.set
51         (
52             patchI,
53             facePointPatch::New(basicBdry[patchI], *this).ptr()
54         );
55     }
59 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
61 void Foam::pointBoundaryMesh::calcGeometry()
63     PstreamBuffers pBufs(Pstream::defaultCommsType);
65     if
66     (
67         Pstream::defaultCommsType == Pstream::blocking
68      || Pstream::defaultCommsType == Pstream::nonBlocking
69     )
70     {
71         forAll(*this, patchi)
72         {
73             operator[](patchi).initGeometry(pBufs);
74         }
76         pBufs.finishedSends();
78         forAll(*this, patchi)
79         {
80             operator[](patchi).calcGeometry(pBufs);
81         }
82     }
83     else if (Pstream::defaultCommsType == Pstream::scheduled)
84     {
85         const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
87         // Dummy.
88         pBufs.finishedSends();
90         forAll(patchSchedule, patchEvali)
91         {
92             label patchi = patchSchedule[patchEvali].patch;
94             if (patchSchedule[patchEvali].init)
95             {
96                 operator[](patchi).initGeometry(pBufs);
97             }
98             else
99             {
100                 operator[](patchi).calcGeometry(pBufs);
101             }
102         }
103     }
107 void Foam::pointBoundaryMesh::movePoints(const pointField& p)
109     PstreamBuffers pBufs(Pstream::defaultCommsType);
111     if
112     (
113         Pstream::defaultCommsType == Pstream::blocking
114      || Pstream::defaultCommsType == Pstream::nonBlocking
115     )
116     {
117         forAll(*this, patchi)
118         {
119             operator[](patchi).initMovePoints(pBufs, p);
120         }
122         pBufs.finishedSends();
124         forAll(*this, patchi)
125         {
126             operator[](patchi).movePoints(pBufs, p);
127         }
128     }
129     else if (Pstream::defaultCommsType == Pstream::scheduled)
130     {
131         const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
133         // Dummy.
134         pBufs.finishedSends();
136         forAll(patchSchedule, patchEvali)
137         {
138             label patchi = patchSchedule[patchEvali].patch;
140             if (patchSchedule[patchEvali].init)
141             {
142                 operator[](patchi).initMovePoints(pBufs, p);
143             }
144             else
145             {
146                 operator[](patchi).movePoints(pBufs, p);
147             }
148         }
149     }
153 void Foam::pointBoundaryMesh::updateMesh()
155     PstreamBuffers pBufs(Pstream::defaultCommsType);
157     if
158     (
159         Pstream::defaultCommsType == Pstream::blocking
160      || Pstream::defaultCommsType == Pstream::nonBlocking
161     )
162     {
163         forAll(*this, patchi)
164         {
165             operator[](patchi).initUpdateMesh(pBufs);
166         }
168         pBufs.finishedSends();
170         forAll(*this, patchi)
171         {
172             operator[](patchi).updateMesh(pBufs);
173         }
174     }
175     else if (Pstream::defaultCommsType == Pstream::scheduled)
176     {
177         const lduSchedule& patchSchedule = mesh().globalData().patchSchedule();
179         // Dummy.
180         pBufs.finishedSends();
182         forAll(patchSchedule, patchEvali)
183         {
184             label patchi = patchSchedule[patchEvali].patch;
186             if (patchSchedule[patchEvali].init)
187             {
188                 operator[](patchi).initUpdateMesh(pBufs);
189             }
190             else
191             {
192                 operator[](patchi).updateMesh(pBufs);
193             }
194         }
195     }
199 // ************************************************************************* //