vtkPVFoamVolFields.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 vtkPVFoamVolFields_H
30 #define vtkPVFoamVolFields_H
31 
32 #include "vtkPVFoam.H"
33 #include "vtkPVFoamReader.h"
34 #include "vtkPVFoamSurfaceField.H"
35 #include "vtkPVFoamPatchField.H"
36 #include "vtkOpenFOAMTupleRemap.H"
37 
38 // OpenFOAM includes
39 #include "domainDecomposition.H"
40 #include "emptyFvPatchField.H"
41 #include "faceSet.H"
42 #include "volFields.H"
43 #include "volPointInterpolation.H"
44 #include "fvFieldReconstructor.H"
45 
46 // VTK includes
47 #include "vtkFloatArray.h"
48 #include "vtkUnstructuredGrid.h"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 template<class Type>
53 void Foam::vtkPVFoam::convertVolFields
54 (
55  const PtrList<PrimitivePatchInterpolation<primitivePatch>>& ppInterpList,
56  const IOobjectList& objects,
57  const bool interpFields,
58  vtkMultiBlockDataSet* output
59 )
60 {
61  const polyBoundaryMesh& patches =
62  procMeshesPtr_->completeMesh().boundaryMesh();
63 
64  forAllConstIter(IOobjectList, objects, iter)
65  {
66  // Restrict to GeometricField<Type, ...>
67  if (iter()->headerClassName() != VolField<Type>::typeName) continue;
68 
69  // Load the field
70  tmp<VolField<Type>> ttf;
71  if (reader_->GetDecomposedCase())
72  {
73  if (!fvReconstructorPtr_.valid())
74  {
75  fvReconstructorPtr_.set
76  (
77  new fvFieldReconstructor
78  (
79  procMeshesPtr_->completeMesh(),
80  procMeshesPtr_->procMeshes(),
81  procMeshesPtr_->procFaceAddressing(),
82  procMeshesPtr_->procCellAddressing(),
83  procMeshesPtr_->procFaceAddressingBf()
84  )
85  );
86  }
87 
88  ttf =
89  fvReconstructorPtr_
90  ->reconstructVolField<Type>(*iter());
91  }
92  else
93  {
94  ttf =
95  new VolField<Type>
96  (
97  *iter(),
98  procMeshesPtr_->completeMesh()
99  );
100  }
101  const VolField<Type>& tf = ttf();
102 
103  // Interpolated field (demand driven)
104  autoPtr<PointField<Type>> ptfPtr;
105  if (interpFields)
106  {
108  << "Interpolating field " << tf.name() << endl;
109 
110  ptfPtr.reset
111  (
113  );
114  }
115 
116  // Convert activated internalMesh regions
117  convertVolFieldBlock
118  (
119  tf,
120  ptfPtr,
121  output,
122  arrayRangeVolume_,
123  regionPolyDecomp_
124  );
125 
126  // Convert activated cellZones
127  convertVolFieldBlock
128  (
129  tf,
130  ptfPtr,
131  output,
132  arrayRangeCellZones_,
133  zonePolyDecomp_
134  );
135 
136  // Convert activated cellSets
137  convertVolFieldBlock
138  (
139  tf,
140  ptfPtr,
141  output,
142  arrayRangeCellSets_,
143  setPolyDecomp_
144  );
145 
146 
147  // Convert patches - if activated
148  for
149  (
150  int partId = arrayRangePatches_.start();
151  partId < arrayRangePatches_.end();
152  ++partId
153  )
154  {
155  const word patchName = getPartName(partId);
156  const label datasetNo = partDataset_[partId];
157  const label patchId = patches.findIndex(patchName);
158 
159  if (!partStatus_[partId] || datasetNo < 0 || patchId < 0) continue;
160 
161  const fvPatchField<Type>& ptf = tf.boundaryField()[patchId];
162 
163  if
164  (
165  isType<emptyFvPatchField<Type>>(ptf)
166  ||
167  (
168  reader_->GetExtrapolatePatches()
170  )
171  )
172  {
173  fvPatch p(ptf.patch().patch(), tf.mesh().boundary());
174 
175  tmp<Field<Type>> tpptf
176  (
177  fvPatchField<Type>(p, tf).patchInternalField()
178  );
179 
180  convertPatchField
181  (
182  tf.name(),
183  tpptf(),
184  output,
185  arrayRangePatches_,
186  datasetNo
187  );
188 
189  if (interpFields)
190  {
191  convertPatchPointField
192  (
193  tf.name(),
194  ppInterpList[patchId].faceToPointInterpolate(tpptf)(),
195  output,
196  arrayRangePatches_,
197  datasetNo
198  );
199  }
200  }
201  else
202  {
203  convertPatchField
204  (
205  tf.name(),
206  ptf,
207  output,
208  arrayRangePatches_,
209  datasetNo
210  );
211 
212  if (interpFields)
213  {
214  convertPatchPointField
215  (
216  tf.name(),
217  ppInterpList[patchId].faceToPointInterpolate(ptf)(),
218  output,
219  arrayRangePatches_,
220  datasetNo
221  );
222  }
223  }
224  }
225 
226  // Convert face zones - if activated
227  for
228  (
229  int partId = arrayRangeFaceZones_.start();
230  partId < arrayRangeFaceZones_.end();
231  ++partId
232  )
233  {
234  const word zoneName = getPartName(partId);
235  const label datasetNo = partDataset_[partId];
236 
237  if (!partStatus_[partId] || datasetNo < 0) continue;
238 
239  const faceZoneList& zMesh = tf.mesh().faceZones();
240  const label zoneId = zMesh.findIndex(zoneName);
241 
242  if (zoneId < 0) continue;
243 
244  convertSurfaceField
245  (
246  tf,
247  output,
248  arrayRangeFaceZones_,
249  datasetNo,
250  tf.mesh(),
251  zMesh[zoneId]
252  );
253  }
254 
255  // Convert face sets - if activated
256  for
257  (
258  int partId = arrayRangeFaceSets_.start();
259  partId < arrayRangeFaceSets_.end();
260  ++partId
261  )
262  {
263  const word selectName = getPartName(partId);
264  const label datasetNo = partDataset_[partId];
265 
266  if (!partStatus_[partId] || datasetNo < 0) continue;
267 
268  const autoPtr<faceSet> fSetPtr =
269  reader_->GetDecomposedCase()
270  ? procMeshesPtr_->reconstructSet<faceSet>(selectName)
271  : autoPtr<faceSet>(new faceSet(tf.mesh(), selectName));
272 
273  convertSurfaceField
274  (
275  tf,
276  output,
277  arrayRangeFaceSets_,
278  datasetNo,
279  tf.mesh(),
280  fSetPtr().toc()
281  );
282  }
283  }
284 }
285 
286 
287 template<class Type>
288 void Foam::vtkPVFoam::convertVolInternalFields
289 (
290  const IOobjectList& objects,
291  vtkMultiBlockDataSet* output
292 )
293 {
294  forAllConstIter(IOobjectList, objects, iter)
295  {
296  // Restrict to GeometricField<Type, ...>::Internal
297  if (iter()->headerClassName() != VolInternalField<Type>::typeName)
298  continue;
299 
300  // Load the field
301  tmp<VolInternalField<Type>> ttf;
302  if (reader_->GetDecomposedCase())
303  {
304  if (!fvReconstructorPtr_.valid())
305  {
306  fvReconstructorPtr_.set
307  (
308  new fvFieldReconstructor
309  (
310  procMeshesPtr_->completeMesh(),
311  procMeshesPtr_->procMeshes(),
312  procMeshesPtr_->procFaceAddressing(),
313  procMeshesPtr_->procCellAddressing(),
314  procMeshesPtr_->procFaceAddressingBf()
315  )
316  );
317  }
318 
319  ttf =
320  fvReconstructorPtr_
321  ->reconstructVolInternalField<Type>(*iter());
322  }
323  else
324  {
325  ttf =
326  new VolInternalField<Type>
327  (
328  *iter(),
329  procMeshesPtr_->completeMesh()
330  );
331  }
332  const VolInternalField<Type>& tf = ttf();
333 
334  // Convert activated internalMesh regions
335  convertVolInternalFieldBlock<Type>
336  (
337  tf,
338  output,
339  arrayRangeVolume_,
340  regionPolyDecomp_
341  );
342 
343  // Convert activated cellZones
344  convertVolInternalFieldBlock<Type>
345  (
346  tf,
347  output,
348  arrayRangeCellZones_,
349  zonePolyDecomp_
350  );
351 
352  // Convert activated cellSets
353  convertVolInternalFieldBlock<Type>
354  (
355  tf,
356  output,
357  arrayRangeCellSets_,
358  setPolyDecomp_
359  );
360  }
361 }
362 
363 
364 template<class Type>
365 void Foam::vtkPVFoam::convertVolFieldBlock
366 (
367  const VolField<Type>& tf,
368  autoPtr<PointField<Type>>& ptfPtr,
369  vtkMultiBlockDataSet* output,
370  const arrayRange& range,
371  const List<polyDecomp>& decompLst
372 )
373 {
374  for (int partId = range.start(); partId < range.end(); ++partId)
375  {
376  const label datasetNo = partDataset_[partId];
377 
378  if (datasetNo >= 0 && partStatus_[partId])
379  {
380  convertVolInternalField<Type>
381  (
382  tf,
383  output,
384  range,
385  datasetNo,
386  decompLst[datasetNo]
387  );
388 
389  if (ptfPtr.valid())
390  {
391  convertPointField
392  (
393  ptfPtr(),
394  tf,
395  output,
396  range,
397  datasetNo,
398  decompLst[datasetNo]
399  );
400  }
401  }
402  }
403 }
404 
405 
406 template<class Type>
407 void Foam::vtkPVFoam::convertVolInternalFieldBlock
408 (
409  const VolInternalField<Type>& tf,
410  vtkMultiBlockDataSet* output,
411  const arrayRange& range,
412  const List<polyDecomp>& decompLst
413 )
414 {
415  for (int partId = range.start(); partId < range.end(); ++partId)
416  {
417  const label datasetNo = partDataset_[partId];
418 
419  if (datasetNo >= 0 && partStatus_[partId])
420  {
421  convertVolInternalField<Type>
422  (
423  tf,
424  output,
425  range,
426  datasetNo,
427  decompLst[datasetNo]
428  );
429  }
430  }
431 }
432 
433 
434 template<class Type>
435 void Foam::vtkPVFoam::convertVolInternalField
436 (
437  const VolInternalField<Type>& tf,
438  vtkMultiBlockDataSet* output,
439  const arrayRange& range,
440  const label datasetNo,
441  const polyDecomp& decompInfo
442 )
443 {
444  const label nComp = pTraits<Type>::nComponents;
445  const labelList& superCells = decompInfo.superCells();
446 
447  vtkFloatArray* celldata = vtkFloatArray::New();
448  celldata->SetNumberOfTuples(superCells.size());
449  celldata->SetNumberOfComponents(nComp);
450  celldata->Allocate(nComp*superCells.size());
451  celldata->SetName(tf.name().c_str());
452 
454  << "Converting Vol field: " << tf.name()
455  << " size=" << tf.size() << " (" << superCells.size()
456  << "), nComp=" << nComp << endl;
457 
458  float vec[nComp];
459  forAll(superCells, i)
460  {
461  const Type& t = tf[superCells[i]];
462  for (direction d=0; d<nComp; ++d)
463  {
464  vec[d] = component(t, d);
465  }
466  vtkOpenFOAMTupleRemap<Type>(vec);
467 
468  celldata->InsertTuple(i, vec);
469  }
470 
471  vtkUnstructuredGrid::SafeDownCast
472  (
473  GetDataSetFromBlock(output, range, datasetNo)
474  ) ->GetCellData()
475  ->AddArray(celldata);
476 
477  celldata->Delete();
478 }
479 
480 
481 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
482 
483 #endif
484 
485 // ************************************************************************* //
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
static volPointInterpolation & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Traits class for primitives.
Definition: pTraits.H:53
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
const tensorField & tf
label patchId(-1)
const fvPatchList & patches
#define DebugInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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)
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:59
uint8_t direction
Definition: direction.H:45
objects
volScalarField & p