vtkPVFoamPatchField.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 vtkPVFoamPatchField_H
30 #define vtkPVFoamPatchField_H
31 
32 #include "vtkPVFoam.H"
33 #include "vtkOpenFOAMTupleRemap.H"
34 
35 // VTK includes
36 #include "vtkCellData.h"
37 #include "vtkFloatArray.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkPointData.h"
40 #include "vtkPolyData.h"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type>
45 void Foam::vtkPVFoam::convertPatchField
46 (
47  const word& name,
48  const Field<Type>& ptf,
49  vtkMultiBlockDataSet* output,
50  const arrayRange& range,
51  const label datasetNo
52 )
53 {
54  const label nComp = pTraits<Type>::nComponents;
55 
56  vtkFloatArray* cellData = vtkFloatArray::New();
57  cellData->SetNumberOfTuples(ptf.size());
58  cellData->SetNumberOfComponents(nComp);
59  cellData->Allocate(nComp*ptf.size());
60  cellData->SetName(name.c_str());
61 
62  float vec[nComp];
63  forAll(ptf, i)
64  {
65  const Type& t = ptf[i];
66  for (direction d=0; d<nComp; ++d)
67  {
68  vec[d] = component(t, d);
69  }
70  vtkOpenFOAMTupleRemap<Type>(vec);
71 
72  cellData->InsertTuple(i, vec);
73  }
74 
75  vtkPolyData::SafeDownCast
76  (
77  GetDataSetFromBlock(output, range, datasetNo)
78  ) ->GetCellData()
79  ->AddArray(cellData);
80 
81  cellData->Delete();
82 }
83 
84 
85 template<class Type>
86 void Foam::vtkPVFoam::convertPatchPointField
87 (
88  const word& name,
89  const Field<Type>& pptf,
90  vtkMultiBlockDataSet* output,
91  const arrayRange& range,
92  const label datasetNo
93 )
94 {
95  const label nComp = pTraits<Type>::nComponents;
96 
97  vtkFloatArray* pointData = vtkFloatArray::New();
98  pointData->SetNumberOfTuples(pptf.size());
99  pointData->SetNumberOfComponents(nComp);
100  pointData->Allocate(nComp*pptf.size());
101  pointData->SetName(name.c_str());
102 
103  float vec[nComp];
104  forAll(pptf, i)
105  {
106  const Type& t = pptf[i];
107  for (direction d=0; d<nComp; ++d)
108  {
109  vec[d] = component(t, d);
110  }
111  vtkOpenFOAMTupleRemap<Type>(vec);
112 
113  pointData->InsertTuple(i, vec);
114  }
115 
116  vtkPolyData::SafeDownCast
117  (
118  GetDataSetFromBlock(output, range, datasetNo)
119  ) ->GetPointData()
120  ->AddArray(pointData);
121 
122  pointData->Delete();
123 }
124 
125 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126 
127 #endif
128 
129 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
uint8_t direction
Definition: direction.H:45