geometryPatchesGather.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) 2019 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // OpenFOAM includes
29 #include "geometryPatches.H"
30 #include "fvMesh.H"
31 #include "volFields.H"
32 #include "emptyPolyPatch.H"
33 #include "processorPolyPatch.H"
34 #include "foamVtkTools.H"
35 #include "runTimePostProcessing.H"
36 
37 // VTK includes
38 #include "vtkCellArray.h"
39 #include "vtkCellData.h"
40 #include "vtkCompositeDataSet.h"
41 #include "vtkMultiPieceDataSet.h"
42 #include "vtkPointData.h"
43 #include "vtkPolyData.h"
44 #include "vtkSmartPointer.h"
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
50 (
51  const labelListList& patchIds
52 ) const
53 {
54  const polyBoundaryMesh& patches = parent().mesh().boundaryMesh();
55 
56  label nPieces = 0;
57  for (const labelList& ids : patchIds)
58  {
59  nPieces += ids.size();
60  }
61 
63  multiPiece->SetNumberOfPieces(nPieces);
64 
65  label pieceId = 0;
66 
67  if (!needsCollective())
68  {
69  // Simple case
70 
71  for (int proci=0; proci < Pstream::myProcNo(); ++proci)
72  {
73  pieceId += patchIds[proci].size();
74  }
75 
76  for (const label patchId : patchIds[Pstream::myProcNo()])
77  {
78  const polyPatch& pp = patches[patchId];
79 
80  multiPiece->SetPiece
81  (
82  pieceId,
84  );
85 
86  ++pieceId;
87  }
88  }
89  else if (Pstream::master())
90  {
91  // Gather pieces on master
92 
93  // Add myself
94  for (const label patchId : patchIds[Pstream::myProcNo()])
95  {
96  const polyPatch& pp = patches[patchId];
97 
98  multiPiece->SetPiece
99  (
100  pieceId,
102  );
103 
104  ++pieceId;
105  }
106 
107  // Receive surfaces
109  faceList faces;
110 
111  for
112  (
113  int slave=Pstream::firstSlave();
114  slave<=Pstream::lastSlave();
115  ++slave
116  )
117  {
118  const label nSlavePatches = patchIds[slave].size();
119 
120  if (!nSlavePatches)
121  {
122  continue;
123  }
124 
125  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
126 
127  for (label recvi=0; recvi < nSlavePatches; ++recvi)
128  {
129  points.clear();
130  faces.clear();
131 
132  fromSlave
133  >> points >> faces;
134 
135  multiPiece->SetPiece
136  (
137  pieceId,
139  );
140 
141  ++pieceId;
142  }
143  }
144  }
145  else
146  {
147  // Slave - send surfaces
148  const labelList& slavePatchIds = patchIds[Pstream::myProcNo()];
149 
150  if (slavePatchIds.size())
151  {
152  OPstream toMaster
153  (
154  Pstream::commsTypes::scheduled,
155  Pstream::masterNo()
156  );
157 
158  for (const label patchId : patchIds[Pstream::myProcNo()])
159  {
160  const polyPatch& pp = patches[patchId];
161 
162  toMaster
163  << pp.localPoints() << pp.localFaces();
164  }
165  }
166  }
167 
168  return multiPiece;
169 }
170 
171 
174 (
175  const labelListList& patchIds
176 ) const
177 {
178  const polyBoundaryMesh& patches = parent().mesh().boundaryMesh();
179 
180  label nPieces = 0;
181  for (const labelList& ids : patchIds)
182  {
183  nPieces += ids.size();
184  }
185 
186  auto multiPiece = vtkSmartPointer<vtkMultiPieceDataSet>::New();
187  multiPiece->SetNumberOfPieces(nPieces);
188 
189  label pieceId = 0;
190 
191  if (!needsCollective())
192  {
193  // Simple case
194 
195  for (int proci=0; proci < Pstream::myProcNo(); ++proci)
196  {
197  pieceId += patchIds[proci].size();
198  }
199 
200  for (const label patchId : patchIds[Pstream::myProcNo()])
201  {
202  const polyPatch& pp = patches[patchId];
203 
204  auto geom = vtkSmartPointer<vtkPolyData>::New();
205 
206  geom->SetPoints(Foam::vtk::Tools::Patch::faceCentres(pp));
207  geom->SetVerts(Foam::vtk::Tools::identityVertices(pp.size()));
208 
209  multiPiece->SetPiece(pieceId, geom);
210 
211  ++pieceId;
212  }
213  }
214  else if (Pstream::master())
215  {
216  // Gather pieces (face centres) on master
217 
218  // Add myself
219  for (const label patchId : patchIds[Pstream::myProcNo()])
220  {
221  const polyPatch& pp = patches[patchId];
222 
223  auto geom = vtkSmartPointer<vtkPolyData>::New();
224 
225  geom->SetPoints(Foam::vtk::Tools::Patch::faceCentres(pp));
226  geom->SetVerts(Foam::vtk::Tools::identityVertices(pp.size()));
227 
228  multiPiece->SetPiece(pieceId, geom);
229 
230  ++pieceId;
231  }
232 
233  // Receive points
235 
236  for
237  (
238  int slave=Pstream::firstSlave();
239  slave<=Pstream::lastSlave();
240  ++slave
241  )
242  {
243  const label nSlavePatches = patchIds[slave].size();
244 
245  if (!nSlavePatches)
246  {
247  continue;
248  }
249 
250  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
251 
252  for (label recvi=0; recvi < nSlavePatches; ++recvi)
253  {
254  points.clear();
255 
256  fromSlave >> points;
257 
258  multiPiece->SetPiece
259  (
260  pieceId,
262  );
263 
264  ++pieceId;
265  }
266  }
267  }
268  else
269  {
270  // Slave - send face centres
271 
272  const labelList& slavePatchIds = patchIds[Pstream::myProcNo()];
273 
274  if (slavePatchIds.size())
275  {
276  OPstream toMaster
277  (
278  Pstream::commsTypes::scheduled,
279  Pstream::masterNo()
280  );
281 
282  for (const label patchId : patchIds[Pstream::myProcNo()])
283  {
284  const polyPatch& pp = patches[patchId];
285 
286  toMaster << pp.faceCentres();
287  }
288  }
289  }
290 
291  return multiPiece;
292 }
293 
294 
295 // ************************************************************************* //
volFields.H
Foam::vtk::Tools::Vertices
vtkSmartPointer< vtkPolyData > Vertices(const UList< point > &pts)
Return vtkPolyData of vertices for each point.
Definition: foamVtkToolsI.H:104
Foam::functionObjects::runTimePostPro::geometryPatches::gatherPatchPieces
vtkSmartPointer< vtkMultiPieceDataSet > gatherPatchPieces(const labelListList &patchIds) const
Definition: geometryPatchesGather.C:50
foamVtkTools.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::vtk::Tools::Patch::mesh
static vtkSmartPointer< vtkPolyData > mesh(const PatchType &p)
Convert patch points/faces to vtkPolyData.
Definition: foamVtkToolsTemplates.C:92
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:458
Foam::vtk::Tools::Patch::faceCentres
static vtkSmartPointer< vtkPoints > faceCentres(const PatchType &p)
Return patch face centres as vtkPoints.
Definition: foamVtkToolsTemplates.C:156
Foam::functionObjects::runTimePostPro::geometryPatches::gatherPatchFaceCentres
vtkSmartPointer< vtkMultiPieceDataSet > gatherPatchFaceCentres(const labelListList &patchIds) const
Definition: geometryPatchesGather.C:174
vtkSmartPointer
Definition: runTimePostProcessing.H:148
geometryPatches.H
patchIds
labelList patchIds
Definition: convertProcessorPatches.H:67
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::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::vtk::Tools::identityVertices
vtkSmartPointer< vtkCellArray > identityVertices(const label size)
An identity list of VTK_VERTEX.
Definition: foamVtkToolsI.H:149
processorPolyPatch.H
fvMesh.H
emptyPolyPatch.H
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
runTimePostProcessing.H
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:311
Foam::List< labelList >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:398
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
patchId
label patchId(-1)