vtkPVFoamSurfaceField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 InClass
25  vtkPVFoam
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #ifndef vtkPVFoamSurfaceField_H
30 #define vtkPVFoamSurfaceField_H
31 
32 #include "surfaceFields.H"
33 #include "emptyFvsPatchField.H"
34 
35 // VTK includes
36 #include "vtkCellData.h"
37 #include "vtkFloatArray.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkPolyData.h"
40 
41 #include "vtkOpenFOAMTupleRemap.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 template<class Type>
46 void Foam::vtkPVFoam::convertSurfaceField
47 (
48  const GeometricField<Type, fvPatchField, volMesh>& tf,
49  vtkMultiBlockDataSet* output,
50  const arrayRange& range,
51  const label datasetNo,
52  const fvMesh& mesh,
53  const labelList& faceLabels
54 )
55 {
56  const label nComp = pTraits<Type>::nComponents;
57  const label nInternalFaces = mesh.nInternalFaces();
58  const labelList& faceOwner = mesh.faceOwner();
59  const labelList& faceNeigh = mesh.faceNeighbour();
60 
61  vtkFloatArray* cellData = vtkFloatArray::New();
62  cellData->SetNumberOfTuples(faceLabels.size());
63  cellData->SetNumberOfComponents(nComp);
64  cellData->Allocate(nComp*faceLabels.size());
65  cellData->SetName(tf.name().c_str());
66 
67  if (debug)
68  {
69  Info<< "convert convertSurfaceField: "
70  << tf.name()
71  << " size = " << tf.size()
72  << " nComp=" << nComp
73  << " nTuples = " << faceLabels.size() << endl;
74  }
75 
76  float vec[nComp];
77 
78  // For interior faces: average owner/neighbour
79  // For boundary faces: owner
80  forAll(faceLabels, facei)
81  {
82  const label faceNo = faceLabels[facei];
83  if (faceNo < nInternalFaces)
84  {
85  Type t = 0.5*(tf[faceOwner[faceNo]] + tf[faceNeigh[faceNo]]);
86 
87  for (direction d=0; d<nComp; ++d)
88  {
89  vec[d] = component(t, d);
90  }
91  }
92  else
93  {
94  const Type& t = tf[faceOwner[faceNo]];
95  for (direction d=0; d<nComp; ++d)
96  {
97  vec[d] = component(t, d);
98  }
99  }
100  vtkOpenFOAMTupleRemap<Type>(vec);
101 
102  cellData->InsertTuple(facei, vec);
103  }
104 
105 
106  vtkPolyData::SafeDownCast
107  (
108  GetDataSetFromBlock(output, range, datasetNo)
109  ) ->GetCellData()
110  ->AddArray(cellData);
111 
112  cellData->Delete();
113 }
114 
115 
116 template<class Type>
117 void Foam::vtkPVFoam::convertSurfaceField
118 (
119  const GeometricField<Type, fvsPatchField, surfaceMesh>& tf,
120  vtkMultiBlockDataSet* output,
121  const arrayRange& range,
122  const label datasetNo,
123  const fvMesh& mesh,
124  const labelList& faceLabels
125 )
126 {
127  const label nComp = pTraits<Type>::nComponents;
128  const label nInternalFaces = mesh.nInternalFaces();
129 
130  vtkFloatArray* cellData = vtkFloatArray::New();
131  cellData->SetNumberOfTuples(faceLabels.size());
132  cellData->SetNumberOfComponents(nComp);
133  cellData->Allocate(nComp*faceLabels.size());
134  cellData->SetName(tf.name().c_str());
135 
136  if (debug)
137  {
138  Info<< "convert convertSurfaceField: "
139  << tf.name()
140  << " size = " << tf.size()
141  << " nComp=" << nComp
142  << " nTuples = " << faceLabels.size() << endl;
143  }
144 
145  // To avoid whichPatch first flatten the field
146  Field<Type> flatFld(mesh.nFaces(), Zero);
147  SubField<Type>(flatFld, nInternalFaces) = tf.internalField();
148  forAll(tf.boundaryField(), patchi)
149  {
150  const fvsPatchField<Type>& fvs = tf.boundaryField()[patchi];
151 
152  SubField<Type>
153  (
154  flatFld,
155  fvs.size(),
156  fvs.patch().start()
157  ) = fvs;
158  }
159 
160 
161  // For interior faces: average owner/neighbour
162  // For boundary faces: owner
163  forAll(faceLabels, facei)
164  {
165  const label faceNo = faceLabels[facei];
166 
167  float vec[nComp];
168  for (direction d=0; d<nComp; ++d)
169  {
170  vec[d] = component(flatFld[faceNo], d);
171  }
172  vtkOpenFOAMTupleRemap<Type>(vec);
173 
174  cellData->InsertTuple(facei, vec);
175  }
176 
177  vtkPolyData::SafeDownCast
178  (
179  GetDataSetFromBlock(output, range, datasetNo)
180  ) ->GetCellData()
181  ->AddArray(cellData);
182 
183  cellData->Delete();
184 }
185 
186 
187 template<class Type>
188 void Foam::vtkPVFoam::convertSurfaceFields
189 (
190  const fvMesh& mesh,
191  const IOobjectList& objects,
192  vtkMultiBlockDataSet* output
193 )
194 {
195  forAllConstIter(IOobjectList, objects, iter)
196  {
197  // Restrict to GeometricField<Type, ...>
198  if
199  (
200  iter()->headerClassName()
202  )
203  {
204  continue;
205  }
206 
207  // Load field
208  const GeometricField<Type, fvsPatchField, surfaceMesh> tf
209  (
210  *iter(),
211  mesh
212  );
213 
214  // Convert patches - if activated
215  for
216  (
217  int partId = arrayRangePatches_.start();
218  partId < arrayRangePatches_.end();
219  ++partId
220  )
221  {
222  const word patchName = getPartName(partId);
223  const label datasetNo = partDataset_[partId];
224  const label patchId = mesh.boundaryMesh().findPatchID(patchName);
225 
226  if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
227  {
228  continue;
229  }
230 
231  const fvsPatchField<Type>& ptf = tf.boundaryField()[patchId];
232 
233  if (!isType<emptyFvsPatchField<Type>>(ptf))
234  {
235  convertPatchField
236  (
237  tf.name(),
238  ptf,
239  output,
240  arrayRangePatches_,
241  datasetNo
242  );
243  }
244  }
245 
246 
247  // Convert face zones - if activated
248  for
249  (
250  int partId = arrayRangeFaceZones_.start();
251  partId < arrayRangeFaceZones_.end();
252  ++partId
253  )
254  {
255  const word zoneName = getPartName(partId);
256  const label datasetNo = partDataset_[partId];
257 
258  if (!partStatus_[partId] || datasetNo < 0)
259  {
260  continue;
261  }
262 
263  const meshFaceZones& zMesh = mesh.faceZones();
264  const label zoneId = zMesh.findZoneID(zoneName);
265 
266  if (zoneId < 0)
267  {
268  continue;
269  }
270 
271  convertSurfaceField
272  (
273  tf,
274  output,
275  arrayRangeFaceZones_,
276  datasetNo,
277  mesh,
278  zMesh[zoneId]
279  );
280  }
281 
282  // Convert face sets - if activated
283  for
284  (
285  int partId = arrayRangeFaceSets_.start();
286  partId < arrayRangeFaceSets_.end();
287  ++partId
288  )
289  {
290  const word selectName = getPartName(partId);
291  const label datasetNo = partDataset_[partId];
292 
293  if (!partStatus_[partId] || datasetNo < 0)
294  {
295  continue;
296  }
297 
298  const faceSet fSet(mesh, selectName);
299 
300  convertSurfaceField
301  (
302  tf,
303  output,
304  arrayRangeFaceSets_,
305  datasetNo,
306  mesh,
307  fSet.toc()
308  );
309  }
310  }
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 #endif
317 
318 // ************************************************************************* //
Foam::surfaceFields.
label patchId(-1)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static const char *const typeName
Definition: Field.H:105
uint8_t direction
Definition: direction.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const tensorField & tf
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:126
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
MeshZones< faceZone, polyMesh > meshFaceZones
A MeshZones with the type faceZone.
label patchi
messageStream Info
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)