This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/utilities/postProcessing/graphics/PV3FoamReader/vtkPV3Foam/vtkPV3FoamLagrangianFields.H

115 lines
3.2 KiB
C
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright held by original author
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM 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.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
InClass
vtkPV3Foam
\*---------------------------------------------------------------------------*/
2010-08-26 14:22:03 +00:00
#ifndef vtkPV3FoamLagrangianFields_H
#define vtkPV3FoamLagrangianFields_H
2010-08-26 14:22:03 +00:00
#include "Cloud.H"
#include "vtkOpenFOAMTupleRemap.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type>
2010-08-26 14:22:03 +00:00
void Foam::vtkPV3Foam::convertLagrangianFields
(
2010-08-26 14:22:03 +00:00
const IOobjectList& objects,
vtkMultiBlockDataSet* output,
const label datasetNo
)
{
2010-08-26 14:22:03 +00:00
const partInfo& selector = partInfoLagrangian_;
2010-08-26 14:22:03 +00:00
forAllConstIter(IOobjectList, objects, iter)
{
// restrict to this IOField<Type>
if (iter()->headerClassName() == IOField<Type>::typeName)
{
IOField<Type> tf(*iter());
convertLagrangianField(tf, output, selector, datasetNo);
}
}
}
template<class Type>
void Foam::vtkPV3Foam::convertLagrangianField
(
const IOField<Type>& tf,
vtkMultiBlockDataSet* output,
const partInfo& selector,
const label datasetNo
)
{
const label nComp = pTraits<Type>::nComponents;
vtkFloatArray *pointData = vtkFloatArray::New();
2010-08-26 14:22:03 +00:00
pointData->SetNumberOfTuples( tf.size() );
pointData->SetNumberOfComponents( nComp );
pointData->Allocate( nComp*tf.size() );
pointData->SetName( tf.name().c_str() );
2010-08-26 14:22:03 +00:00
if (debug)
{
Info<< "convert LagrangianField: "
<< tf.name()
<< " size = " << tf.size()
<< " nComp=" << nComp
<< " nTuples = " << tf.size() << endl;
}
2010-08-26 14:22:03 +00:00
float vec[nComp];
forAll(tf, i)
{
2010-08-26 14:22:03 +00:00
const Type& t = tf[i];
for (direction d=0; d<nComp; d++)
{
2010-08-26 14:22:03 +00:00
vec[d] = component(t, d);
}
2010-08-26 14:22:03 +00:00
vtkOpenFOAMTupleRemap<Type>(vec);
pointData->InsertTuple(i, vec);
}
2010-08-26 14:22:03 +00:00
vtkPolyData::SafeDownCast
(
GetDataSetFromBlock(output, selector, datasetNo)
) ->GetPointData()
->AddArray(pointData);
pointData->Delete();
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //