geometryPatchesTemplates.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 "foamVtkTools.H"
33 
34 // VTK includes
35 #include "vtkCellData.h"
36 #include "vtkMultiPieceDataSet.h"
37 #include "vtkPointData.h"
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class Type>
43 (
44  vtkMultiPieceDataSet* multiPiece,
45  const labelListList& patchIds,
47  const word& fieldName
48 ) const
49 {
50  if (!multiPiece || !fldptr)
51  {
52  return 0;
53  }
54 
55  const int nCmpt(pTraits<Type>::nComponents);
56 
57  const auto& bf = fldptr->boundaryField();
58 
59  label pieceId = 0;
60 
61  if (!needsCollective())
62  {
63  // Simple case (serial-serial, parallel-parallel)
64 
65  for (int proci=0; proci < Pstream::myProcNo(); ++proci)
66  {
67  pieceId += patchIds[proci].size();
68  }
69 
70  for (const label patchId : patchIds[Pstream::myProcNo()])
71  {
72  const auto& pf = bf[patchId];
73 
74  auto vtkfield =
75  (
76  nearCellValue_
77  ? Foam::vtk::Tools::convertFieldToVTK<Type>
78  (
79  fieldName,
80  pf.patchInternalField()()
81  )
82  : Foam::vtk::Tools::convertFieldToVTK<Type>
83  (
84  fieldName,
85  pf
86  )
87  );
88 
89  auto piece = multiPiece->GetPiece(pieceId);
90 
91  if (piece->GetNumberOfCells() == piece->GetNumberOfPoints())
92  {
93  // Only has verts
94  piece->GetPointData()->AddArray(vtkfield);
95  }
96  else
97  {
98  piece->GetCellData()->AddArray(vtkfield);
99  }
100 
101  ++pieceId;
102  }
103  }
104  else if (Pstream::master())
105  {
106  // Gather pieces on master
107 
108  // Add myself
109  for (const label patchId : patchIds[Pstream::myProcNo()])
110  {
111  const auto& pf = bf[patchId];
112 
113  auto vtkfield =
114  (
115  nearCellValue_
116  ? Foam::vtk::Tools::convertFieldToVTK<Type>
117  (
118  fieldName,
119  pf.patchInternalField()()
120  )
121  : Foam::vtk::Tools::convertFieldToVTK<Type>
122  (
123  fieldName,
124  pf
125  )
126  );
127 
128  auto piece = multiPiece->GetPiece(pieceId);
129 
130  if (piece->GetNumberOfCells() == piece->GetNumberOfPoints())
131  {
132  // Only has verts
133  piece->GetPointData()->AddArray(vtkfield);
134  }
135  else
136  {
137  piece->GetCellData()->AddArray(vtkfield);
138  }
139 
140  ++pieceId;
141  }
142 
143  // Receive field
144  Field<Type> recv;
145 
146  for
147  (
148  int slave=Pstream::firstSlave();
149  slave<=Pstream::lastSlave();
150  ++slave
151  )
152  {
153  const label nSlavePatches = patchIds[slave].size();
154 
155  if (!nSlavePatches)
156  {
157  continue;
158  }
159 
160  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
161 
162  for (label recvi=0; recvi < nSlavePatches; ++recvi)
163  {
164  recv.clear();
165 
166  fromSlave >> recv;
167 
168  auto vtkfield = Foam::vtk::Tools::convertFieldToVTK<Type>
169  (
170  fieldName,
171  recv
172  );
173 
174  auto piece = multiPiece->GetPiece(pieceId);
175 
176  if (piece->GetNumberOfCells() == piece->GetNumberOfPoints())
177  {
178  // Only has verts
179  piece->GetPointData()->AddArray(vtkfield);
180  }
181  else
182  {
183  piece->GetCellData()->AddArray(vtkfield);
184  }
185 
186  ++pieceId;
187  }
188  }
189  }
190  else
191  {
192  // Slave - send fields
193  const labelList& slavePatchIds = patchIds[Pstream::myProcNo()];
194 
195  if (slavePatchIds.size())
196  {
197  OPstream toMaster
198  (
199  Pstream::commsTypes::scheduled,
200  Pstream::masterNo()
201  );
202 
203  for (const label patchId : patchIds[Pstream::myProcNo()])
204  {
205  const auto& pf = bf[patchId];
206 
207  if (nearCellValue_)
208  {
209  toMaster << pf.patchInternalField()();
210  }
211  else
212  {
213  toMaster << pf;
214  }
215  }
216  }
217  }
218 
219  return nCmpt;
220 }
221 
222 
223 // ************************************************************************* //
volFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
foamVtkTools.H
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
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
Generic templated field type.
Definition: Field.H:63
Foam::functionObjects::runTimePostPro::geometryPatches::addPatchField
int addPatchField(vtkMultiPieceDataSet *multiPiece, const labelListList &patchIds, const GeometricField< Type, fvPatchField, volMesh > *fldptr, const word &fieldName) const
Add patch values.
Definition: geometryPatchesTemplates.C:43
fvMesh.H
Foam::List< labelList >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
patchId
label patchId(-1)
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62