BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / graphics / PV3Readers / PV3FoamReader / vtkPV3Foam / vtkPV3FoamFaceField.H
blob030c9aca05a5407dfede8c02354fc7f5530ed384
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 InClass
25     vtkPV3Foam
27 \*---------------------------------------------------------------------------*/
29 #ifndef vtkPV3FoamFaceField_H
30 #define vtkPV3FoamFaceField_H
32 // VTK includes
33 #include "vtkCellData.h"
34 #include "vtkFloatArray.h"
35 #include "vtkMultiBlockDataSet.h"
36 #include "vtkPolyData.h"
38 #include "vtkOpenFOAMTupleRemap.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 template<class Type>
43 void Foam::vtkPV3Foam::convertFaceField
45     const GeometricField<Type, fvPatchField, volMesh>& tf,
46     vtkMultiBlockDataSet* output,
47     const arrayRange& range,
48     const label datasetNo,
49     const fvMesh& mesh,
50     const labelList& faceLabels
53     const label nComp = pTraits<Type>::nComponents;
54     const label nInternalFaces = mesh.nInternalFaces();
55     const labelList& faceOwner = mesh.faceOwner();
56     const labelList& faceNeigh = mesh.faceNeighbour();
58     vtkFloatArray* cellData = vtkFloatArray::New();
59     cellData->SetNumberOfTuples(faceLabels.size());
60     cellData->SetNumberOfComponents(nComp);
61     cellData->Allocate(nComp*faceLabels.size());
62     cellData->SetName(tf.name().c_str());
64     if (debug)
65     {
66         Info<< "convert convertFaceField: "
67             << tf.name()
68             << " size = " << tf.size()
69             << " nComp=" << nComp
70             << " nTuples = " << faceLabels.size() <<  endl;
71     }
73     float vec[nComp];
75     // for interior faces: average owner/neighbour
76     // for boundary faces: owner
77     forAll(faceLabels, faceI)
78     {
79         const label faceNo = faceLabels[faceI];
80         if (faceNo < nInternalFaces)
81         {
82             Type t = 0.5*(tf[faceOwner[faceNo]] + tf[faceNeigh[faceNo]]);
84             for (direction d=0; d<nComp; ++d)
85             {
86                 vec[d] = component(t, d);
87             }
88         }
89         else
90         {
91             const Type& t = tf[faceOwner[faceNo]];
92             for (direction d=0; d<nComp; ++d)
93             {
94                 vec[d] = component(t, d);
95             }
96         }
97         vtkOpenFOAMTupleRemap<Type>(vec);
99         cellData->InsertTuple(faceI, vec);
100     }
103     vtkPolyData::SafeDownCast
104     (
105         GetDataSetFromBlock(output, range, datasetNo)
106     )   ->GetCellData()
107         ->AddArray(cellData);
109     cellData->Delete();
113 template<class Type>
114 void Foam::vtkPV3Foam::convertFaceField
116     const GeometricField<Type, fvPatchField, volMesh>& tf,
117     vtkMultiBlockDataSet* output,
118     const arrayRange& range,
119     const label datasetNo,
120     const fvMesh& mesh,
121     const faceSet& fSet
124     const label nComp = pTraits<Type>::nComponents;
125     const label nInternalFaces = mesh.nInternalFaces();
126     const labelList& faceOwner = mesh.faceOwner();
127     const labelList& faceNeigh = mesh.faceNeighbour();
129     vtkFloatArray* cellData = vtkFloatArray::New();
130     cellData->SetNumberOfTuples(fSet.size());
131     cellData->SetNumberOfComponents(nComp);
132     cellData->Allocate(nComp*fSet.size());
133     cellData->SetName(tf.name().c_str());
135     if (debug)
136     {
137         Info<< "convert convertFaceField: "
138             << tf.name()
139             << " size = " << tf.size()
140             << " nComp=" << nComp
141             << " nTuples = " << fSet.size() <<  endl;
142     }
144     float vec[nComp];
146     // for interior faces: average owner/neighbour
147     // for boundary faces: owner
148     label faceI = 0;
149     forAllConstIter(faceSet, fSet, iter)
150     {
151         const label faceNo = iter.key();
153         if (faceNo < nInternalFaces)
154         {
155             Type t = 0.5*(tf[faceOwner[faceNo]] + tf[faceNeigh[faceNo]]);
157             for (direction d=0; d<nComp; ++d)
158             {
159                 vec[d] = component(t, d);
160             }
161         }
162         else
163         {
164             const Type& t = tf[faceOwner[faceNo]];
165             for (direction d=0; d<nComp; ++d)
166             {
167                 vec[d] = component(t, d);
168             }
169         }
170         vtkOpenFOAMTupleRemap<Type>(vec);
172         cellData->InsertTuple(faceI, vec);
173         ++faceI;
174     }
177     vtkPolyData::SafeDownCast
178     (
179         GetDataSetFromBlock(output, range, datasetNo)
180     )   ->GetCellData()
181         ->AddArray(cellData);
183     cellData->Delete();
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 #endif
190 // ************************************************************************* //