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 / vtkPV3FoamMesh.C
blob348f436ce0e7a874155047abb7c799e262fff7e5
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 Description
26 \*---------------------------------------------------------------------------*/
28 #include "vtkPV3Foam.H"
30 // Foam includes
31 #include "cellSet.H"
32 #include "faceSet.H"
33 #include "pointSet.H"
34 #include "fvMeshSubset.H"
35 #include "vtkPV3FoamReader.h"
37 // VTK includes
38 #include "vtkDataArraySelection.h"
39 #include "vtkMultiBlockDataSet.h"
40 #include "vtkPolyData.h"
41 #include "vtkUnstructuredGrid.h"
43 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
45 void Foam::vtkPV3Foam::convertMeshVolume
47     vtkMultiBlockDataSet* output,
48     int& blockNo
51     partInfo& selector = partInfoVolume_;
52     selector.block(blockNo);   // set output block
53     label datasetNo = 0;       // restart at dataset 0
54     const fvMesh& mesh = *meshPtr_;
56     // resize for decomposed polyhedra
57     regionPolyDecomp_.setSize(selector.size());
59     if (debug)
60     {
61         Info<< "<beg> Foam::vtkPV3Foam::convertMeshVolume" << endl;
62         printMemory();
63     }
65     // Convert the internalMesh
66     // this looks like more than one part, but it isn't
67     for (int partId = selector.start(); partId < selector.end(); ++partId)
68     {
69         const word partName = "internalMesh";
71         if (!partStatus_[partId])
72         {
73             continue;
74         }
76         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
77         (
78             mesh,
79             regionPolyDecomp_[datasetNo]
80         );
82         if (vtkmesh)
83         {
84             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
85             vtkmesh->Delete();
87             partDataset_[partId] = datasetNo++;
88         }
89     }
91     // anything added?
92     if (datasetNo)
93     {
94         ++blockNo;
95     }
97     if (debug)
98     {
99         Info<< "<end> Foam::vtkPV3Foam::convertMeshVolume" << endl;
100         printMemory();
101     }
105 void Foam::vtkPV3Foam::convertMeshLagrangian
107     vtkMultiBlockDataSet* output,
108     int& blockNo
111     partInfo& selector = partInfoLagrangian_;
112     selector.block(blockNo);   // set output block
113     label datasetNo = 0;       // restart at dataset 0
114     const fvMesh& mesh = *meshPtr_;
116     if (debug)
117     {
118         Info<< "<beg> Foam::vtkPV3Foam::convertMeshLagrangian" << endl;
119         printMemory();
120     }
122     for (int partId = selector.start(); partId < selector.end(); ++partId)
123     {
124         const word cloudName = getPartName(partId);
126         if (!partStatus_[partId])
127         {
128             continue;
129         }
131         vtkPolyData* vtkmesh = lagrangianVTKMesh(mesh, cloudName);
133         if (vtkmesh)
134         {
135             AddToBlock(output, vtkmesh, selector, datasetNo, cloudName);
136             vtkmesh->Delete();
138             partDataset_[partId] = datasetNo++;
139         }
140     }
142     // anything added?
143     if (datasetNo)
144     {
145         ++blockNo;
146     }
148     if (debug)
149     {
150         Info<< "<end> Foam::vtkPV3Foam::convertMeshLagrangian" << endl;
151         printMemory();
152     }
156 void Foam::vtkPV3Foam::convertMeshPatches
158     vtkMultiBlockDataSet* output,
159     int& blockNo
162     partInfo& selector = partInfoPatches_;
163     selector.block(blockNo);   // set output block
164     label datasetNo = 0;       // restart at dataset 0
165     const fvMesh& mesh = *meshPtr_;
166     const polyBoundaryMesh& patches = mesh.boundaryMesh();
168     if (debug)
169     {
170         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPatches" << endl;
171         printMemory();
172     }
174     for (int partId = selector.start(); partId < selector.end(); ++partId)
175     {
176         const word patchName = getPartName(partId);
177         const label  patchId = patches.findPatchID(patchName);
179         if (!partStatus_[partId] || patchId < 0)
180         {
181             continue;
182         }
184         if (debug)
185         {
186             Info<< "Creating VTK mesh for patch[" << patchId <<"] "
187                 << patchName  << endl;
188         }
190         vtkPolyData* vtkmesh = patchVTKMesh(patches[patchId]);
192         if (vtkmesh)
193         {
194             AddToBlock(output, vtkmesh, selector, datasetNo, patchName);
195             vtkmesh->Delete();
197             partDataset_[partId] = datasetNo++;
198         }
199     }
201     // anything added?
202     if (datasetNo)
203     {
204         ++blockNo;
205     }
207     if (debug)
208     {
209         Info<< "<end> Foam::vtkPV3Foam::convertMeshPatches" << endl;
210         printMemory();
211     }
215 void Foam::vtkPV3Foam::convertMeshCellZones
217     vtkMultiBlockDataSet* output,
218     int& blockNo
221     partInfo& selector = partInfoCellZones_;
222     selector.block(blockNo);   // set output block
223     label datasetNo = 0;       // restart at dataset 0
224     const fvMesh& mesh = *meshPtr_;
226     // resize for decomposed polyhedra
227     zonePolyDecomp_.setSize(selector.size());
229     if (!selector.size())
230     {
231         return;
232     }
234     if (debug)
235     {
236         Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
237         printMemory();
238     }
240     const cellZoneMesh& zMesh = mesh.cellZones();
241     for (int partId = selector.start(); partId < selector.end(); ++partId)
242     {
243         const word zoneName = getPartName(partId);
244         const label  zoneId = zMesh.findZoneID(zoneName);
246         if (!partStatus_[partId] || zoneId < 0)
247         {
248             continue;
249         }
251         if (debug)
252         {
253             Info<< "Creating VTK mesh for cellZone[" << zoneId << "] "
254                 << zoneName << endl;
255         }
257         fvMeshSubset subsetMesh
258         (
259             IOobject
260             (
261                 "set",
262                 mesh.time().constant(),
263                 mesh,
264                 IOobject::NO_READ,
265                 IOobject::NO_WRITE
266             ),
267             mesh
268         );
270         subsetMesh.setLargeCellSubset(labelHashSet(zMesh[zoneId]));
272         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
273         (
274             subsetMesh.subMesh(),
275             zonePolyDecomp_[datasetNo]
276         );
278         if (vtkmesh)
279         {
280             // superCells + addPointCellLabels must contain global cell ids
281             inplaceRenumber
282             (
283                 subsetMesh.cellMap(),
284                 zonePolyDecomp_[datasetNo].superCells()
285             );
286             inplaceRenumber
287             (
288                 subsetMesh.cellMap(),
289                 zonePolyDecomp_[datasetNo].addPointCellLabels()
290             );
292             // copy pointMap as well, otherwise pointFields fail
293             zonePolyDecomp_[datasetNo].pointMap() = subsetMesh.pointMap();
295             AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
296             vtkmesh->Delete();
298             partDataset_[partId] = datasetNo++;
299         }
300     }
302     // anything added?
303     if (datasetNo)
304     {
305         ++blockNo;
306     }
308     if (debug)
309     {
310         Info<< "<end> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
311         printMemory();
312     }
316 void Foam::vtkPV3Foam::convertMeshCellSets
318     vtkMultiBlockDataSet* output,
319     int& blockNo
322     partInfo& selector = partInfoCellSets_;
323     selector.block(blockNo);   // set output block
324     label datasetNo = 0;       // restart at dataset 0
325     const fvMesh& mesh = *meshPtr_;
327     // resize for decomposed polyhedra
328     csetPolyDecomp_.setSize(selector.size());
330     if (debug)
331     {
332         Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
333         printMemory();
334     }
336     for (int partId = selector.start(); partId < selector.end(); ++partId)
337     {
338         const word partName = getPartName(partId);
340         if (!partStatus_[partId])
341         {
342             continue;
343         }
345         if (debug)
346         {
347             Info<< "Creating VTK mesh for cellSet=" << partName << endl;
348         }
350         const cellSet cSet(mesh, partName);
352         fvMeshSubset subsetMesh
353         (
354             IOobject
355             (
356                 "set",
357                 mesh.time().constant(),
358                 mesh,
359                 IOobject::NO_READ,
360                 IOobject::NO_WRITE
361             ),
362             mesh
363         );
365         subsetMesh.setLargeCellSubset(cSet);
367         vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
368         (
369             subsetMesh.subMesh(),
370             csetPolyDecomp_[datasetNo]
371         );
373         if (vtkmesh)
374         {
375             // superCells + addPointCellLabels must contain global cell ids
376             inplaceRenumber
377             (
378                 subsetMesh.cellMap(),
379                 csetPolyDecomp_[datasetNo].superCells()
380             );
381             inplaceRenumber
382             (
383                 subsetMesh.cellMap(),
384                 csetPolyDecomp_[datasetNo].addPointCellLabels()
385             );
387             // copy pointMap as well, otherwise pointFields fail
388             csetPolyDecomp_[datasetNo].pointMap() = subsetMesh.pointMap();
390             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
391             vtkmesh->Delete();
393             partDataset_[partId] = datasetNo++;
394         }
395     }
397     // anything added?
398     if (datasetNo)
399     {
400         ++blockNo;
401     }
403     if (debug)
404     {
405         Info<< "<end> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
406         printMemory();
407     }
411 void Foam::vtkPV3Foam::convertMeshFaceZones
413     vtkMultiBlockDataSet* output,
414     int& blockNo
417     partInfo& selector = partInfoFaceZones_;
418     selector.block(blockNo);   // set output block
419     label datasetNo = 0;       // restart at dataset 0
420     const fvMesh& mesh = *meshPtr_;
422     if (!selector.size())
423     {
424         return;
425     }
427     if (debug)
428     {
429         Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
430         printMemory();
431     }
433     const faceZoneMesh& zMesh = mesh.faceZones();
434     for (int partId = selector.start(); partId < selector.end(); ++partId)
435     {
436         const word zoneName = getPartName(partId);
437         const label  zoneId = zMesh.findZoneID(zoneName);
439         if (!partStatus_[partId] || zoneId < 0)
440         {
441             continue;
442         }
444         if (debug)
445         {
446             Info<< "Creating VTKmesh for faceZone[" << zoneId << "] "
447                 << zoneName << endl;
448         }
450         vtkPolyData* vtkmesh = faceZoneVTKMesh(mesh, zMesh[zoneId]);
451         if (vtkmesh)
452         {
453             AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
454             vtkmesh->Delete();
456             partDataset_[partId] = datasetNo++;
457         }
458     }
460     // anything added?
461     if (datasetNo)
462     {
463         ++blockNo;
464     }
466     if (debug)
467     {
468         Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
469         printMemory();
470     }
474 void Foam::vtkPV3Foam::convertMeshFaceSets
476     vtkMultiBlockDataSet* output,
477     int& blockNo
480     partInfo& selector = partInfoFaceSets_;
481     selector.block(blockNo);   // set output block
482     label datasetNo = 0;       // restart at dataset 0
483     const fvMesh& mesh = *meshPtr_;
485     if (debug)
486     {
487         Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
488         printMemory();
489     }
491     for (int partId = selector.start(); partId < selector.end(); ++partId)
492     {
493         const word partName = getPartName(partId);
495         if (!partStatus_[partId])
496         {
497             continue;
498         }
500         if (debug)
501         {
502             Info<< "Creating VTK mesh for faceSet=" << partName << endl;
503         }
505         const faceSet fSet(mesh, partName);
507         vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSet);
508         if (vtkmesh)
509         {
510             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
511             vtkmesh->Delete();
513             partDataset_[partId] = datasetNo++;
514         }
515     }
517     // anything added?
518     if (datasetNo)
519     {
520         ++blockNo;
521     }
523     if (debug)
524     {
525         Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
526         printMemory();
527     }
531 void Foam::vtkPV3Foam::convertMeshPointZones
533     vtkMultiBlockDataSet* output,
534     int& blockNo
537     partInfo& selector = partInfoPointZones_;
538     selector.block(blockNo);   // set output block
539     label datasetNo = 0;       // restart at dataset 0
540     const fvMesh& mesh = *meshPtr_;
542     if (debug)
543     {
544         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
545         printMemory();
546     }
548     if (selector.size())
549     {
550         const pointZoneMesh& zMesh = mesh.pointZones();
551         for (int partId = selector.start(); partId < selector.end(); ++partId)
552         {
553             word zoneName = getPartName(partId);
554             label zoneId = zMesh.findZoneID(zoneName);
556             if (!partStatus_[partId] || zoneId < 0)
557             {
558                 continue;
559             }
561             vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
562             if (vtkmesh)
563             {
564                 AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
565                 vtkmesh->Delete();
567                 partDataset_[partId] = datasetNo++;
568             }
569         }
570     }
572     // anything added?
573     if (datasetNo)
574     {
575         ++blockNo;
576     }
578     if (debug)
579     {
580         Info<< "<end> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
581         printMemory();
582     }
587 void Foam::vtkPV3Foam::convertMeshPointSets
589     vtkMultiBlockDataSet* output,
590     int& blockNo
593     partInfo& selector = partInfoPointSets_;
594     selector.block(blockNo);   // set output block
595     label datasetNo = 0;       // restart at dataset 0
596     const fvMesh& mesh = *meshPtr_;
598     if (debug)
599     {
600         Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
601         printMemory();
602     }
604     for (int partId = selector.start(); partId < selector.end(); ++partId)
605     {
606         word partName = getPartName(partId);
608         if (!partStatus_[partId])
609         {
610             continue;
611         }
613         if (debug)
614         {
615             Info<< "Creating VTK mesh for pointSet=" << partName << endl;
616         }
618         const pointSet pSet(mesh, partName);
620         vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSet);
621         if (vtkmesh)
622         {
623             AddToBlock(output, vtkmesh, selector, datasetNo, partName);
624             vtkmesh->Delete();
626             partDataset_[partId] = datasetNo++;
627         }
628     }
630     // anything added?
631     if (datasetNo)
632     {
633         ++blockNo;
634     }
636     if (debug)
637     {
638         Info<< "<end> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
639         printMemory();
640     }
643 // ************************************************************************* //