bugfix: use ptf.name() on point fields
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamPointFields.H
blob3a609b2bde50f28d1946ae23fb598175d5720046
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 vtkPV3FoamPointFields_H
31 #define vtkPV3FoamPointFields_H
33 // Foam includes
34 #include "interpolatePointToCell.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 template<class Type>
39 void Foam::vtkPV3Foam::convertPointFields
41     const fvMesh& mesh,
42     const pointMesh& pMesh,
43     const IOobjectList& objects,
44     vtkMultiBlockDataSet* output
47     const polyBoundaryMesh& patches = mesh.boundaryMesh();
49     forAllConstIter(IOobjectList, objects, iter)
50     {
51         const word& fieldName = iter()->name();
52         // restrict to this GeometricField<Type, ...>
53         if
54         (
55             iter()->headerClassName()
56          != GeometricField<Type, pointPatchField, pointMesh>::typeName
57         )
58         {
59             continue;
60         }
62         if (debug)
63         {
64             Info<< "Foam::vtkPV3Foam::convertPointFields : "
65                 << fieldName << endl;
66         }
68         GeometricField<Type, pointPatchField, pointMesh> ptf
69         (
70             *iter(),
71             pMesh
72         );
75         // Convert activated internalMesh regions
76         convertPointFieldBlock
77         (
78             ptf,
79             output,
80             partInfoVolume_,
81             regionPolyDecomp_
82         );
84         // Convert activated cellZones
85         convertPointFieldBlock
86         (
87             ptf,
88             output,
89             partInfoCellZones_,
90             zonePolyDecomp_
91         );
93         // Convert activated cellSets
94         convertPointFieldBlock
95         (
96             ptf,
97             output,
98             partInfoCellSets_,
99             csetPolyDecomp_
100         );
103         //
104         // Convert patches - if activated
105         //
106         for
107         (
108             int partId = partInfoPatches_.start();
109             partId < partInfoPatches_.end();
110             ++partId
111         )
112         {
113             const word  patchName = getPartName(partId);
114             const label datasetNo = partDataset_[partId];
115             const label   patchId = patches.findPatchID(patchName);
117             if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
118             {
119                 continue;
120             }
122             convertPatchPointField
123             (
124                 fieldName,
125                 ptf.boundaryField()[patchId].patchInternalField()(),
126                 output,
127                 partInfoPatches_,
128                 datasetNo
129             );
130         }
131     }
135 template<class Type>
136 void Foam::vtkPV3Foam::convertPointFieldBlock
138     const GeometricField<Type, pointPatchField, pointMesh>& ptf,
139     vtkMultiBlockDataSet* output,
140     const partInfo& selector,
141     const List<polyDecomp>& decompLst
144    for (int partId = selector.start(); partId < selector.end(); ++partId)
145    {
146        const label datasetNo = partDataset_[partId];
148        if (datasetNo >= 0 && partStatus_[partId])
149        {
150            convertPointField
151            (
152                ptf,
153                GeometricField<Type, fvPatchField, volMesh>::null(),
154                output,
155                selector,
156                datasetNo,
157                decompLst[datasetNo]
158            );
159        }
160    }
164 template<class Type>
165 void Foam::vtkPV3Foam::convertPointField
167     const GeometricField<Type, pointPatchField, pointMesh>& ptf,
168     const GeometricField<Type, fvPatchField, volMesh>& tf,
169     vtkMultiBlockDataSet* output,
170     const partInfo& selector,
171     const label datasetNo,
172     const polyDecomp& decomp
175     const label nComp = pTraits<Type>::nComponents;
176     const labelList& addPointCellLabels = decomp.addPointCellLabels();
177     const labelList& pointMap = decomp.pointMap();
179     // use a pointMap or address directly into mesh
180     label nPoints;
181     if (pointMap.size())
182     {
183         nPoints = pointMap.size();
184     }
185     else
186     {
187         nPoints = ptf.size();
188     }
190     vtkFloatArray *pointData = vtkFloatArray::New();
191     pointData->SetNumberOfTuples(nPoints + addPointCellLabels.size());
192     pointData->SetNumberOfComponents(nComp);
193     pointData->Allocate(nComp*(nPoints + addPointCellLabels.size()));
194     pointData->SetName(ptf.name().c_str());
196     if (debug)
197     {
198         Info<< "convert convertPointField: "
199             << ptf.name()
200             << " size = " << nPoints
201             << " nComp=" << nComp
202             << " nTuples = " << (nPoints + addPointCellLabels.size())
203             <<  endl;
204     }
206     float vec[nComp];
208     if (pointMap.size())
209     {
210         forAll(pointMap, i)
211         {
212             const Type& t = ptf[pointMap[i]];
213             for (direction d=0; d<nComp; d++)
214             {
215                 vec[d] = component(t, d);
216             }
217             pointData->InsertTuple(i, vec);
218         }
219     }
220     else
221     {
222         forAll(ptf, i)
223         {
224             const Type& t = ptf[i];
225             for (direction d=0; d<nComp; d++)
226             {
227                 vec[d] = component(t, d);
228             }
229             pointData->InsertTuple(i, vec);
230         }
231     }
233     // continue insertion from here
234     label i = nPoints;
236     if (&tf != &GeometricField<Type, fvPatchField, volMesh>::null())
237     {
238         forAll(addPointCellLabels, apI)
239         {
240             const Type& t = tf[addPointCellLabels[apI]];
241             for (direction d=0; d<nComp; d++)
242             {
243                 vec[d] = component(t, d);
244             }
245             pointData->InsertTuple(i++, vec);
246         }
247     }
248     else
249     {
250         forAll(addPointCellLabels, apI)
251         {
252             Type t = interpolatePointToCell(ptf, addPointCellLabels[apI]);
253             for (direction d=0; d<nComp; d++)
254             {
255                 vec[d] = component(t, d);
256             }
257             pointData->InsertTuple(i++, vec);
258         }
259     }
261     vtkUnstructuredGrid::SafeDownCast
262     (
263         GetDataSetFromBlock(output, selector, datasetNo)
264     )   ->GetPointData()
265         ->AddArray(pointData);
267     pointData->Delete();
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 #endif
275 // ************************************************************************* //