ParaFoam Reader: Corrected ordering of the components of the symmetric tensor to...
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamPointFields.H
blob254635802a4ec69c1c44206b43ed8dcecbaa4f5f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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 #include "vtkOpenFOAMTupleRemap.H"
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 template<class Type>
41 void Foam::vtkPV3Foam::convertPointFields
43     const fvMesh& mesh,
44     const pointMesh& pMesh,
45     const IOobjectList& objects,
46     vtkMultiBlockDataSet* output
49     const polyBoundaryMesh& patches = mesh.boundaryMesh();
51     forAllConstIter(IOobjectList, objects, iter)
52     {
53         const word& fieldName = iter()->name();
54         // restrict to this GeometricField<Type, ...>
55         if
56         (
57             iter()->headerClassName()
58          != GeometricField<Type, pointPatchField, pointMesh>::typeName
59         )
60         {
61             continue;
62         }
64         if (debug)
65         {
66             Info<< "Foam::vtkPV3Foam::convertPointFields : "
67                 << fieldName << endl;
68         }
70         GeometricField<Type, pointPatchField, pointMesh> ptf
71         (
72             *iter(),
73             pMesh
74         );
77         // Convert activated internalMesh regions
78         convertPointFieldBlock
79         (
80             ptf,
81             output,
82             partInfoVolume_,
83             regionPolyDecomp_
84         );
86         // Convert activated cellZones
87         convertPointFieldBlock
88         (
89             ptf,
90             output,
91             partInfoCellZones_,
92             zonePolyDecomp_
93         );
95         // Convert activated cellSets
96         convertPointFieldBlock
97         (
98             ptf,
99             output,
100             partInfoCellSets_,
101             csetPolyDecomp_
102         );
105         //
106         // Convert patches - if activated
107         //
108         for
109         (
110             int partId = partInfoPatches_.start();
111             partId < partInfoPatches_.end();
112             ++partId
113         )
114         {
115             const word  patchName = getPartName(partId);
116             const label datasetNo = partDataset_[partId];
117             const label   patchId = patches.findPatchID(patchName);
119             if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
120             {
121                 continue;
122             }
124             convertPatchPointField
125             (
126                 fieldName,
127                 ptf.boundaryField()[patchId].patchInternalField()(),
128                 output,
129                 partInfoPatches_,
130                 datasetNo
131             );
132         }
133     }
137 template<class Type>
138 void Foam::vtkPV3Foam::convertPointFieldBlock
140     const GeometricField<Type, pointPatchField, pointMesh>& ptf,
141     vtkMultiBlockDataSet* output,
142     const partInfo& selector,
143     const List<polyDecomp>& decompLst
146    for (int partId = selector.start(); partId < selector.end(); ++partId)
147    {
148        const label datasetNo = partDataset_[partId];
150        if (datasetNo >= 0 && partStatus_[partId])
151        {
152            convertPointField
153            (
154                ptf,
155                GeometricField<Type, fvPatchField, volMesh>::null(),
156                output,
157                selector,
158                datasetNo,
159                decompLst[datasetNo]
160            );
161        }
162    }
166 template<class Type>
167 void Foam::vtkPV3Foam::convertPointField
169     const GeometricField<Type, pointPatchField, pointMesh>& ptf,
170     const GeometricField<Type, fvPatchField, volMesh>& tf,
171     vtkMultiBlockDataSet* output,
172     const partInfo& selector,
173     const label datasetNo,
174     const polyDecomp& decomp
177     const label nComp = pTraits<Type>::nComponents;
178     const labelList& addPointCellLabels = decomp.addPointCellLabels();
179     const labelList& pointMap = decomp.pointMap();
181     // use a pointMap or address directly into mesh
182     label nPoints;
183     if (pointMap.size())
184     {
185         nPoints = pointMap.size();
186     }
187     else
188     {
189         nPoints = ptf.size();
190     }
192     vtkFloatArray *pointData = vtkFloatArray::New();
193     pointData->SetNumberOfTuples(nPoints + addPointCellLabels.size());
194     pointData->SetNumberOfComponents(nComp);
195     pointData->Allocate(nComp*(nPoints + addPointCellLabels.size()));
196     pointData->SetName(ptf.name().c_str());
198     if (debug)
199     {
200         Info<< "convert convertPointField: "
201             << ptf.name()
202             << " size = " << nPoints
203             << " nComp=" << nComp
204             << " nTuples = " << (nPoints + addPointCellLabels.size())
205             <<  endl;
206     }
208     float vec[nComp];
210     if (pointMap.size())
211     {
212         forAll(pointMap, i)
213         {
214             const Type& t = ptf[pointMap[i]];
215             for (direction d=0; d<nComp; d++)
216             {
217                 vec[d] = component(t, d);
218             }
219             vtkOpenFOAMTupleRemap<Type>(vec);
221             pointData->InsertTuple(i, vec);
222         }
223     }
224     else
225     {
226         forAll(ptf, i)
227         {
228             const Type& t = ptf[i];
229             for (direction d=0; d<nComp; d++)
230             {
231                 vec[d] = component(t, d);
232             }
233             vtkOpenFOAMTupleRemap<Type>(vec);
235             pointData->InsertTuple(i, vec);
236         }
237     }
239     // continue insertion from here
240     label i = nPoints;
242     if (&tf != &GeometricField<Type, fvPatchField, volMesh>::null())
243     {
244         forAll(addPointCellLabels, apI)
245         {
246             const Type& t = tf[addPointCellLabels[apI]];
247             for (direction d=0; d<nComp; d++)
248             {
249                 vec[d] = component(t, d);
250             }
251             vtkOpenFOAMTupleRemap<Type>(vec);
253             pointData->InsertTuple(i++, vec);
254         }
255     }
256     else
257     {
258         forAll(addPointCellLabels, apI)
259         {
260             Type t = interpolatePointToCell(ptf, addPointCellLabels[apI]);
261             for (direction d=0; d<nComp; d++)
262             {
263                 vec[d] = component(t, d);
264             }
265             vtkOpenFOAMTupleRemap<Type>(vec);
267             pointData->InsertTuple(i++, vec);
268         }
269     }
271     vtkUnstructuredGrid::SafeDownCast
272     (
273         GetDataSetFromBlock(output, selector, datasetNo)
274     )   ->GetPointData()
275         ->AddArray(pointData);
277     pointData->Delete();
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 #endif
285 // ************************************************************************* //