1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ 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 3 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, see <http://www.gnu.org/licenses/>.
InClass
vtkPVFoam
\*---------------------------------------------------------------------------*/
#ifndef vtkPVFoamLagrangianFields_H
#define vtkPVFoamLagrangianFields_H
#include "Cloud.H"
#include "vtkOpenFOAMTupleRemap.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type>
void Foam::vtkPVFoam::convertLagrangianFields
(
const IOobjectList& objects,
vtkMultiBlockDataSet* output,
const label datasetNo
)
{
const arrayRange& range = arrayRangeLagrangian_;
forAllConstIter(IOobjectList, objects, iter)
{
// restrict to this IOField<Type>
if (iter()->headerClassName() == IOField<Type>::typeName)
{
IOField<Type> tf(*iter());
convertLagrangianField(tf, output, range, datasetNo);
}
}
}
template<class Type>
void Foam::vtkPVFoam::convertLagrangianField
(
const IOField<Type>& tf,
vtkMultiBlockDataSet* output,
const arrayRange& range,
const label datasetNo
)
{
const label nComp = pTraits<Type>::nComponents;
vtkFloatArray* pointData = vtkFloatArray::New();
pointData->SetNumberOfTuples(tf.size());
pointData->SetNumberOfComponents(nComp);
pointData->Allocate(nComp*tf.size());
pointData->SetName(tf.name().c_str());
if (debug)
{
Info<< "convert LagrangianField: "
<< tf.name()
<< " size = " << tf.size()
<< " nComp=" << nComp
<< " nTuples = " << tf.size() << endl;
}
float vec[nComp];
forAll(tf, i)
{
const Type& t = tf[i];
for (direction d=0; d<nComp; ++d)
{
vec[d] = component(t, d);
}
vtkOpenFOAMTupleRemap<Type>(vec);
pointData->InsertTuple(i, vec);
}
vtkPolyData::SafeDownCast
(
GetDataSetFromBlock(output, range, datasetNo)
) ->GetPointData()
->AddArray(pointData);
pointData->Delete();
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
|