BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamLagrangianFields.H
blob7ae588b25cc3435d963e3311ce54dec40c919b7e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 InClass
26     vtkPV3Foam
28 \*---------------------------------------------------------------------------*/
30 #ifndef vtkPV3FoamLagrangianFields_H
31 #define vtkPV3FoamLagrangianFields_H
33 #include "Cloud.H"
35 #include "vtkOpenFOAMTupleRemap.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 template<class Type>
40 void Foam::vtkPV3Foam::convertLagrangianFields
42     const IOobjectList& objects,
43     vtkMultiBlockDataSet* output,
44     const label datasetNo
47     const partInfo& selector = partInfoLagrangian_;
49     forAllConstIter(IOobjectList, objects, iter)
50     {
51         // restrict to this IOField<Type>
52         if (iter()->headerClassName() == IOField<Type>::typeName)
53         {
54             IOField<Type> tf(*iter());
55             convertLagrangianField(tf, output, selector, datasetNo);
56         }
57     }
61 template<class Type>
62 void Foam::vtkPV3Foam::convertLagrangianField
64     const IOField<Type>& tf,
65     vtkMultiBlockDataSet* output,
66     const partInfo& selector,
67     const label datasetNo
70     const label nComp = pTraits<Type>::nComponents;
72     vtkFloatArray *pointData = vtkFloatArray::New();
73     pointData->SetNumberOfTuples( tf.size() );
74     pointData->SetNumberOfComponents( nComp );
75     pointData->Allocate( nComp*tf.size() );
76     pointData->SetName( tf.name().c_str() );
78     if (debug)
79     {
80         Info<< "convert LagrangianField: "
81             << tf.name()
82             << " size = " << tf.size()
83             << " nComp=" << nComp
84             << " nTuples = " << tf.size() <<  endl;
85     }
87     float vec[nComp];
88     forAll(tf, i)
89     {
90         const Type& t = tf[i];
91         for (direction d=0; d<nComp; d++)
92         {
93             vec[d] = component(t, d);
94         }
95         vtkOpenFOAMTupleRemap<Type>(vec);
97         pointData->InsertTuple(i, vec);
98     }
101     vtkPolyData::SafeDownCast
102     (
103         GetDataSetFromBlock(output, selector, datasetNo)
104     )   ->GetPointData()
105         ->AddArray(pointData);
107     pointData->Delete();
110 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 #endif
114 // ************************************************************************* //