BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMesh.C
blob867ddf099d9a35b4dc59ff572e74cc99cca9af4e
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 Description
27 \*---------------------------------------------------------------------------*/
29 #include "vtkPV3Foam.H"
31 // Foam includes
32 #include "cellSet.H"
33 #include "faceSet.H"
34 #include "pointSet.H"
35 #include "fvMeshSubset.H"
36 #include "vtkPV3FoamReader.h"
38 // VTK includes
39 #include "vtkDataArraySelection.h"
40 #include "vtkMultiBlockDataSet.h"
41 #include "vtkPolyData.h"
42 #include "vtkUnstructuredGrid.h"
44 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
46 void Foam::vtkPV3Foam::convertMeshVolume
48     vtkMultiBlockDataSet* output,
49     int& blockNo
52     partInfo& selector = partInfoVolume_;
53     selector.block(blockNo);   // set output block
54     label datasetNo = 0;       // restart at dataset 0
55     const fvMesh& mesh = *meshPtr_;
57     // resize for decomposed polyhedra
58     regionPolyDecomp_.setSize(selector.size());
60     if (debug)
61     {
62         Info<< "<beg> Foam::vtkPV3Foam::convertMeshVolume" << endl;
63         printMemory();
64     }
66     // Convert the internalMesh
67     // this looks like more than one part, but it isn't
68     for (int partId = selector.start(); partId < selector.end(); ++partId)
69     {
70         const word partName = "internalMesh";
72         if (!partStatus_[partId])
73         {
74             continue;
75         }
77         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
78         (
79             mesh,
80             regionPolyDecomp_[datasetNo]
81         );
83         if (vtkmesh)
84         {
85             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
86             vtkmesh->Delete();
88             partDataset_[partId] = datasetNo++;
89         }
90     }
92     // anything added?
93     if (datasetNo)
94     {
95         ++blockNo;
96     }
98     if (debug)
99     {
100         Info<< "<end> Foam::vtkPV3Foam::convertMeshVolume" << endl;
101         printMemory();
102     }
106 void Foam::vtkPV3Foam::convertMeshLagrangian
108     vtkMultiBlockDataSet* output,
109     int& blockNo
112     partInfo& selector = partInfoLagrangian_;
113     selector.block(blockNo);   // set output block
114     label datasetNo = 0;       // restart at dataset 0
115     const fvMesh& mesh = *meshPtr_;
117     if (debug)
118     {
119         Info<< "<beg> Foam::vtkPV3Foam::convertMeshLagrangian" << endl;
120         printMemory();
121     }
123     for (int partId = selector.start(); partId < selector.end(); ++partId)
124     {
125         const word cloudName = getPartName(partId);
127         if (!partStatus_[partId])
128         {
129             continue;
130         }
132         vtkPolyData* vtkmesh = lagrangianVTKMesh(mesh, cloudName);
134         if (vtkmesh)
135         {
136             AddToBlock(output, vtkmesh, selector, datasetNo, cloudName);
137             vtkmesh->Delete();
139             partDataset_[partId] = datasetNo++;
140         }
141     }
143     // anything added?
144     if (datasetNo)
145     {
146         ++blockNo;
147     }
149     if (debug)
150     {
151         Info<< "<end> Foam::vtkPV3Foam::convertMeshLagrangian" << endl;
152         printMemory();
153     }
157 void Foam::vtkPV3Foam::convertMeshPatches
159     vtkMultiBlockDataSet* output,
160     int& blockNo
163     partInfo& selector = partInfoPatches_;
164     selector.block(blockNo);   // set output block
165     label datasetNo = 0;       // restart at dataset 0
166     const fvMesh& mesh = *meshPtr_;
167     const polyBoundaryMesh& patches = mesh.boundaryMesh();
169     if (debug)
170     {
171         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPatches" << endl;
172         printMemory();
173     }
175     for (int partId = selector.start(); partId < selector.end(); ++partId)
176     {
177         const word patchName = getPartName(partId);
178         const label  patchId = patches.findPatchID(patchName);
180         if (!partStatus_[partId] || patchId < 0)
181         {
182             continue;
183         }
185         if (debug)
186         {
187             Info<< "Creating VTK mesh for patch[" << patchId <<"] "
188                 << patchName  << endl;
189         }
191         vtkPolyData* vtkmesh = patchVTKMesh(patches[patchId]);
193         if (vtkmesh)
194         {
195             AddToBlock(output, vtkmesh, selector, datasetNo, patchName);
196             vtkmesh->Delete();
198             partDataset_[partId] = datasetNo++;
199         }
200     }
202     // anything added?
203     if (datasetNo)
204     {
205         ++blockNo;
206     }
208     if (debug)
209     {
210         Info<< "<end> Foam::vtkPV3Foam::convertMeshPatches" << endl;
211         printMemory();
212     }
216 void Foam::vtkPV3Foam::convertMeshCellZones
218     vtkMultiBlockDataSet* output,
219     int& blockNo
222     partInfo& selector = partInfoCellZones_;
223     selector.block(blockNo);   // set output block
224     label datasetNo = 0;       // restart at dataset 0
225     const fvMesh& mesh = *meshPtr_;
227     // resize for decomposed polyhedra
228     zonePolyDecomp_.setSize(selector.size());
230     if (!selector.size())
231     {
232         return;
233     }
235     if (debug)
236     {
237         Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
238         printMemory();
239     }
241     const cellZoneMesh& zMesh = mesh.cellZones();
242     for (int partId = selector.start(); partId < selector.end(); ++partId)
243     {
244         const word zoneName = getPartName(partId);
245         const label  zoneId = zMesh.findZoneID(zoneName);
247         if (!partStatus_[partId] || zoneId < 0)
248         {
249             continue;
250         }
252         if (debug)
253         {
254             Info<< "Creating VTK mesh for cellZone[" << zoneId << "] "
255                 << zoneName << endl;
256         }
258         fvMeshSubset subsetMesh
259         (
260             IOobject
261             (
262                 "set",
263                 mesh.time().constant(),
264                 mesh,
265                 IOobject::NO_READ,
266                 IOobject::NO_WRITE
267             ),
268             mesh
269         );
271         subsetMesh.setLargeCellSubset(labelHashSet(zMesh[zoneId]));
273         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
274         (
275             subsetMesh.subMesh(),
276             zonePolyDecomp_[datasetNo]
277         );
279         if (vtkmesh)
280         {
281             // superCells + addPointCellLabels must contain global cell ids
282             inplaceRenumber
283             (
284                 subsetMesh.cellMap(),
285                 zonePolyDecomp_[datasetNo].superCells()
286             );
287             inplaceRenumber
288             (
289                 subsetMesh.cellMap(),
290                 zonePolyDecomp_[datasetNo].addPointCellLabels()
291             );
293             // copy pointMap as well, otherwise pointFields fail
294             zonePolyDecomp_[datasetNo].pointMap() = subsetMesh.pointMap();
296             AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
297             vtkmesh->Delete();
299             partDataset_[partId] = datasetNo++;
300         }
301     }
303     // anything added?
304     if (datasetNo)
305     {
306         ++blockNo;
307     }
309     if (debug)
310     {
311         Info<< "<end> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
312         printMemory();
313     }
317 void Foam::vtkPV3Foam::convertMeshCellSets
319     vtkMultiBlockDataSet* output,
320     int& blockNo
323     partInfo& selector = partInfoCellSets_;
324     selector.block(blockNo);   // set output block
325     label datasetNo = 0;       // restart at dataset 0
326     const fvMesh& mesh = *meshPtr_;
328     // resize for decomposed polyhedra
329     csetPolyDecomp_.setSize(selector.size());
331     if (debug)
332     {
333         Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
334         printMemory();
335     }
337     for (int partId = selector.start(); partId < selector.end(); ++partId)
338     {
339         const word partName = getPartName(partId);
341         if (!partStatus_[partId])
342         {
343             continue;
344         }
346         if (debug)
347         {
348             Info<< "Creating VTK mesh for cellSet=" << partName << endl;
349         }
351         const cellSet cSet(mesh, partName);
353         fvMeshSubset subsetMesh
354         (
355             IOobject
356             (
357                 "set",
358                 mesh.time().constant(),
359                 mesh,
360                 IOobject::NO_READ,
361                 IOobject::NO_WRITE
362             ),
363             mesh
364         );
366         subsetMesh.setLargeCellSubset(cSet);
368         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
369         (
370             subsetMesh.subMesh(),
371             csetPolyDecomp_[datasetNo]
372         );
374         if (vtkmesh)
375         {
376             // superCells + addPointCellLabels must contain global cell ids
377             inplaceRenumber
378             (
379                 subsetMesh.cellMap(),
380                 csetPolyDecomp_[datasetNo].superCells()
381             );
382             inplaceRenumber
383             (
384                 subsetMesh.cellMap(),
385                 csetPolyDecomp_[datasetNo].addPointCellLabels()
386             );
388             // copy pointMap as well, otherwise pointFields fail
389             csetPolyDecomp_[datasetNo].pointMap() = subsetMesh.pointMap();
391             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
392             vtkmesh->Delete();
394             partDataset_[partId] = datasetNo++;
395         }
396     }
398     // anything added?
399     if (datasetNo)
400     {
401         ++blockNo;
402     }
404     if (debug)
405     {
406         Info<< "<end> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
407         printMemory();
408     }
412 void Foam::vtkPV3Foam::convertMeshFaceZones
414     vtkMultiBlockDataSet* output,
415     int& blockNo
418     partInfo& selector = partInfoFaceZones_;
419     selector.block(blockNo);   // set output block
420     label datasetNo = 0;       // restart at dataset 0
421     const fvMesh& mesh = *meshPtr_;
423     if (!selector.size())
424     {
425         return;
426     }
428     if (debug)
429     {
430         Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
431         printMemory();
432     }
434     const faceZoneMesh& zMesh = mesh.faceZones();
435     for (int partId = selector.start(); partId < selector.end(); ++partId)
436     {
437         const word zoneName = getPartName(partId);
438         const label  zoneId = zMesh.findZoneID(zoneName);
440         if (!partStatus_[partId] || zoneId < 0)
441         {
442             continue;
443         }
445         if (debug)
446         {
447             Info<< "Creating VTKmesh for faceZone[" << zoneId << "] "
448                 << zoneName << endl;
449         }
451         vtkPolyData* vtkmesh = faceZoneVTKMesh(mesh, zMesh[zoneId]);
452         if (vtkmesh)
453         {
454             AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
455             vtkmesh->Delete();
457             partDataset_[partId] = datasetNo++;
458         }
459     }
461     // anything added?
462     if (datasetNo)
463     {
464         ++blockNo;
465     }
467     if (debug)
468     {
469         Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
470         printMemory();
471     }
475 void Foam::vtkPV3Foam::convertMeshFaceSets
477     vtkMultiBlockDataSet* output,
478     int& blockNo
481     partInfo& selector = partInfoFaceSets_;
482     selector.block(blockNo);   // set output block
483     label datasetNo = 0;       // restart at dataset 0
484     const fvMesh& mesh = *meshPtr_;
486     if (debug)
487     {
488         Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
489         printMemory();
490     }
492     for (int partId = selector.start(); partId < selector.end(); ++partId)
493     {
494         const word partName = getPartName(partId);
496         if (!partStatus_[partId])
497         {
498             continue;
499         }
501         if (debug)
502         {
503             Info<< "Creating VTK mesh for faceSet=" << partName << endl;
504         }
506         const faceSet fSet(mesh, partName);
508         vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSet);
509         if (vtkmesh)
510         {
511             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
512             vtkmesh->Delete();
514             partDataset_[partId] = datasetNo++;
515         }
516     }
518     // anything added?
519     if (datasetNo)
520     {
521         ++blockNo;
522     }
524     if (debug)
525     {
526         Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
527         printMemory();
528     }
532 void Foam::vtkPV3Foam::convertMeshPointZones
534     vtkMultiBlockDataSet* output,
535     int& blockNo
538     partInfo& selector = partInfoPointZones_;
539     selector.block(blockNo);   // set output block
540     label datasetNo = 0;       // restart at dataset 0
541     const fvMesh& mesh = *meshPtr_;
543     if (debug)
544     {
545         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
546         printMemory();
547     }
549     if (selector.size())
550     {
551         const pointZoneMesh& zMesh = mesh.pointZones();
552         for (int partId = selector.start(); partId < selector.end(); ++partId)
553         {
554             word zoneName = getPartName(partId);
555             label zoneId = zMesh.findZoneID(zoneName);
557             if (!partStatus_[partId] || zoneId < 0)
558             {
559                 continue;
560             }
562             vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
563             if (vtkmesh)
564             {
565                 AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
566                 vtkmesh->Delete();
568                 partDataset_[partId] = datasetNo++;
569             }
570         }
571     }
573     // anything added?
574     if (datasetNo)
575     {
576         ++blockNo;
577     }
579     if (debug)
580     {
581         Info<< "<end> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
582         printMemory();
583     }
588 void Foam::vtkPV3Foam::convertMeshPointSets
590     vtkMultiBlockDataSet* output,
591     int& blockNo
594     partInfo& selector = partInfoPointSets_;
595     selector.block(blockNo);   // set output block
596     label datasetNo = 0;       // restart at dataset 0
597     const fvMesh& mesh = *meshPtr_;
599     if (debug)
600     {
601         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
602         printMemory();
603     }
605     for (int partId = selector.start(); partId < selector.end(); ++partId)
606     {
607         word partName = getPartName(partId);
609         if (!partStatus_[partId])
610         {
611             continue;
612         }
614         if (debug)
615         {
616             Info<< "Creating VTK mesh for pointSet=" << partName << endl;
617         }
619         const pointSet pSet(mesh, partName);
621         vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSet);
622         if (vtkmesh)
623         {
624             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
625             vtkmesh->Delete();
627             partDataset_[partId] = datasetNo++;
628         }
629     }
631     // anything added?
632     if (datasetNo)
633     {
634         ++blockNo;
635     }
637     if (debug)
638     {
639         Info<< "<end> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
640         printMemory();
641     }
644 // ************************************************************************* //