surfaceGather.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 "surface.H"
30 #include "runTimePostProcessing.H"
31 
32 #include "foamVtkTools.H"
33 #include "polySurface.H"
34 
35 // VTK includes
36 #include "vtkMultiPieceDataSet.h"
37 #include "vtkPolyData.h"
38 #include "vtkSmartPointer.h"
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
44 (
45  const polySurface* surf
46 ) const
47 {
49  multiPiece->SetNumberOfPieces(Pstream::nProcs());
50 
51  if (!needsCollective())
52  {
53  // Simple case (serial-serial, parallel-parallel)
54 
55  if (surf)
56  {
57  multiPiece->SetPiece
58  (
59  Pstream::myProcNo(),
61  );
62  }
63  }
64  else if (Pstream::master())
65  {
66  // Gather pieces on master
67 
68  if (surf)
69  {
70  // Add myself
71 
72  multiPiece->SetPiece
73  (
74  Pstream::myProcNo(),
76  );
77  }
78 
79  // Receive surfaces
81  faceList faces;
82 
83  for
84  (
85  int slave=Pstream::firstSlave();
86  slave<=Pstream::lastSlave();
87  ++slave
88  )
89  {
90  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
91 
92  points.clear();
93  faces.clear();
94 
95  fromSlave >> points >> faces;
96 
97  if (points.size())
98  {
99  multiPiece->SetPiece
100  (
101  slave,
103  );
104  }
105  }
106  }
107  else
108  {
109  // Slave - send surfaces
110 
111  OPstream toMaster
112  (
113  Pstream::commsTypes::scheduled,
114  Pstream::masterNo()
115  );
116 
117  if (surf)
118  {
119  toMaster
120  << surf->points() << surf->faces();
121  }
122  else
123  {
124  toMaster
125  << pointField() << faceList();
126  }
127  }
128 
129  return multiPiece;
130 }
131 
132 
135 (
136  const polySurface* surf
137 ) const
138 {
139  auto multiPiece = vtkSmartPointer<vtkMultiPieceDataSet>::New();
140  multiPiece->SetNumberOfPieces(Pstream::nProcs());
141 
142  if (!needsCollective())
143  {
144  // Simple case
145 
146  if (surf)
147  {
148  auto dataset = vtkSmartPointer<vtkPolyData>::New();
149 
150  auto geom = vtkSmartPointer<vtkPolyData>::New();
151 
152  geom->SetPoints(Foam::vtk::Tools::Patch::faceCentres(*surf));
153  geom->SetVerts(Foam::vtk::Tools::identityVertices(surf->nFaces()));
154 
155  multiPiece->SetPiece(Pstream::myProcNo(), geom);
156  }
157  }
158  else if (Pstream::master())
159  {
160  // Gather pieces (face centres) on master
161 
162  if (surf)
163  {
164  // Add myself
165 
166  auto geom = vtkSmartPointer<vtkPolyData>::New();
167 
168  geom->SetPoints(Foam::vtk::Tools::Patch::faceCentres(*surf));
169  geom->SetVerts(Foam::vtk::Tools::identityVertices(surf->nFaces()));
170 
171  multiPiece->SetPiece(Pstream::myProcNo(), geom);
172  }
173 
174  // Receive points
176 
177  for
178  (
179  int slave=Pstream::firstSlave();
180  slave<=Pstream::lastSlave();
181  ++slave
182  )
183  {
184  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
185 
186  points.clear();
187 
188  fromSlave >> points;
189 
190  if (points.size())
191  {
192  multiPiece->SetPiece
193  (
194  slave,
196  );
197  }
198  }
199  }
200  else
201  {
202  // Slave - send face centres
203 
204  OPstream toMaster
205  (
206  Pstream::commsTypes::scheduled,
207  Pstream::masterNo()
208  );
209 
210  if (surf)
211  {
212  toMaster << surf->faceCentres();
213  }
214  else
215  {
216  toMaster << pointField();
217  }
218  }
219 
220  return multiPiece;
221 }
222 
223 
224 // ************************************************************************* //
surface.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::vtk::Tools::Vertices
vtkSmartPointer< vtkPolyData > Vertices(const UList< point > &pts)
Return vtkPolyData of vertices for each point.
Definition: foamVtkToolsI.H:104
foamVtkTools.H
Foam::functionObjects::runTimePostPro::surface::gatherFaceCentres
vtkSmartPointer< vtkMultiPieceDataSet > gatherFaceCentres(const polySurface *surf) const
Definition: surfaceGather.C:135
polySurface.H
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::vtk::Tools::Patch::faceCentres
static vtkSmartPointer< vtkPoints > faceCentres(const PatchType &p)
Return patch face centres as vtkPoints.
Definition: foamVtkToolsTemplates.C:156
vtkSmartPointer
Definition: runTimePostProcessing.H:148
Foam::polySurface
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:67
Foam::Field< vector >
Foam::polySurface::nFaces
virtual label nFaces() const
Return the number of faces.
Definition: polySurface.H:225
Foam::polySurface::points
virtual const pointField & points() const
Return points.
Definition: polySurface.H:238
Foam::vtk::Tools::identityVertices
vtkSmartPointer< vtkCellArray > identityVertices(const label size)
An identity list of VTK_VERTEX.
Definition: foamVtkToolsI.H:149
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
Foam::functionObjects::runTimePostPro::surface::gatherSurfacePieces
vtkSmartPointer< vtkMultiPieceDataSet > gatherSurfacePieces(const polySurface *surf) const
Definition: surfaceGather.C:44
runTimePostProcessing.H
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::List< face >
Foam::PrimitivePatch::faceCentres
const Field< PointType > & faceCentres() const
Return face centres for patch.
Definition: PrimitivePatch.C:517
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
Foam::polySurface::faces
virtual const faceList & faces() const
Return faces.
Definition: polySurface.H:244
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52