Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamPointFields.H
blob4adff917ff334b4db271c5b13759ae4619c43395
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 InClass
25     vtkPV3Foam
27 \*---------------------------------------------------------------------------*/
29 #ifndef vtkPV3FoamPointFields_H
30 #define vtkPV3FoamPointFields_H
32 // Foam includes
33 #include "interpolatePointToCell.H"
35 #include "vtkFOAMTupleRemap.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 template<class Type>
40 void Foam::vtkPV3Foam::convertPointFields
42     const fvMesh& mesh,
43     const pointMesh& pMesh,
44     const IOobjectList& objects,
45     vtkMultiBlockDataSet* output
48     const polyBoundaryMesh& patches = mesh.boundaryMesh();
50     forAllConstIter(IOobjectList, objects, iter)
51     {
52         const word& fieldName = iter()->name();
53         // restrict to this GeometricField<Type, ...>
54         if
55         (
56             iter()->headerClassName()
57          != GeometricField<Type, pointPatchField, pointMesh>::typeName
58         )
59         {
60             continue;
61         }
63         if (debug)
64         {
65             Info<< "Foam::vtkPV3Foam::convertPointFields : "
66                 << fieldName << endl;
67         }
69         GeometricField<Type, pointPatchField, pointMesh> ptf
70         (
71             *iter(),
72             pMesh
73         );
76         // Convert activated internalMesh regions
77         convertPointFieldBlock
78         (
79             ptf,
80             output,
81             partInfoVolume_,
82             regionPolyDecomp_
83         );
85         // Convert activated cellZones
86         convertPointFieldBlock
87         (
88             ptf,
89             output,
90             partInfoCellZones_,
91             zonePolyDecomp_
92         );
94         // Convert activated cellSets
95         convertPointFieldBlock
96         (
97             ptf,
98             output,
99             partInfoCellSets_,
100             csetPolyDecomp_
101         );
104         //
105         // Convert patches - if activated
106         //
107         for
108         (
109             int partId = partInfoPatches_.start();
110             partId < partInfoPatches_.end();
111             ++partId
112         )
113         {
114             const word  patchName = getPartName(partId);
115             const label datasetNo = partDataset_[partId];
116             const label   patchId = patches.findPatchID(patchName);
118             if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
119             {
120                 continue;
121             }
123             convertPatchPointField
124             (
125                 fieldName,
126                 ptf.boundaryField()[patchId].patchInternalField()(),
127                 output,
128                 partInfoPatches_,
129                 datasetNo
130             );
131         }
132     }
136 template<class Type>
137 void Foam::vtkPV3Foam::convertPointFieldBlock
139     const GeometricField<Type, pointPatchField, pointMesh>& ptf,
140     vtkMultiBlockDataSet* output,
141     const partInfo& selector,
142     const List<polyDecomp>& decompLst
145    for (int partId = selector.start(); partId < selector.end(); ++partId)
146    {
147        const label datasetNo = partDataset_[partId];
149        if (datasetNo >= 0 && partStatus_[partId])
150        {
151            convertPointField
152            (
153                ptf,
154                GeometricField<Type, fvPatchField, volMesh>::null(),
155                output,
156                selector,
157                datasetNo,
158                decompLst[datasetNo]
159            );
160        }
161    }
165 template<class Type>
166 void Foam::vtkPV3Foam::convertPointField
168     const GeometricField<Type, pointPatchField, pointMesh>& ptf,
169     const GeometricField<Type, fvPatchField, volMesh>& tf,
170     vtkMultiBlockDataSet* output,
171     const partInfo& selector,
172     const label datasetNo,
173     const polyDecomp& decomp
176     const label nComp = pTraits<Type>::nComponents;
177     const labelList& addPointCellLabels = decomp.addPointCellLabels();
178     const labelList& pointMap = decomp.pointMap();
180     // use a pointMap or address directly into mesh
181     label nPoints;
182     if (pointMap.size())
183     {
184         nPoints = pointMap.size();
185     }
186     else
187     {
188         nPoints = ptf.size();
189     }
191     vtkFloatArray *pointData = vtkFloatArray::New();
192     pointData->SetNumberOfTuples(nPoints + addPointCellLabels.size());
193     pointData->SetNumberOfComponents(nComp);
194     pointData->Allocate(nComp*(nPoints + addPointCellLabels.size()));
195     pointData->SetName(ptf.name().c_str());
197     if (debug)
198     {
199         Info<< "convert convertPointField: "
200             << ptf.name()
201             << " size = " << nPoints
202             << " nComp=" << nComp
203             << " nTuples = " << (nPoints + addPointCellLabels.size())
204             <<  endl;
205     }
207     float vec[nComp];
209     if (pointMap.size())
210     {
211         forAll(pointMap, i)
212         {
213             const Type& t = ptf[pointMap[i]];
214             for (direction d=0; d<nComp; d++)
215             {
216                 vec[d] = component(t, d);
217             }
218             vtkFOAMTupleRemap<Type>(vec);
220             pointData->InsertTuple(i, vec);
221         }
222     }
223     else
224     {
225         forAll(ptf, i)
226         {
227             const Type& t = ptf[i];
228             for (direction d=0; d<nComp; d++)
229             {
230                 vec[d] = component(t, d);
231             }
232             vtkFOAMTupleRemap<Type>(vec);
234             pointData->InsertTuple(i, vec);
235         }
236     }
238     // continue insertion from here
239     label i = nPoints;
241     if (&tf != &GeometricField<Type, fvPatchField, volMesh>::null())
242     {
243         forAll(addPointCellLabels, apI)
244         {
245             const Type& t = tf[addPointCellLabels[apI]];
246             for (direction d=0; d<nComp; d++)
247             {
248                 vec[d] = component(t, d);
249             }
250             vtkFOAMTupleRemap<Type>(vec);
252             pointData->InsertTuple(i++, vec);
253         }
254     }
255     else
256     {
257         forAll(addPointCellLabels, apI)
258         {
259             Type t = interpolatePointToCell(ptf, addPointCellLabels[apI]);
260             for (direction d=0; d<nComp; d++)
261             {
262                 vec[d] = component(t, d);
263             }
264             vtkFOAMTupleRemap<Type>(vec);
266             pointData->InsertTuple(i++, vec);
267         }
268     }
270     vtkUnstructuredGrid::SafeDownCast
271     (
272         GetDataSetFromBlock(output, selector, datasetNo)
273     )   ->GetPointData()
274         ->AddArray(pointData);
276     pointData->Delete();
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 #endif
284 // ************************************************************************* //