29 #ifndef vtkPVFoamSurfaceField_H
30 #define vtkPVFoamSurfaceField_H
33 #include "vtkPVFoamReader.h"
44 #include "vtkCellData.h"
45 #include "vtkFloatArray.h"
46 #include "vtkMultiBlockDataSet.h"
47 #include "vtkPolyData.h"
52 void Foam::vtkPVFoam::convertSurfaceField
54 const VolField<Type>&
tf,
55 vtkMultiBlockDataSet* output,
56 const arrayRange&
range,
57 const label datasetNo,
62 const label nComp = pTraits<Type>::nComponents;
68 cellData->SetNumberOfTuples(faceLabels.size());
69 cellData->SetNumberOfComponents(nComp);
70 cellData->Allocate(nComp*faceLabels.size());
71 cellData->SetName(
tf.name().c_str());
74 <<
"Converting Surface field: " <<
tf.name()
75 <<
" size=" <<
tf.size() <<
" (" << faceLabels.size()
76 <<
"), nComp=" << nComp <<
endl;
84 const label facei = faceLabels[i];
85 if (facei < nInternalFaces)
87 const Type t = 0.5*(
tf[faceOwner[facei]] +
tf[faceNeigh[facei]]);
98 const Type& t =
tf.boundaryField()[
patchi][pFacei];
105 vtkOpenFOAMTupleRemap<Type>(vec);
107 cellData->InsertTuple(i, vec);
111 vtkPolyData::SafeDownCast
113 GetDataSetFromBlock(output,
range, datasetNo)
115 ->AddArray(cellData);
122 void Foam::vtkPVFoam::convertSurfaceField
124 const SurfaceField<Type>&
tf,
125 vtkMultiBlockDataSet* output,
126 const arrayRange&
range,
127 const label datasetNo,
132 const label nComp = pTraits<Type>::nComponents;
136 cellData->SetNumberOfTuples(faceLabels.size());
137 cellData->SetNumberOfComponents(nComp);
138 cellData->Allocate(nComp*faceLabels.size());
139 cellData->SetName(
tf.name().c_str());
142 <<
"Converting Surface field: " <<
tf.name()
143 <<
" size=" <<
tf.
size() <<
" (" << faceLabels.size()
144 <<
"), nComp=" << nComp <<
endl;
148 SubField<Type>(flatFld, nInternalFaces) =
tf.internalField();
151 const fvsPatchField<Type>& fvs =
tf.boundaryField()[
patchi];
163 const label facei = faceLabels[i];
170 vtkOpenFOAMTupleRemap<Type>(vec);
172 cellData->InsertTuple(i, vec);
175 vtkPolyData::SafeDownCast
177 GetDataSetFromBlock(output,
range, datasetNo)
179 ->AddArray(cellData);
186 void Foam::vtkPVFoam::convertSurfaceFields
189 vtkMultiBlockDataSet* output
195 if (iter()->headerClassName() != SurfaceField<Type>::typeName)
201 tmp<SurfaceField<Type>> ttf;
202 if (reader_->GetDecomposedCase())
204 if (!fvReconstructorPtr_.valid())
206 fvReconstructorPtr_.set
208 new fvFieldReconstructor
210 procMeshesPtr_->completeMesh(),
211 procMeshesPtr_->procMeshes(),
212 procMeshesPtr_->procFaceAddressing(),
213 procMeshesPtr_->procCellAddressing(),
214 procMeshesPtr_->procFaceAddressingBf()
221 ->reconstructFvSurfaceField<Type>(*iter());
226 new SurfaceField<Type>
229 procMeshesPtr_->completeMesh()
232 const SurfaceField<Type>&
tf = ttf();
237 int partId = arrayRangePatches_.start();
238 partId < arrayRangePatches_.end();
242 const word patchName = getPartName(partId);
243 const label datasetNo = partDataset_[partId];
244 const label patchId =
tf.mesh().boundaryMesh().findIndex(patchName);
246 if (!partStatus_[partId] || datasetNo < 0 ||
patchId < 0)
251 const fvsPatchField<Type>& ptf =
tf.boundaryField()[
patchId];
253 if (!
isType<emptyFvsPatchField<Type>>(ptf))
269 int partId = arrayRangeFaceZones_.start();
270 partId < arrayRangeFaceZones_.end();
274 const word zoneName = getPartName(partId);
275 const label datasetNo = partDataset_[partId];
277 if (!partStatus_[partId] || datasetNo < 0)
282 const faceZoneList& zMesh =
tf.mesh().faceZones();
283 const label zoneId = zMesh.findIndex(zoneName);
294 arrayRangeFaceZones_,
304 int partId = arrayRangeFaceSets_.start();
305 partId < arrayRangeFaceSets_.end();
309 const word selectName = getPartName(partId);
310 const label datasetNo = partDataset_[partId];
312 if (!partStatus_[partId] || datasetNo < 0)
317 const autoPtr<faceSet> fSetPtr =
318 reader_->GetDecomposedCase()
319 ? procMeshesPtr_->reconstructSet<faceSet>(selectName)
320 : autoPtr<faceSet>(
new faceSet(
tf.mesh(), selectName));
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
void size(const label)
Override size to be inconsistent with allocated storage.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
virtual const labelList & faceOwner() const
Return face owner.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
virtual const labelList & faceNeighbour() const
Return face neighbour.
label nInternalFaces() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define DebugInFunction
Report an information message using Foam::Info.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
bool isType(const Type &t)
Check the typeid.
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)