29 #ifndef vtkPVFoamVolFields_H
30 #define vtkPVFoamVolFields_H
33 #include "vtkPVFoamReader.h"
47 #include "vtkFloatArray.h"
48 #include "vtkUnstructuredGrid.h"
53 void Foam::vtkPVFoam::convertVolFields
55 const PtrList<PrimitivePatchInterpolation<primitivePatch>>& ppInterpList,
57 const bool interpFields,
58 vtkMultiBlockDataSet* output
61 const polyBoundaryMesh&
patches =
62 procMeshesPtr_->completeMesh().boundaryMesh();
67 if (iter()->headerClassName() != VolField<Type>::typeName)
continue;
70 tmp<VolField<Type>> ttf;
71 if (reader_->GetDecomposedCase())
73 if (!fvReconstructorPtr_.valid())
75 fvReconstructorPtr_.set
77 new fvFieldReconstructor
79 procMeshesPtr_->completeMesh(),
80 procMeshesPtr_->procMeshes(),
81 procMeshesPtr_->procFaceAddressing(),
82 procMeshesPtr_->procCellAddressing(),
83 procMeshesPtr_->procFaceAddressingBf()
90 ->reconstructVolField<Type>(*iter());
98 procMeshesPtr_->completeMesh()
101 const VolField<Type>&
tf = ttf();
104 autoPtr<PointField<Type>> ptfPtr;
108 <<
"Interpolating field " <<
tf.name() <<
endl;
132 arrayRangeCellZones_,
150 int partId = arrayRangePatches_.start();
151 partId < arrayRangePatches_.end();
155 const word patchName = getPartName(partId);
156 const label datasetNo = partDataset_[partId];
159 if (!partStatus_[partId] || datasetNo < 0 ||
patchId < 0)
continue;
161 const fvPatchField<Type>& ptf =
tf.boundaryField()[
patchId];
165 isType<emptyFvPatchField<Type>>(ptf)
168 reader_->GetExtrapolatePatches()
173 fvPatch
p(ptf.patch().patch(),
tf.mesh().boundary());
175 tmp<Field<Type>> tpptf
177 fvPatchField<Type>(
p,
tf).patchInternalField()
191 convertPatchPointField
194 ppInterpList[
patchId].faceToPointInterpolate(tpptf)(),
214 convertPatchPointField
217 ppInterpList[
patchId].faceToPointInterpolate(ptf)(),
229 int partId = arrayRangeFaceZones_.start();
230 partId < arrayRangeFaceZones_.end();
234 const word zoneName = getPartName(partId);
235 const label datasetNo = partDataset_[partId];
237 if (!partStatus_[partId] || datasetNo < 0)
continue;
239 const faceZoneList& zMesh =
tf.mesh().faceZones();
240 const label zoneId = zMesh.findIndex(zoneName);
242 if (zoneId < 0)
continue;
248 arrayRangeFaceZones_,
258 int partId = arrayRangeFaceSets_.start();
259 partId < arrayRangeFaceSets_.end();
263 const word selectName = getPartName(partId);
264 const label datasetNo = partDataset_[partId];
266 if (!partStatus_[partId] || datasetNo < 0)
continue;
268 const autoPtr<faceSet> fSetPtr =
269 reader_->GetDecomposedCase()
270 ? procMeshesPtr_->reconstructSet<faceSet>(selectName)
271 : autoPtr<faceSet>(
new faceSet(
tf.mesh(), selectName));
288 void Foam::vtkPVFoam::convertVolInternalFields
291 vtkMultiBlockDataSet* output
297 if (iter()->headerClassName() != VolInternalField<Type>::typeName)
301 tmp<VolInternalField<Type>> ttf;
302 if (reader_->GetDecomposedCase())
304 if (!fvReconstructorPtr_.valid())
306 fvReconstructorPtr_.set
308 new fvFieldReconstructor
310 procMeshesPtr_->completeMesh(),
311 procMeshesPtr_->procMeshes(),
312 procMeshesPtr_->procFaceAddressing(),
313 procMeshesPtr_->procCellAddressing(),
314 procMeshesPtr_->procFaceAddressingBf()
321 ->reconstructVolInternalField<Type>(*iter());
326 new VolInternalField<Type>
329 procMeshesPtr_->completeMesh()
332 const VolInternalField<Type>&
tf = ttf();
335 convertVolInternalFieldBlock<Type>
344 convertVolInternalFieldBlock<Type>
348 arrayRangeCellZones_,
353 convertVolInternalFieldBlock<Type>
365 void Foam::vtkPVFoam::convertVolFieldBlock
367 const VolField<Type>&
tf,
368 autoPtr<PointField<Type>>& ptfPtr,
369 vtkMultiBlockDataSet* output,
370 const arrayRange&
range,
371 const List<polyDecomp>& decompLst
374 for (
int partId =
range.start(); partId <
range.end(); ++partId)
376 const label datasetNo = partDataset_[partId];
378 if (datasetNo >= 0 && partStatus_[partId])
380 convertVolInternalField<Type>
407 void Foam::vtkPVFoam::convertVolInternalFieldBlock
410 vtkMultiBlockDataSet* output,
411 const arrayRange&
range,
415 for (
int partId =
range.start(); partId <
range.end(); ++partId)
417 const label datasetNo = partDataset_[partId];
419 if (datasetNo >= 0 && partStatus_[partId])
421 convertVolInternalField<Type>
435 void Foam::vtkPVFoam::convertVolInternalField
438 vtkMultiBlockDataSet* output,
439 const arrayRange&
range,
440 const label datasetNo,
441 const polyDecomp& decompInfo
445 const labelList& superCells = decompInfo.superCells();
448 celldata->SetNumberOfTuples(superCells.
size());
449 celldata->SetNumberOfComponents(nComp);
450 celldata->Allocate(nComp*superCells.
size());
451 celldata->SetName(
tf.name().c_str());
454 <<
"Converting Vol field: " <<
tf.name()
455 <<
" size=" <<
tf.
size() <<
" (" << superCells.
size()
456 <<
"), nComp=" << nComp <<
endl;
461 const Type& t =
tf[superCells[i]];
466 vtkOpenFOAMTupleRemap<Type>(vec);
468 celldata->InsertTuple(i, vec);
471 vtkUnstructuredGrid::SafeDownCast
473 GetDataSetFromBlock(output,
range, datasetNo)
475 ->AddArray(celldata);
#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.
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.
Traits class for primitives.
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
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.
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)
typename VolField< Type >::Internal VolInternalField