29 #ifndef vtkPVFoamPatchField_H
30 #define vtkPVFoamPatchField_H
36 #include "vtkCellData.h"
37 #include "vtkFloatArray.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkPointData.h"
40 #include "vtkPolyData.h"
45 void Foam::vtkPVFoam::convertPatchField
48 const Field<Type>& ptf,
49 vtkMultiBlockDataSet* output,
50 const arrayRange&
range,
54 const label nComp = pTraits<Type>::nComponents;
57 cellData->SetNumberOfTuples(ptf.size());
58 cellData->SetNumberOfComponents(nComp);
59 cellData->Allocate(nComp*ptf.size());
60 cellData->SetName(
name.c_str());
65 const Type& t = ptf[i];
70 vtkOpenFOAMTupleRemap<Type>(vec);
72 cellData->InsertTuple(i, vec);
75 vtkPolyData::SafeDownCast
77 GetDataSetFromBlock(output,
range, datasetNo)
86 void Foam::vtkPVFoam::convertPatchPointField
89 const Field<Type>& pptf,
90 vtkMultiBlockDataSet* output,
91 const arrayRange&
range,
95 const label nComp = pTraits<Type>::nComponents;
98 pointData->SetNumberOfTuples(pptf.size());
99 pointData->SetNumberOfComponents(nComp);
100 pointData->Allocate(nComp*pptf.size());
101 pointData->SetName(
name.c_str());
106 const Type& t = pptf[i];
111 vtkOpenFOAMTupleRemap<Type>(vec);
113 pointData->InsertTuple(i, vec);
116 vtkPolyData::SafeDownCast
118 GetDataSetFromBlock(output,
range, datasetNo)
120 ->AddArray(pointData);
#define forAll(list, i)
Loop across all elements in list.
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.
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.