vtkPVFoamMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "vtkPVFoam.H"
27 
28 // OpenFOAM includes
29 #include "cellSet.H"
30 #include "faceSet.H"
31 #include "pointSet.H"
32 #include "fvMeshSubset.H"
33 #include "vtkPVFoamReader.h"
35 
36 // VTK includes
37 #include "vtkDataArraySelection.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkPolyData.h"
40 #include "vtkUnstructuredGrid.h"
41 
42 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
43 
44 void Foam::vtkPVFoam::convertMeshVolume
45 (
46  vtkMultiBlockDataSet* output,
47  int& blockNo
48 )
49 {
50  arrayRange& range = arrayRangeVolume_;
51  range.block(blockNo); // set output block
52  label datasetNo = 0; // restart at dataset 0
53  const fvMesh& mesh = *meshPtr_;
54 
55  // resize for decomposed polyhedra
56  regionPolyDecomp_.setSize(range.size());
57 
58  if (debug)
59  {
60  Info<< "<beg> Foam::vtkPVFoam::convertMeshVolume" << endl;
61  printMemory();
62  }
63 
64  // Convert the internalMesh
65  // this looks like more than one part, but it isn't
66  for (int partId = range.start(); partId < range.end(); ++partId)
67  {
68  const word partName = "internalMesh";
69 
70  if (!partStatus_[partId])
71  {
72  continue;
73  }
74 
75  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
76  (
77  mesh,
78  regionPolyDecomp_[datasetNo]
79  );
80 
81  if (vtkmesh)
82  {
83  AddToBlock(output, vtkmesh, range, datasetNo, partName);
84  vtkmesh->Delete();
85 
86  partDataset_[partId] = datasetNo++;
87  }
88  }
89 
90  // anything added?
91  if (datasetNo)
92  {
93  ++blockNo;
94  }
95 
96  if (debug)
97  {
98  Info<< "<end> Foam::vtkPVFoam::convertMeshVolume" << endl;
99  printMemory();
100  }
101 }
102 
103 
104 void Foam::vtkPVFoam::convertMeshLagrangian
105 (
106  vtkMultiBlockDataSet* output,
107  int& blockNo
108 )
109 {
110  arrayRange& range = arrayRangeLagrangian_;
111  range.block(blockNo); // set output block
112  label datasetNo = 0; // restart at dataset 0
113  const fvMesh& mesh = *meshPtr_;
114 
115  if (debug)
116  {
117  Info<< "<beg> Foam::vtkPVFoam::convertMeshLagrangian" << endl;
118  printMemory();
119  }
120 
121  for (int partId = range.start(); partId < range.end(); ++partId)
122  {
123  const word cloudName = getPartName(partId);
124 
125  if (!partStatus_[partId])
126  {
127  continue;
128  }
129 
130  vtkPolyData* vtkmesh = lagrangianVTKMesh(mesh, cloudName);
131 
132  if (vtkmesh)
133  {
134  AddToBlock(output, vtkmesh, range, datasetNo, cloudName);
135  vtkmesh->Delete();
136 
137  partDataset_[partId] = datasetNo++;
138  }
139  }
140 
141  // anything added?
142  if (datasetNo)
143  {
144  ++blockNo;
145  }
146 
147  if (debug)
148  {
149  Info<< "<end> Foam::vtkPVFoam::convertMeshLagrangian" << endl;
150  printMemory();
151  }
152 }
153 
154 
155 void Foam::vtkPVFoam::convertMeshPatches
156 (
157  vtkMultiBlockDataSet* output,
158  int& blockNo
159 )
160 {
161  arrayRange& range = arrayRangePatches_;
162  range.block(blockNo); // set output block
163  label datasetNo = 0; // restart at dataset 0
164  const fvMesh& mesh = *meshPtr_;
165  const polyBoundaryMesh& patches = mesh.boundaryMesh();
166 
167  if (debug)
168  {
169  Info<< "<beg> Foam::vtkPVFoam::convertMeshPatches" << endl;
170  printMemory();
171  }
172 
173  for (int partId = range.start(); partId < range.end(); ++partId)
174  {
175  if (!partStatus_[partId])
176  {
177  continue;
178  }
179 
180  const word patchName = getPartName(partId);
181 
183  patchIds(patches.patchSet(List<wordRe>(1, wordRe(patchName))));
184 
185  if (debug)
186  {
187  Info<< "Creating VTK mesh for patches [" << patchIds <<"] "
188  << patchName << endl;
189  }
190 
191  vtkPolyData* vtkmesh = nullptr;
192  if (patchIds.size() == 1)
193  {
194  vtkmesh = patchVTKMesh(patchName, patches[patchIds.begin().key()]);
195  }
196  else
197  {
198  // Patch group. Collect patch faces.
199  label sz = 0;
200  forAllConstIter(labelHashSet, patchIds, iter)
201  {
202  sz += patches[iter.key()].size();
203  }
204  labelList meshFaceLabels(sz);
205  sz = 0;
206  forAllConstIter(labelHashSet, patchIds, iter)
207  {
208  const polyPatch& pp = patches[iter.key()];
209  forAll(pp, i)
210  {
211  meshFaceLabels[sz++] = pp.start()+i;
212  }
213  }
214  UIndirectList<face> fcs(mesh.faces(), meshFaceLabels);
215  uindirectPrimitivePatch pp(fcs, mesh.points());
216 
217  vtkmesh = patchVTKMesh(patchName, pp);
218  }
219 
220 
221  if (vtkmesh)
222  {
223  AddToBlock(output, vtkmesh, range, datasetNo, patchName);
224  vtkmesh->Delete();
225 
226  partDataset_[partId] = datasetNo++;
227  }
228  }
229 
230  // anything added?
231  if (datasetNo)
232  {
233  ++blockNo;
234  }
235 
236  if (debug)
237  {
238  Info<< "<end> Foam::vtkPVFoam::convertMeshPatches" << endl;
239  printMemory();
240  }
241 }
242 
243 
244 void Foam::vtkPVFoam::convertMeshCellZones
245 (
246  vtkMultiBlockDataSet* output,
247  int& blockNo
248 )
249 {
250  arrayRange& range = arrayRangeCellZones_;
251  range.block(blockNo); // set output block
252  label datasetNo = 0; // restart at dataset 0
253  const fvMesh& mesh = *meshPtr_;
254 
255  // resize for decomposed polyhedra
256  zonePolyDecomp_.setSize(range.size());
257 
258  if (range.empty())
259  {
260  return;
261  }
262 
263  if (debug)
264  {
265  Info<< "<beg> Foam::vtkPVFoam::convertMeshCellZones" << endl;
266  printMemory();
267  }
268 
269  const cellZoneMesh& zMesh = mesh.cellZones();
270  for (int partId = range.start(); partId < range.end(); ++partId)
271  {
272  const word zoneName = getPartName(partId);
273  const label zoneId = zMesh.findZoneID(zoneName);
274 
275  if (!partStatus_[partId] || zoneId < 0)
276  {
277  continue;
278  }
279 
280  if (debug)
281  {
282  Info<< "Creating VTK mesh for cellZone[" << zoneId << "] "
283  << zoneName << endl;
284  }
285 
286  fvMeshSubset subsetter(mesh);
287  subsetter.setLargeCellSubset(zMesh[zoneId]);
288 
289  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
290  (
291  subsetter.subMesh(),
292  zonePolyDecomp_[datasetNo]
293  );
294 
295  if (vtkmesh)
296  {
297  // superCells + addPointCellLabels must contain global cell ids
299  (
300  subsetter.cellMap(),
301  zonePolyDecomp_[datasetNo].superCells()
302  );
304  (
305  subsetter.cellMap(),
306  zonePolyDecomp_[datasetNo].addPointCellLabels()
307  );
308 
309  // copy pointMap as well, otherwise pointFields fail
310  zonePolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
311 
312  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
313  vtkmesh->Delete();
314 
315  partDataset_[partId] = datasetNo++;
316  }
317  }
318 
319  // anything added?
320  if (datasetNo)
321  {
322  ++blockNo;
323  }
324 
325  if (debug)
326  {
327  Info<< "<end> Foam::vtkPVFoam::convertMeshCellZones" << endl;
328  printMemory();
329  }
330 }
331 
332 
333 void Foam::vtkPVFoam::convertMeshCellSets
334 (
335  vtkMultiBlockDataSet* output,
336  int& blockNo
337 )
338 {
339  arrayRange& range = arrayRangeCellSets_;
340  range.block(blockNo); // set output block
341  label datasetNo = 0; // restart at dataset 0
342  const fvMesh& mesh = *meshPtr_;
343 
344  // resize for decomposed polyhedra
345  csetPolyDecomp_.setSize(range.size());
346 
347  if (debug)
348  {
349  Info<< "<beg> Foam::vtkPVFoam::convertMeshCellSets" << endl;
350  printMemory();
351  }
352 
353  for (int partId = range.start(); partId < range.end(); ++partId)
354  {
355  const word partName = getPartName(partId);
356 
357  if (!partStatus_[partId])
358  {
359  continue;
360  }
361 
362  if (debug)
363  {
364  Info<< "Creating VTK mesh for cellSet=" << partName << endl;
365  }
366 
367  const cellSet cSet(mesh, partName);
368  fvMeshSubset subsetter(mesh);
369  subsetter.setLargeCellSubset(cSet);
370 
371  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
372  (
373  subsetter.subMesh(),
374  csetPolyDecomp_[datasetNo]
375  );
376 
377  if (vtkmesh)
378  {
379  // superCells + addPointCellLabels must contain global cell ids
381  (
382  subsetter.cellMap(),
383  csetPolyDecomp_[datasetNo].superCells()
384  );
386  (
387  subsetter.cellMap(),
388  csetPolyDecomp_[datasetNo].addPointCellLabels()
389  );
390 
391  // copy pointMap as well, otherwise pointFields fail
392  csetPolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
393 
394  AddToBlock(output, vtkmesh, range, datasetNo, partName);
395  vtkmesh->Delete();
396 
397  partDataset_[partId] = datasetNo++;
398  }
399  }
400 
401  // anything added?
402  if (datasetNo)
403  {
404  ++blockNo;
405  }
406 
407  if (debug)
408  {
409  Info<< "<end> Foam::vtkPVFoam::convertMeshCellSets" << endl;
410  printMemory();
411  }
412 }
413 
414 
415 void Foam::vtkPVFoam::convertMeshFaceZones
416 (
417  vtkMultiBlockDataSet* output,
418  int& blockNo
419 )
420 {
421  arrayRange& range = arrayRangeFaceZones_;
422  range.block(blockNo); // set output block
423  label datasetNo = 0; // restart at dataset 0
424  const fvMesh& mesh = *meshPtr_;
425 
426  if (range.empty())
427  {
428  return;
429  }
430 
431  if (debug)
432  {
433  Info<< "<beg> Foam::vtkPVFoam::convertMeshFaceZones" << endl;
434  printMemory();
435  }
436 
437  const faceZoneMesh& zMesh = mesh.faceZones();
438  for (int partId = range.start(); partId < range.end(); ++partId)
439  {
440  const word zoneName = getPartName(partId);
441  const label zoneId = zMesh.findZoneID(zoneName);
442 
443  if (!partStatus_[partId] || zoneId < 0)
444  {
445  continue;
446  }
447 
448  if (debug)
449  {
450  Info<< "Creating VTKmesh for faceZone[" << zoneId << "] "
451  << zoneName << endl;
452  }
453 
454  vtkPolyData* vtkmesh = patchVTKMesh(zoneName, zMesh[zoneId]());
455 
456  if (vtkmesh)
457  {
458  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
459  vtkmesh->Delete();
460 
461  partDataset_[partId] = datasetNo++;
462  }
463  }
464 
465  // anything added?
466  if (datasetNo)
467  {
468  ++blockNo;
469  }
470 
471  if (debug)
472  {
473  Info<< "<end> Foam::vtkPVFoam::convertMeshFaceZones" << endl;
474  printMemory();
475  }
476 }
477 
478 
479 void Foam::vtkPVFoam::convertMeshFaceSets
480 (
481  vtkMultiBlockDataSet* output,
482  int& blockNo
483 )
484 {
485  arrayRange& range = arrayRangeFaceSets_;
486  range.block(blockNo); // set output block
487  label datasetNo = 0; // restart at dataset 0
488  const fvMesh& mesh = *meshPtr_;
489 
490  if (debug)
491  {
492  Info<< "<beg> Foam::vtkPVFoam::convertMeshFaceSets" << endl;
493  printMemory();
494  }
495 
496  for (int partId = range.start(); partId < range.end(); ++partId)
497  {
498  const word partName = getPartName(partId);
499 
500  if (!partStatus_[partId])
501  {
502  continue;
503  }
504 
505  if (debug)
506  {
507  Info<< "Creating VTK mesh for faceSet=" << partName << endl;
508  }
509 
510  const faceSet fSet(mesh, partName);
511 
512  vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSet);
513  if (vtkmesh)
514  {
515  AddToBlock(output, vtkmesh, range, datasetNo, partName);
516  vtkmesh->Delete();
517 
518  partDataset_[partId] = datasetNo++;
519  }
520  }
521 
522  // anything added?
523  if (datasetNo)
524  {
525  ++blockNo;
526  }
527 
528  if (debug)
529  {
530  Info<< "<end> Foam::vtkPVFoam::convertMeshFaceSets" << endl;
531  printMemory();
532  }
533 }
534 
535 
536 void Foam::vtkPVFoam::convertMeshPointZones
537 (
538  vtkMultiBlockDataSet* output,
539  int& blockNo
540 )
541 {
542  arrayRange& range = arrayRangePointZones_;
543  range.block(blockNo); // set output block
544  label datasetNo = 0; // restart at dataset 0
545  const fvMesh& mesh = *meshPtr_;
546 
547  if (debug)
548  {
549  Info<< "<beg> Foam::vtkPVFoam::convertMeshPointZones" << endl;
550  printMemory();
551  }
552 
553  if (range.size())
554  {
555  const pointZoneMesh& zMesh = mesh.pointZones();
556  for (int partId = range.start(); partId < range.end(); ++partId)
557  {
558  word zoneName = getPartName(partId);
559  label zoneId = zMesh.findZoneID(zoneName);
560 
561  if (!partStatus_[partId] || zoneId < 0)
562  {
563  continue;
564  }
565 
566  vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
567  if (vtkmesh)
568  {
569  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
570  vtkmesh->Delete();
571 
572  partDataset_[partId] = datasetNo++;
573  }
574  }
575  }
576 
577  // anything added?
578  if (datasetNo)
579  {
580  ++blockNo;
581  }
582 
583  if (debug)
584  {
585  Info<< "<end> Foam::vtkPVFoam::convertMeshPointZones" << endl;
586  printMemory();
587  }
588 }
589 
590 
591 
592 void Foam::vtkPVFoam::convertMeshPointSets
593 (
594  vtkMultiBlockDataSet* output,
595  int& blockNo
596 )
597 {
598  arrayRange& range = arrayRangePointSets_;
599  range.block(blockNo); // set output block
600  label datasetNo = 0; // restart at dataset 0
601  const fvMesh& mesh = *meshPtr_;
602 
603  if (debug)
604  {
605  Info<< "<beg> Foam::vtkPVFoam::convertMeshPointSets" << endl;
606  printMemory();
607  }
608 
609  for (int partId = range.start(); partId < range.end(); ++partId)
610  {
611  word partName = getPartName(partId);
612 
613  if (!partStatus_[partId])
614  {
615  continue;
616  }
617 
618  if (debug)
619  {
620  Info<< "Creating VTK mesh for pointSet=" << partName << endl;
621  }
622 
623  const pointSet pSet(mesh, partName);
624 
625  vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSet);
626  if (vtkmesh)
627  {
628  AddToBlock(output, vtkmesh, range, datasetNo, partName);
629  vtkmesh->Delete();
630 
631  partDataset_[partId] = datasetNo++;
632  }
633  }
634 
635  // anything added?
636  if (datasetNo)
637  {
638  ++blockNo;
639  }
640 
641  if (debug)
642  {
643  Info<< "<end> Foam::vtkPVFoam::convertMeshPointSets" << endl;
644  printMemory();
645  }
646 }
647 
648 
649 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static void printMemory()
Simple memory used debugging information.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:208
List< label > labelList
A List of labels.
Definition: labelList.H:56
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
Foam::uindirectPrimitivePatch.
messageStream Info
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.