initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMesh.C
blob6550a960e14aac0b594aa80215e33f4deb536bd1
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 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 subsetter(mesh);
259         subsetter.setLargeCellSubset(zMesh[zoneId]);
261         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
262         (
263             subsetter.subMesh(),
264             zonePolyDecomp_[datasetNo]
265         );
267         if (vtkmesh)
268         {
269             // superCells + addPointCellLabels must contain global cell ids
270             inplaceRenumber
271             (
272                 subsetter.cellMap(),
273                 zonePolyDecomp_[datasetNo].superCells()
274             );
275             inplaceRenumber
276             (
277                 subsetter.cellMap(),
278                 zonePolyDecomp_[datasetNo].addPointCellLabels()
279             );
281             // copy pointMap as well, otherwise pointFields fail
282             zonePolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
284             AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
285             vtkmesh->Delete();
287             partDataset_[partId] = datasetNo++;
288         }
289     }
291     // anything added?
292     if (datasetNo)
293     {
294         ++blockNo;
295     }
297     if (debug)
298     {
299         Info<< "<end> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
300         printMemory();
301     }
305 void Foam::vtkPV3Foam::convertMeshCellSets
307     vtkMultiBlockDataSet* output,
308     int& blockNo
311     partInfo& selector = partInfoCellSets_;
312     selector.block(blockNo);   // set output block
313     label datasetNo = 0;       // restart at dataset 0
314     const fvMesh& mesh = *meshPtr_;
316     // resize for decomposed polyhedra
317     csetPolyDecomp_.setSize(selector.size());
319     if (debug)
320     {
321         Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
322         printMemory();
323     }
325     for (int partId = selector.start(); partId < selector.end(); ++partId)
326     {
327         const word partName = getPartName(partId);
329         if (!partStatus_[partId])
330         {
331             continue;
332         }
334         if (debug)
335         {
336             Info<< "Creating VTK mesh for cellSet=" << partName << endl;
337         }
339         const cellSet cSet(mesh, partName);
340         fvMeshSubset subsetter(mesh);
341         subsetter.setLargeCellSubset(cSet);
343         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
344         (
345             subsetter.subMesh(),
346             csetPolyDecomp_[datasetNo]
347         );
349         if (vtkmesh)
350         {
351             // superCells + addPointCellLabels must contain global cell ids
352             inplaceRenumber
353             (
354                 subsetter.cellMap(),
355                 csetPolyDecomp_[datasetNo].superCells()
356             );
357             inplaceRenumber
358             (
359                 subsetter.cellMap(),
360                 csetPolyDecomp_[datasetNo].addPointCellLabels()
361             );
363             // copy pointMap as well, otherwise pointFields fail
364             csetPolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
366             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
367             vtkmesh->Delete();
369             partDataset_[partId] = datasetNo++;
370         }
371     }
373     // anything added?
374     if (datasetNo)
375     {
376         ++blockNo;
377     }
379     if (debug)
380     {
381         Info<< "<end> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
382         printMemory();
383     }
387 void Foam::vtkPV3Foam::convertMeshFaceZones
389     vtkMultiBlockDataSet* output,
390     int& blockNo
393     partInfo& selector = partInfoFaceZones_;
394     selector.block(blockNo);   // set output block
395     label datasetNo = 0;       // restart at dataset 0
396     const fvMesh& mesh = *meshPtr_;
398     if (!selector.size())
399     {
400         return;
401     }
403     if (debug)
404     {
405         Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
406         printMemory();
407     }
409     const faceZoneMesh& zMesh = mesh.faceZones();
410     for (int partId = selector.start(); partId < selector.end(); ++partId)
411     {
412         const word zoneName = getPartName(partId);
413         const label  zoneId = zMesh.findZoneID(zoneName);
415         if (!partStatus_[partId] || zoneId < 0)
416         {
417             continue;
418         }
420         if (debug)
421         {
422             Info<< "Creating VTKmesh for faceZone[" << zoneId << "] "
423                 << zoneName << endl;
424         }
426         vtkPolyData* vtkmesh = faceZoneVTKMesh(mesh, zMesh[zoneId]);
427         if (vtkmesh)
428         {
429             AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
430             vtkmesh->Delete();
432             partDataset_[partId] = datasetNo++;
433         }
434     }
436     // anything added?
437     if (datasetNo)
438     {
439         ++blockNo;
440     }
442     if (debug)
443     {
444         Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
445         printMemory();
446     }
450 void Foam::vtkPV3Foam::convertMeshFaceSets
452     vtkMultiBlockDataSet* output,
453     int& blockNo
456     partInfo& selector = partInfoFaceSets_;
457     selector.block(blockNo);   // set output block
458     label datasetNo = 0;       // restart at dataset 0
459     const fvMesh& mesh = *meshPtr_;
461     if (debug)
462     {
463         Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
464         printMemory();
465     }
467     for (int partId = selector.start(); partId < selector.end(); ++partId)
468     {
469         const word partName = getPartName(partId);
471         if (!partStatus_[partId])
472         {
473             continue;
474         }
476         if (debug)
477         {
478             Info<< "Creating VTK mesh for faceSet=" << partName << endl;
479         }
481         const faceSet fSet(mesh, partName);
483         vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSet);
484         if (vtkmesh)
485         {
486             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
487             vtkmesh->Delete();
489             partDataset_[partId] = datasetNo++;
490         }
491     }
493     // anything added?
494     if (datasetNo)
495     {
496         ++blockNo;
497     }
499     if (debug)
500     {
501         Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
502         printMemory();
503     }
507 void Foam::vtkPV3Foam::convertMeshPointZones
509     vtkMultiBlockDataSet* output,
510     int& blockNo
513     partInfo& selector = partInfoPointZones_;
514     selector.block(blockNo);   // set output block
515     label datasetNo = 0;       // restart at dataset 0
516     const fvMesh& mesh = *meshPtr_;
518     if (debug)
519     {
520         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
521         printMemory();
522     }
524     if (selector.size())
525     {
526         const pointZoneMesh& zMesh = mesh.pointZones();
527         for (int partId = selector.start(); partId < selector.end(); ++partId)
528         {
529             word zoneName = getPartName(partId);
530             label zoneId = zMesh.findZoneID(zoneName);
532             if (!partStatus_[partId] || zoneId < 0)
533             {
534                 continue;
535             }
537             vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
538             if (vtkmesh)
539             {
540                 AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
541                 vtkmesh->Delete();
543                 partDataset_[partId] = datasetNo++;
544             }
545         }
546     }
548     // anything added?
549     if (datasetNo)
550     {
551         ++blockNo;
552     }
554     if (debug)
555     {
556         Info<< "<end> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
557         printMemory();
558     }
563 void Foam::vtkPV3Foam::convertMeshPointSets
565     vtkMultiBlockDataSet* output,
566     int& blockNo
569     partInfo& selector = partInfoPointSets_;
570     selector.block(blockNo);   // set output block
571     label datasetNo = 0;       // restart at dataset 0
572     const fvMesh& mesh = *meshPtr_;
574     if (debug)
575     {
576         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
577         printMemory();
578     }
580     for (int partId = selector.start(); partId < selector.end(); ++partId)
581     {
582         word partName = getPartName(partId);
584         if (!partStatus_[partId])
585         {
586             continue;
587         }
589         if (debug)
590         {
591             Info<< "Creating VTK mesh for pointSet=" << partName << endl;
592         }
594         const pointSet pSet(mesh, partName);
596         vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSet);
597         if (vtkmesh)
598         {
599             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
600             vtkmesh->Delete();
602             partDataset_[partId] = datasetNo++;
603         }
604     }
606     // anything added?
607     if (datasetNo)
608     {
609         ++blockNo;
610     }
612     if (debug)
613     {
614         Info<< "<end> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
615         printMemory();
616     }
619 // ************************************************************************* //