1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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/>.
26 \*---------------------------------------------------------------------------*/
28 #include "vtkPV3Foam.H"
34 #include "fvMeshSubset.H"
35 #include "vtkPV3FoamReader.h"
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,
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());
61 Info<< "<beg> Foam::vtkPV3Foam::convertMeshVolume" << endl;
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)
69 const word partName = "internalMesh";
71 if (!partStatus_[partId])
76 vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
79 regionPolyDecomp_[datasetNo]
84 AddToBlock(output, vtkmesh, selector, datasetNo, partName);
87 partDataset_[partId] = datasetNo++;
99 Info<< "<end> Foam::vtkPV3Foam::convertMeshVolume" << endl;
105 void Foam::vtkPV3Foam::convertMeshLagrangian
107 vtkMultiBlockDataSet* output,
111 partInfo& selector = partInfoLagrangian_;
112 selector.block(blockNo); // set output block
113 label datasetNo = 0; // restart at dataset 0
114 const fvMesh& mesh = *meshPtr_;
118 Info<< "<beg> Foam::vtkPV3Foam::convertMeshLagrangian" << endl;
122 for (int partId = selector.start(); partId < selector.end(); ++partId)
124 const word cloudName = getPartName(partId);
126 if (!partStatus_[partId])
131 vtkPolyData* vtkmesh = lagrangianVTKMesh(mesh, cloudName);
135 AddToBlock(output, vtkmesh, selector, datasetNo, cloudName);
138 partDataset_[partId] = datasetNo++;
150 Info<< "<end> Foam::vtkPV3Foam::convertMeshLagrangian" << endl;
156 void Foam::vtkPV3Foam::convertMeshPatches
158 vtkMultiBlockDataSet* output,
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();
170 Info<< "<beg> Foam::vtkPV3Foam::convertMeshPatches" << endl;
174 for (int partId = selector.start(); partId < selector.end(); ++partId)
176 const word patchName = getPartName(partId);
177 const label patchId = patches.findPatchID(patchName);
179 if (!partStatus_[partId] || patchId < 0)
186 Info<< "Creating VTK mesh for patch[" << patchId <<"] "
187 << patchName << endl;
190 vtkPolyData* vtkmesh = patchVTKMesh(patches[patchId]);
194 AddToBlock(output, vtkmesh, selector, datasetNo, patchName);
197 partDataset_[partId] = datasetNo++;
209 Info<< "<end> Foam::vtkPV3Foam::convertMeshPatches" << endl;
215 void Foam::vtkPV3Foam::convertMeshCellZones
217 vtkMultiBlockDataSet* output,
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())
236 Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
240 const cellZoneMesh& zMesh = mesh.cellZones();
241 for (int partId = selector.start(); partId < selector.end(); ++partId)
243 const word zoneName = getPartName(partId);
244 const label zoneId = zMesh.findZoneID(zoneName);
246 if (!partStatus_[partId] || zoneId < 0)
253 Info<< "Creating VTK mesh for cellZone[" << zoneId << "] "
257 fvMeshSubset subsetMesh
262 mesh.time().constant(),
270 subsetMesh.setLargeCellSubset(labelHashSet(zMesh[zoneId]));
272 vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
274 subsetMesh.subMesh(),
275 zonePolyDecomp_[datasetNo]
280 // superCells + addPointCellLabels must contain global cell ids
283 subsetMesh.cellMap(),
284 zonePolyDecomp_[datasetNo].superCells()
288 subsetMesh.cellMap(),
289 zonePolyDecomp_[datasetNo].addPointCellLabels()
292 // copy pointMap as well, otherwise pointFields fail
293 zonePolyDecomp_[datasetNo].pointMap() = subsetMesh.pointMap();
295 AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
298 partDataset_[partId] = datasetNo++;
310 Info<< "<end> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
316 void Foam::vtkPV3Foam::convertMeshCellSets
318 vtkMultiBlockDataSet* output,
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());
332 Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
336 for (int partId = selector.start(); partId < selector.end(); ++partId)
338 const word partName = getPartName(partId);
340 if (!partStatus_[partId])
347 Info<< "Creating VTK mesh for cellSet=" << partName << endl;
350 const cellSet cSet(mesh, partName);
352 fvMeshSubset subsetMesh
357 mesh.time().constant(),
365 subsetMesh.setLargeCellSubset(cSet);
367 vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
369 subsetMesh.subMesh(),
370 csetPolyDecomp_[datasetNo]
375 // superCells + addPointCellLabels must contain global cell ids
378 subsetMesh.cellMap(),
379 csetPolyDecomp_[datasetNo].superCells()
383 subsetMesh.cellMap(),
384 csetPolyDecomp_[datasetNo].addPointCellLabels()
387 // copy pointMap as well, otherwise pointFields fail
388 csetPolyDecomp_[datasetNo].pointMap() = subsetMesh.pointMap();
390 AddToBlock(output, vtkmesh, selector, datasetNo, partName);
393 partDataset_[partId] = datasetNo++;
405 Info<< "<end> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
411 void Foam::vtkPV3Foam::convertMeshFaceZones
413 vtkMultiBlockDataSet* output,
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())
429 Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
433 const faceZoneMesh& zMesh = mesh.faceZones();
434 for (int partId = selector.start(); partId < selector.end(); ++partId)
436 const word zoneName = getPartName(partId);
437 const label zoneId = zMesh.findZoneID(zoneName);
439 if (!partStatus_[partId] || zoneId < 0)
446 Info<< "Creating VTKmesh for faceZone[" << zoneId << "] "
450 vtkPolyData* vtkmesh = faceZoneVTKMesh(mesh, zMesh[zoneId]);
453 AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
456 partDataset_[partId] = datasetNo++;
468 Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
474 void Foam::vtkPV3Foam::convertMeshFaceSets
476 vtkMultiBlockDataSet* output,
480 partInfo& selector = partInfoFaceSets_;
481 selector.block(blockNo); // set output block
482 label datasetNo = 0; // restart at dataset 0
483 const fvMesh& mesh = *meshPtr_;
487 Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
491 for (int partId = selector.start(); partId < selector.end(); ++partId)
493 const word partName = getPartName(partId);
495 if (!partStatus_[partId])
502 Info<< "Creating VTK mesh for faceSet=" << partName << endl;
505 const faceSet fSet(mesh, partName);
507 vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSet);
510 AddToBlock(output, vtkmesh, selector, datasetNo, partName);
513 partDataset_[partId] = datasetNo++;
525 Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
531 void Foam::vtkPV3Foam::convertMeshPointZones
533 vtkMultiBlockDataSet* output,
537 partInfo& selector = partInfoPointZones_;
538 selector.block(blockNo); // set output block
539 label datasetNo = 0; // restart at dataset 0
540 const fvMesh& mesh = *meshPtr_;
544 Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
550 const pointZoneMesh& zMesh = mesh.pointZones();
551 for (int partId = selector.start(); partId < selector.end(); ++partId)
553 word zoneName = getPartName(partId);
554 label zoneId = zMesh.findZoneID(zoneName);
556 if (!partStatus_[partId] || zoneId < 0)
561 vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
564 AddToBlock(output, vtkmesh, selector, datasetNo, zoneName);
567 partDataset_[partId] = datasetNo++;
580 Info<< "<end> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
587 void Foam::vtkPV3Foam::convertMeshPointSets
589 vtkMultiBlockDataSet* output,
593 partInfo& selector = partInfoPointSets_;
594 selector.block(blockNo); // set output block
595 label datasetNo = 0; // restart at dataset 0
596 const fvMesh& mesh = *meshPtr_;
600 Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
604 for (int partId = selector.start(); partId < selector.end(); ++partId)
606 word partName = getPartName(partId);
608 if (!partStatus_[partId])
615 Info<< "Creating VTK mesh for pointSet=" << partName << endl;
618 const pointSet pSet(mesh, partName);
620 vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSet);
623 AddToBlock(output, vtkmesh, selector, datasetNo, partName);
626 partDataset_[partId] = datasetNo++;
638 Info<< "<end> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
643 // ************************************************************************* //