BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamVolFields.H
blob35e4093c89721d1b9adcea573ea7094fdbd5cd23
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 vtkPV3FoamVolFields_H
31 #define vtkPV3FoamVolFields_H
33 // Foam includes
34 #include "emptyFvPatchField.H"
35 #include "wallPolyPatch.H"
36 #include "faceSet.H"
37 #include "volPointInterpolation.H"
39 #include "vtkPV3FoamFaceField.H"
40 #include "vtkPV3FoamPatchField.H"
42 #include "vtkOpenFOAMTupleRemap.H"
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 template<class Type>
47 void Foam::vtkPV3Foam::convertVolFields
49     const fvMesh& mesh,
50     const PtrList<PrimitivePatchInterpolation<primitivePatch> >& ppInterpList,
51     const IOobjectList& objects,
52     vtkMultiBlockDataSet* output
55     const polyBoundaryMesh& patches = mesh.boundaryMesh();
57     forAllConstIter(IOobjectList, objects, iter)
58     {
59         // restrict to GeometricField<Type, ...>
60         if
61         (
62             iter()->headerClassName()
63          != GeometricField<Type, fvPatchField, volMesh>::typeName
64         )
65         {
66             continue;
67         }
69         // Load field
70         GeometricField<Type, fvPatchField, volMesh> tf
71         (
72             *iter(),
73             mesh
74         );
76         // Interpolated field (demand driven)
77         autoPtr<GeometricField<Type, pointPatchField, pointMesh> > ptfPtr;
80         // Convert activated internalMesh regions
81         convertVolFieldBlock
82         (
83             tf,
84             ptfPtr,
85             output,
86             partInfoVolume_,
87             regionPolyDecomp_
88         );
90         // Convert activated cellZones
91         convertVolFieldBlock
92         (
93             tf,
94             ptfPtr,
95             output,
96             partInfoCellZones_,
97             zonePolyDecomp_
98         );
100         // Convert activated cellSets
101         convertVolFieldBlock
102         (
103             tf,
104             ptfPtr,
105             output,
106             partInfoCellSets_,
107             csetPolyDecomp_
108         );
111         //
112         // Convert patches - if activated
113         //
115         // The name for the interpolated patch point field must be consistent
116         // with the interpolated volume point field.
117         // This could be done better.
118         const word pointFldName = "volPointInterpolate(" + tf.name() + ')';
120         for
121         (
122             int partId = partInfoPatches_.start();
123             partId < partInfoPatches_.end();
124             ++partId
125         )
126         {
127             const word patchName = getPartName(partId);
128             const label datasetNo = partDataset_[partId];
129             const label patchId = patches.findPatchID(patchName);
131             if (!partStatus_[partId] || datasetNo < 0 || patchId < 0)
132             {
133                 continue;
134             }
136             const fvPatchField<Type>& ptf = tf.boundaryField()[patchId];
138             if
139             (
140                 isType<emptyFvPatchField<Type> >(ptf)
141              ||
142                 (
143                     reader_->GetExtrapolatePatches()
144                 && !polyPatch::constraintType(patches[patchId].type())
145                 )
146             )
147             {
148                 fvPatch p(ptf.patch().patch(), tf.mesh().boundary());
150                 tmp<Field<Type> > tpptf
151                 (
152                     fvPatchField<Type>(p, tf).patchInternalField()
153                 );
155                 convertPatchField
156                 (
157                     tf.name(),
158                     tpptf(),
159                     output,
160                     partInfoPatches_,
161                     datasetNo
162                 );
164                 convertPatchPointField
165                 (
166                     pointFldName,
167                     ppInterpList[patchId].faceToPointInterpolate(tpptf)(),
168                     output,
169                     partInfoPatches_,
170                     datasetNo
171                 );
172             }
173             else
174             {
175                 convertPatchField
176                 (
177                     tf.name(),
178                     ptf,
179                     output,
180                     partInfoPatches_,
181                     datasetNo
182                 );
184                 convertPatchPointField
185                 (
186                     pointFldName,
187                     ppInterpList[patchId].faceToPointInterpolate(ptf)(),
188                     output,
189                     partInfoPatches_,
190                     datasetNo
191                 );
192             }
193         }
195         //
196         // Convert face zones - if activated
197         //
198         for
199         (
200             int partId = partInfoFaceZones_.start();
201             partId < partInfoFaceZones_.end();
202             ++partId
203         )
204         {
205             const word zoneName = getPartName(partId);
206             const label datasetNo = partDataset_[partId];
208             if (!partStatus_[partId] || datasetNo < 0)
209             {
210                 continue;
211             }
213             const faceZoneMesh& zMesh = mesh.faceZones();
214             const label zoneId = zMesh.findZoneID(zoneName);
216             if (zoneId < 0)
217             {
218                 continue;
219             }
221             convertFaceField
222             (
223                 tf,
224                 output,
225                 partInfoFaceZones_,
226                 datasetNo,
227                 mesh,
228                 zMesh[zoneId]
229             );
231             // TODO: points
232         }
234         //
235         // Convert face sets - if activated
236         //
237         for
238         (
239             int partId = partInfoFaceSets_.start();
240             partId < partInfoFaceSets_.end();
241             ++partId
242         )
243         {
244             const word selectName = getPartName(partId);
245             const label datasetNo = partDataset_[partId];
247             if (!partStatus_[partId] || datasetNo < 0)
248             {
249                 continue;
250             }
252             const faceSet fSet(mesh, selectName);
254             convertFaceField
255             (
256                 tf,
257                 output,
258                 partInfoFaceSets_,
259                 datasetNo,
260                 mesh,
261                 fSet
262             );
264             // TODO: points
265         }
266     }
270 template<class Type>
271 void Foam::vtkPV3Foam::convertVolFieldBlock
273     const GeometricField<Type, fvPatchField, volMesh>& tf,
274     autoPtr<GeometricField<Type, pointPatchField, pointMesh> >& ptfPtr,
275     vtkMultiBlockDataSet* output,
276     const partInfo& selector,
277     const List<polyDecomp>& decompLst
280     for (int partId = selector.start(); partId < selector.end(); ++partId)
281     {
282         const label datasetNo = partDataset_[partId];
284         if (datasetNo >= 0 && partStatus_[partId])
285         {
286             convertVolField
287             (
288                 tf,
289                 output,
290                 selector,
291                 datasetNo,
292                 decompLst[datasetNo]
293             );
295             if (!ptfPtr.valid())
296             {
297                 if (debug)
298                 {
299                     Info<< "convertVolFieldBlock interpolating:" << tf.name()
300                         << endl;
301                 }
303                 ptfPtr.reset
304                 (
305                     volPointInterpolation::New(tf.mesh()).interpolate(tf).ptr()
306                 );
307             }
309             convertPointField
310             (
311                 ptfPtr(),
312                 tf,
313                 output,
314                 selector,
315                 datasetNo,
316                 decompLst[datasetNo]
317             );
318         }
319     }
323 template<class Type>
324 void Foam::vtkPV3Foam::convertVolField
326     const GeometricField<Type, fvPatchField, volMesh>& tf,
327     vtkMultiBlockDataSet* output,
328     const partInfo& selector,
329     const label datasetNo,
330     const polyDecomp& decompInfo
333     const label nComp = pTraits<Type>::nComponents;
334     const labelList& superCells = decompInfo.superCells();
336     vtkFloatArray* celldata = vtkFloatArray::New();
337     celldata->SetNumberOfTuples(superCells.size());
338     celldata->SetNumberOfComponents(nComp);
339     celldata->Allocate(nComp*superCells.size());
340     celldata->SetName(tf.name().c_str());
342     if (debug)
343     {
344         Info<< "convert volField: "
345             << tf.name()
346             << " size = " << tf.size()
347             << " nComp=" << nComp
348             << " nTuples = " << superCells.size() <<  endl;
349     }
351     float vec[nComp];
352     forAll(superCells, i)
353     {
354         const Type& t = tf[superCells[i]];
355         for (direction d=0; d<nComp; d++)
356         {
357             vec[d] = component(t, d);
358         }
359         vtkOpenFOAMTupleRemap<Type>(vec);
361         celldata->InsertTuple(i, vec);
362     }
364     vtkUnstructuredGrid::SafeDownCast
365     (
366         GetDataSetFromBlock(output, selector, datasetNo)
367     )   ->GetCellData()
368         ->AddArray(celldata);
370     celldata->Delete();
374 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
376 #endif
378 // ************************************************************************* //