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-2025 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 "vtkPVFoam.H"
33 #include "vtkPVFoamReader.h"
34 #include "vtkOpenFOAMTupleRemap.H"
35 
36 // OpenFOAM includes
37 #include "domainDecomposition.H"
38 #include "emptyFvsPatchField.H"
39 #include "faceSet.H"
40 #include "surfaceFields.H"
41 #include "fvFieldReconstructor.H"
42 
43 // VTK includes
44 #include "vtkCellData.h"
45 #include "vtkFloatArray.h"
46 #include "vtkMultiBlockDataSet.h"
47 #include "vtkPolyData.h"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 template<class Type>
52 void Foam::vtkPVFoam::convertSurfaceField
53 (
54  const VolField<Type>& tf,
55  vtkMultiBlockDataSet* output,
56  const arrayRange& range,
57  const label datasetNo,
58  const fvMesh& mesh,
59  const labelList& faceLabels
60 )
61 {
62  const label nComp = pTraits<Type>::nComponents;
63  const label nInternalFaces = mesh.nInternalFaces();
64  const labelList& faceOwner = mesh.faceOwner();
65  const labelList& faceNeigh = mesh.faceNeighbour();
66 
67  vtkFloatArray* cellData = vtkFloatArray::New();
68  cellData->SetNumberOfTuples(faceLabels.size());
69  cellData->SetNumberOfComponents(nComp);
70  cellData->Allocate(nComp*faceLabels.size());
71  cellData->SetName(tf.name().c_str());
72 
74  << "Converting Surface field: " << tf.name()
75  << " size=" << tf.size() << " (" << faceLabels.size()
76  << "), nComp=" << nComp << endl;
77 
78  float vec[nComp];
79 
80  // For interior faces: average owner/neighbour
81  // For boundary faces: owner
82  forAll(faceLabels, i)
83  {
84  const label facei = faceLabels[i];
85  if (facei < nInternalFaces)
86  {
87  const Type t = 0.5*(tf[faceOwner[facei]] + tf[faceNeigh[facei]]);
88 
89  for (direction d=0; d<nComp; ++d)
90  {
91  vec[d] = component(t, d);
92  }
93  }
94  else
95  {
96  const label patchi = mesh.boundaryMesh().whichPatch(facei);
97  const label pFacei = mesh.boundaryMesh()[patchi].whichFace(facei);
98  const Type& t = tf.boundaryField()[patchi][pFacei];
99 
100  for (direction d=0; d<nComp; ++d)
101  {
102  vec[d] = component(t, d);
103  }
104  }
105  vtkOpenFOAMTupleRemap<Type>(vec);
106 
107  cellData->InsertTuple(i, vec);
108  }
109 
110 
111  vtkPolyData::SafeDownCast
112  (
113  GetDataSetFromBlock(output, range, datasetNo)
114  ) ->GetCellData()
115  ->AddArray(cellData);
116 
117  cellData->Delete();
118 }
119 
120 
121 template<class Type>
122 void Foam::vtkPVFoam::convertSurfaceField
123 (
124  const SurfaceField<Type>& tf,
125  vtkMultiBlockDataSet* output,
126  const arrayRange& range,
127  const label datasetNo,
128  const fvMesh& mesh,
129  const labelList& faceLabels
130 )
131 {
132  const label nComp = pTraits<Type>::nComponents;
133  const label nInternalFaces = mesh.nInternalFaces();
134 
135  vtkFloatArray* cellData = vtkFloatArray::New();
136  cellData->SetNumberOfTuples(faceLabels.size());
137  cellData->SetNumberOfComponents(nComp);
138  cellData->Allocate(nComp*faceLabels.size());
139  cellData->SetName(tf.name().c_str());
140 
142  << "Converting Surface field: " << tf.name()
143  << " size=" << tf.size() << " (" << faceLabels.size()
144  << "), nComp=" << nComp << endl;
145 
146  // To avoid whichPatch first flatten the field
147  Field<Type> flatFld(mesh.nFaces(), Zero);
148  SubField<Type>(flatFld, nInternalFaces) = tf.internalField();
149  forAll(tf.boundaryField(), patchi)
150  {
151  const fvsPatchField<Type>& fvs = tf.boundaryField()[patchi];
152 
153  SubField<Type>
154  (
155  flatFld,
156  fvs.size(),
157  fvs.patch().start()
158  ) = fvs;
159  }
160 
161  forAll(faceLabels, i)
162  {
163  const label facei = faceLabels[i];
164 
165  float vec[nComp];
166  for (direction d=0; d<nComp; ++d)
167  {
168  vec[d] = component(flatFld[facei], d);
169  }
170  vtkOpenFOAMTupleRemap<Type>(vec);
171 
172  cellData->InsertTuple(i, vec);
173  }
174 
175  vtkPolyData::SafeDownCast
176  (
177  GetDataSetFromBlock(output, range, datasetNo)
178  ) ->GetCellData()
179  ->AddArray(cellData);
180 
181  cellData->Delete();
182 }
183 
184 
185 template<class Type>
186 void Foam::vtkPVFoam::convertSurfaceFields
187 (
188  const IOobjectList& objects,
189  vtkMultiBlockDataSet* output
190 )
191 {
192  forAllConstIter(IOobjectList, objects, iter)
193  {
194  // Restrict to GeometricField<Type, ...>
195  if (iter()->headerClassName() != SurfaceField<Type>::typeName)
196  {
197  continue;
198  }
199 
200  // Load the field
201  tmp<SurfaceField<Type>> ttf;
202  if (reader_->GetDecomposedCase())
203  {
204  if (!fvReconstructorPtr_.valid())
205  {
206  fvReconstructorPtr_.set
207  (
208  new fvFieldReconstructor
209  (
210  procMeshesPtr_->completeMesh(),
211  procMeshesPtr_->procMeshes(),
212  procMeshesPtr_->procFaceAddressing(),
213  procMeshesPtr_->procCellAddressing(),
214  procMeshesPtr_->procFaceAddressingBf()
215  )
216  );
217  }
218 
219  ttf =
220  fvReconstructorPtr_
221  ->reconstructFvSurfaceField<Type>(*iter());
222  }
223  else
224  {
225  ttf =
226  new SurfaceField<Type>
227  (
228  *iter(),
229  procMeshesPtr_->completeMesh()
230  );
231  }
232  const SurfaceField<Type>& tf = ttf();
233 
234  // Convert patches - if activated
235  for
236  (
237  int partId = arrayRangePatches_.start();
238  partId < arrayRangePatches_.end();
239  ++partId
240  )
241  {
242  const word patchName = getPartName(partId);
243  const label datasetNo = partDataset_[partId];
244  const label patchId = tf.mesh().boundaryMesh().findIndex(patchName);
245 
246  if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
247  {
248  continue;
249  }
250 
251  const fvsPatchField<Type>& ptf = tf.boundaryField()[patchId];
252 
253  if (!isType<emptyFvsPatchField<Type>>(ptf))
254  {
255  convertPatchField
256  (
257  tf.name(),
258  ptf,
259  output,
260  arrayRangePatches_,
261  datasetNo
262  );
263  }
264  }
265 
266  // Convert face zones - if activated
267  for
268  (
269  int partId = arrayRangeFaceZones_.start();
270  partId < arrayRangeFaceZones_.end();
271  ++partId
272  )
273  {
274  const word zoneName = getPartName(partId);
275  const label datasetNo = partDataset_[partId];
276 
277  if (!partStatus_[partId] || datasetNo < 0)
278  {
279  continue;
280  }
281 
282  const faceZoneList& zMesh = tf.mesh().faceZones();
283  const label zoneId = zMesh.findIndex(zoneName);
284 
285  if (zoneId < 0)
286  {
287  continue;
288  }
289 
290  convertSurfaceField
291  (
292  tf,
293  output,
294  arrayRangeFaceZones_,
295  datasetNo,
296  tf.mesh(),
297  zMesh[zoneId]
298  );
299  }
300 
301  // Convert face sets - if activated
302  for
303  (
304  int partId = arrayRangeFaceSets_.start();
305  partId < arrayRangeFaceSets_.end();
306  ++partId
307  )
308  {
309  const word selectName = getPartName(partId);
310  const label datasetNo = partDataset_[partId];
311 
312  if (!partStatus_[partId] || datasetNo < 0)
313  {
314  continue;
315  }
316 
317  const autoPtr<faceSet> fSetPtr =
318  reader_->GetDecomposedCase()
319  ? procMeshesPtr_->reconstructSet<faceSet>(selectName)
320  : autoPtr<faceSet>(new faceSet(tf.mesh(), selectName));
321 
322  convertSurfaceField
323  (
324  tf,
325  output,
326  arrayRangeFaceSets_,
327  datasetNo,
328  tf.mesh(),
329  fSetPtr().toc()
330  );
331  }
332  }
333 }
334 
335 
336 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 
338 #endif
339 
340 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1363
label nInternalFaces() const
label nFaces() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
const tensorField & tf
label patchId(-1)
#define DebugInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static const zero Zero
Definition: zero.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:163
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
uint8_t direction
Definition: direction.H:45
objects
Foam::surfaceFields.