2010-08-26 14:22:03 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
2013-12-11 16:09:41 +00:00
|
|
|
\\ / F ield | foam-extend: Open Source CFD
|
2016-06-20 15:00:40 +00:00
|
|
|
\\ / O peration | Version: 4.0
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
2010-08-26 14:22:03 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2013-12-11 16:09:41 +00:00
|
|
|
This file is part of foam-extend.
|
2010-08-26 14:22:03 +00:00
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2010-08-26 14:22:03 +00:00
|
|
|
under the terms of the GNU General Public License as published by the
|
2013-12-11 16:09:41 +00:00
|
|
|
Free Software Foundation, either version 3 of the License, or (at your
|
2010-08-26 14:22:03 +00:00
|
|
|
option) any later version.
|
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is distributed in the hope that it will be useful, but
|
|
|
|
WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
|
|
|
|
General Public License for more details.
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
2013-12-11 16:09:41 +00:00
|
|
|
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
InClass
|
|
|
|
vtkPV3Foam
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#ifndef vtkPV3FoamPointFields_H
|
|
|
|
#define vtkPV3FoamPointFields_H
|
|
|
|
|
|
|
|
// Foam includes
|
|
|
|
#include "interpolatePointToCell.H"
|
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
#include "vtkFOAMTupleRemap.H"
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
template<class Type>
|
|
|
|
void Foam::vtkPV3Foam::convertPointFields
|
|
|
|
(
|
|
|
|
const fvMesh& mesh,
|
|
|
|
const pointMesh& pMesh,
|
|
|
|
const IOobjectList& objects,
|
|
|
|
vtkMultiBlockDataSet* output
|
|
|
|
)
|
|
|
|
{
|
|
|
|
const polyBoundaryMesh& patches = mesh.boundaryMesh();
|
|
|
|
|
|
|
|
forAllConstIter(IOobjectList, objects, iter)
|
|
|
|
{
|
|
|
|
const word& fieldName = iter()->name();
|
|
|
|
// restrict to this GeometricField<Type, ...>
|
|
|
|
if
|
|
|
|
(
|
|
|
|
iter()->headerClassName()
|
|
|
|
!= GeometricField<Type, pointPatchField, pointMesh>::typeName
|
|
|
|
)
|
|
|
|
{
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
|
|
|
|
if (debug)
|
|
|
|
{
|
|
|
|
Info<< "Foam::vtkPV3Foam::convertPointFields : "
|
|
|
|
<< fieldName << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
GeometricField<Type, pointPatchField, pointMesh> ptf
|
|
|
|
(
|
|
|
|
*iter(),
|
|
|
|
pMesh
|
|
|
|
);
|
|
|
|
|
|
|
|
|
|
|
|
// Convert activated internalMesh regions
|
|
|
|
convertPointFieldBlock
|
|
|
|
(
|
|
|
|
ptf,
|
|
|
|
output,
|
|
|
|
partInfoVolume_,
|
|
|
|
regionPolyDecomp_
|
|
|
|
);
|
|
|
|
|
|
|
|
// Convert activated cellZones
|
|
|
|
convertPointFieldBlock
|
|
|
|
(
|
|
|
|
ptf,
|
|
|
|
output,
|
|
|
|
partInfoCellZones_,
|
|
|
|
zonePolyDecomp_
|
|
|
|
);
|
|
|
|
|
|
|
|
// Convert activated cellSets
|
|
|
|
convertPointFieldBlock
|
|
|
|
(
|
|
|
|
ptf,
|
|
|
|
output,
|
|
|
|
partInfoCellSets_,
|
|
|
|
csetPolyDecomp_
|
|
|
|
);
|
|
|
|
|
|
|
|
|
|
|
|
//
|
|
|
|
// Convert patches - if activated
|
|
|
|
//
|
|
|
|
for
|
|
|
|
(
|
|
|
|
int partId = partInfoPatches_.start();
|
|
|
|
partId < partInfoPatches_.end();
|
|
|
|
++partId
|
|
|
|
)
|
|
|
|
{
|
|
|
|
const word patchName = getPartName(partId);
|
|
|
|
const label datasetNo = partDataset_[partId];
|
|
|
|
const label patchId = patches.findPatchID(patchName);
|
|
|
|
|
|
|
|
if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
|
|
|
|
{
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
|
|
|
|
convertPatchPointField
|
|
|
|
(
|
|
|
|
fieldName,
|
|
|
|
ptf.boundaryField()[patchId].patchInternalField()(),
|
|
|
|
output,
|
|
|
|
partInfoPatches_,
|
|
|
|
datasetNo
|
|
|
|
);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
template<class Type>
|
|
|
|
void Foam::vtkPV3Foam::convertPointFieldBlock
|
|
|
|
(
|
|
|
|
const GeometricField<Type, pointPatchField, pointMesh>& ptf,
|
|
|
|
vtkMultiBlockDataSet* output,
|
|
|
|
const partInfo& selector,
|
|
|
|
const List<polyDecomp>& decompLst
|
|
|
|
)
|
|
|
|
{
|
|
|
|
for (int partId = selector.start(); partId < selector.end(); ++partId)
|
|
|
|
{
|
|
|
|
const label datasetNo = partDataset_[partId];
|
|
|
|
|
|
|
|
if (datasetNo >= 0 && partStatus_[partId])
|
|
|
|
{
|
|
|
|
convertPointField
|
|
|
|
(
|
|
|
|
ptf,
|
|
|
|
GeometricField<Type, fvPatchField, volMesh>::null(),
|
|
|
|
output,
|
|
|
|
selector,
|
|
|
|
datasetNo,
|
|
|
|
decompLst[datasetNo]
|
|
|
|
);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
template<class Type>
|
|
|
|
void Foam::vtkPV3Foam::convertPointField
|
|
|
|
(
|
|
|
|
const GeometricField<Type, pointPatchField, pointMesh>& ptf,
|
|
|
|
const GeometricField<Type, fvPatchField, volMesh>& tf,
|
|
|
|
vtkMultiBlockDataSet* output,
|
|
|
|
const partInfo& selector,
|
|
|
|
const label datasetNo,
|
|
|
|
const polyDecomp& decomp
|
|
|
|
)
|
|
|
|
{
|
|
|
|
const label nComp = pTraits<Type>::nComponents;
|
|
|
|
const labelList& addPointCellLabels = decomp.addPointCellLabels();
|
|
|
|
const labelList& pointMap = decomp.pointMap();
|
|
|
|
|
|
|
|
// use a pointMap or address directly into mesh
|
|
|
|
label nPoints;
|
|
|
|
if (pointMap.size())
|
|
|
|
{
|
|
|
|
nPoints = pointMap.size();
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
nPoints = ptf.size();
|
|
|
|
}
|
|
|
|
|
|
|
|
vtkFloatArray *pointData = vtkFloatArray::New();
|
|
|
|
pointData->SetNumberOfTuples(nPoints + addPointCellLabels.size());
|
|
|
|
pointData->SetNumberOfComponents(nComp);
|
|
|
|
pointData->Allocate(nComp*(nPoints + addPointCellLabels.size()));
|
|
|
|
pointData->SetName(ptf.name().c_str());
|
|
|
|
|
|
|
|
if (debug)
|
|
|
|
{
|
|
|
|
Info<< "convert convertPointField: "
|
|
|
|
<< ptf.name()
|
|
|
|
<< " size = " << nPoints
|
|
|
|
<< " nComp=" << nComp
|
|
|
|
<< " nTuples = " << (nPoints + addPointCellLabels.size())
|
|
|
|
<< endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
float vec[nComp];
|
|
|
|
|
|
|
|
if (pointMap.size())
|
|
|
|
{
|
|
|
|
forAll(pointMap, i)
|
|
|
|
{
|
|
|
|
const Type& t = ptf[pointMap[i]];
|
|
|
|
for (direction d=0; d<nComp; d++)
|
|
|
|
{
|
|
|
|
vec[d] = component(t, d);
|
|
|
|
}
|
2013-12-11 16:09:41 +00:00
|
|
|
vtkFOAMTupleRemap<Type>(vec);
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
pointData->InsertTuple(i, vec);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
forAll(ptf, i)
|
|
|
|
{
|
|
|
|
const Type& t = ptf[i];
|
|
|
|
for (direction d=0; d<nComp; d++)
|
|
|
|
{
|
|
|
|
vec[d] = component(t, d);
|
|
|
|
}
|
2013-12-11 16:09:41 +00:00
|
|
|
vtkFOAMTupleRemap<Type>(vec);
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
pointData->InsertTuple(i, vec);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// continue insertion from here
|
|
|
|
label i = nPoints;
|
|
|
|
|
|
|
|
if (&tf != &GeometricField<Type, fvPatchField, volMesh>::null())
|
|
|
|
{
|
|
|
|
forAll(addPointCellLabels, apI)
|
|
|
|
{
|
|
|
|
const Type& t = tf[addPointCellLabels[apI]];
|
|
|
|
for (direction d=0; d<nComp; d++)
|
|
|
|
{
|
|
|
|
vec[d] = component(t, d);
|
|
|
|
}
|
2013-12-11 16:09:41 +00:00
|
|
|
vtkFOAMTupleRemap<Type>(vec);
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
pointData->InsertTuple(i++, vec);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
forAll(addPointCellLabels, apI)
|
|
|
|
{
|
|
|
|
Type t = interpolatePointToCell(ptf, addPointCellLabels[apI]);
|
|
|
|
for (direction d=0; d<nComp; d++)
|
|
|
|
{
|
|
|
|
vec[d] = component(t, d);
|
|
|
|
}
|
2013-12-11 16:09:41 +00:00
|
|
|
vtkFOAMTupleRemap<Type>(vec);
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
pointData->InsertTuple(i++, vec);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
vtkUnstructuredGrid::SafeDownCast
|
|
|
|
(
|
|
|
|
GetDataSetFromBlock(output, selector, datasetNo)
|
|
|
|
) ->GetPointData()
|
|
|
|
->AddArray(pointData);
|
|
|
|
|
|
|
|
pointData->Delete();
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
#endif
|
|
|
|
|
|
|
|
// ************************************************************************* //
|