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-2021 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  {
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  printMemory();
99  }
100 }
101 
102 
103 void Foam::vtkPVFoam::convertMeshLagrangian
104 (
105  vtkMultiBlockDataSet* output,
106  int& blockNo
107 )
108 {
109  arrayRange& range = arrayRangeLagrangian_;
110  range.block(blockNo); // Set output block
111  label datasetNo = 0; // Restart at dataset 0
112  const fvMesh& mesh = *meshPtr_;
113 
114  if (debug)
115  {
116  InfoInFunction << endl;
117  printMemory();
118  }
119 
120  for (int partId = range.start(); partId < range.end(); ++partId)
121  {
122  const word cloudName = getPartName(partId);
123 
124  if (!partStatus_[partId])
125  {
126  continue;
127  }
128 
129  vtkPolyData* vtkmesh = lagrangianVTKMesh(mesh, cloudName);
130 
131  if (vtkmesh)
132  {
133  AddToBlock(output, vtkmesh, range, datasetNo, cloudName);
134  vtkmesh->Delete();
135 
136  partDataset_[partId] = datasetNo++;
137  }
138  }
139 
140  // Anything added?
141  if (datasetNo)
142  {
143  ++blockNo;
144  }
145 
146  if (debug)
147  {
148  printMemory();
149  }
150 }
151 
152 
153 void Foam::vtkPVFoam::convertMeshPatches
154 (
155  vtkMultiBlockDataSet* output,
156  int& blockNo
157 )
158 {
159  arrayRange& range = arrayRangePatches_;
160  range.block(blockNo); // Set output block
161  label datasetNo = 0; // Restart at dataset 0
162  const fvMesh& mesh = *meshPtr_;
163  const polyBoundaryMesh& patches = mesh.boundaryMesh();
164 
165  if (debug)
166  {
167  InfoInFunction << endl;
168  printMemory();
169  }
170 
171  for (int partId = range.start(); partId < range.end(); ++partId)
172  {
173  if (!partStatus_[partId])
174  {
175  continue;
176  }
177 
178  const word patchName = getPartName(partId);
179 
181  patchIds(patches.patchSet(List<wordRe>(1, wordRe(patchName))));
182 
183  if (debug)
184  {
186  << "Creating VTK mesh for patches [" << patchIds <<"] "
187  << patchName << endl;
188  }
189 
190  vtkPolyData* vtkmesh = nullptr;
191  if (patchIds.size() == 1)
192  {
193  vtkmesh = patchVTKMesh(patchName, patches[patchIds.begin().key()]);
194  }
195  else
196  {
197  // Patch group. Collect patch faces.
198  label sz = 0;
199  forAllConstIter(labelHashSet, patchIds, iter)
200  {
201  sz += patches[iter.key()].size();
202  }
203  labelList meshFaceLabels(sz);
204  sz = 0;
205  forAllConstIter(labelHashSet, patchIds, iter)
206  {
207  const polyPatch& pp = patches[iter.key()];
208  forAll(pp, i)
209  {
210  meshFaceLabels[sz++] = pp.start()+i;
211  }
212  }
213  UIndirectList<face> fcs(mesh.faces(), meshFaceLabels);
214  uindirectPrimitivePatch pp(fcs, mesh.points());
215 
216  vtkmesh = patchVTKMesh(patchName, pp);
217  }
218 
219 
220  if (vtkmesh)
221  {
222  AddToBlock(output, vtkmesh, range, datasetNo, patchName);
223  vtkmesh->Delete();
224 
225  partDataset_[partId] = datasetNo++;
226  }
227  }
228 
229  // Anything added?
230  if (datasetNo)
231  {
232  ++blockNo;
233  }
234 
235  if (debug)
236  {
237  printMemory();
238  }
239 }
240 
241 
242 void Foam::vtkPVFoam::convertMeshCellZones
243 (
244  vtkMultiBlockDataSet* output,
245  int& blockNo
246 )
247 {
248  arrayRange& range = arrayRangeCellZones_;
249  range.block(blockNo); // Set output block
250  label datasetNo = 0; // Restart at dataset 0
251  const fvMesh& mesh = *meshPtr_;
252 
253  // Resize for decomposed polyhedra
254  zonePolyDecomp_.setSize(range.size());
255 
256  if (range.empty())
257  {
258  return;
259  }
260 
261  if (debug)
262  {
263  InfoInFunction << endl;
264  printMemory();
265  }
266 
267  const meshCellZones& zMesh = mesh.cellZones();
268  for (int partId = range.start(); partId < range.end(); ++partId)
269  {
270  const word zoneName = getPartName(partId);
271  const label zoneId = zMesh.findZoneID(zoneName);
272 
273  if (!partStatus_[partId] || zoneId < 0)
274  {
275  continue;
276  }
277 
278  if (debug)
279  {
281  << "Creating VTK mesh for cellZone[" << zoneId << "] "
282  << zoneName << endl;
283  }
284 
285  fvMeshSubset subsetter(mesh);
286  subsetter.setLargeCellSubset(zMesh[zoneId]);
287 
288  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
289  (
290  subsetter.subMesh(),
291  zonePolyDecomp_[datasetNo]
292  );
293 
294  if (vtkmesh)
295  {
296  // superCells + addPointCellLabels must contain global cell ids
298  (
299  subsetter.cellMap(),
300  zonePolyDecomp_[datasetNo].superCells()
301  );
303  (
304  subsetter.cellMap(),
305  zonePolyDecomp_[datasetNo].addPointCellLabels()
306  );
307 
308  // Copy pointMap as well, otherwise pointFields fail
309  zonePolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
310 
311  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
312  vtkmesh->Delete();
313 
314  partDataset_[partId] = datasetNo++;
315  }
316  }
317 
318  // Anything added?
319  if (datasetNo)
320  {
321  ++blockNo;
322  }
323 
324  if (debug)
325  {
326  printMemory();
327  }
328 }
329 
330 
331 void Foam::vtkPVFoam::convertMeshCellSets
332 (
333  vtkMultiBlockDataSet* output,
334  int& blockNo
335 )
336 {
337  arrayRange& range = arrayRangeCellSets_;
338  range.block(blockNo); // Set output block
339  label datasetNo = 0; // Restart at dataset 0
340  const fvMesh& mesh = *meshPtr_;
341 
342  // Resize for decomposed polyhedra
343  csetPolyDecomp_.setSize(range.size());
344 
345  if (debug)
346  {
347  InfoInFunction << endl;
348  printMemory();
349  }
350 
351  for (int partId = range.start(); partId < range.end(); ++partId)
352  {
353  const word partName = getPartName(partId);
354 
355  if (!partStatus_[partId])
356  {
357  continue;
358  }
359 
360  if (debug)
361  {
363  << "Creating VTK mesh for cellSet=" << partName << endl;
364  }
365 
366  const cellSet cSet(mesh, partName);
367  fvMeshSubset subsetter(mesh);
368  subsetter.setLargeCellSubset(cSet);
369 
370  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
371  (
372  subsetter.subMesh(),
373  csetPolyDecomp_[datasetNo]
374  );
375 
376  if (vtkmesh)
377  {
378  // superCells + addPointCellLabels must contain global cell ids
380  (
381  subsetter.cellMap(),
382  csetPolyDecomp_[datasetNo].superCells()
383  );
385  (
386  subsetter.cellMap(),
387  csetPolyDecomp_[datasetNo].addPointCellLabels()
388  );
389 
390  // Copy pointMap as well, otherwise pointFields fail
391  csetPolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
392 
393  AddToBlock(output, vtkmesh, range, datasetNo, partName);
394  vtkmesh->Delete();
395 
396  partDataset_[partId] = datasetNo++;
397  }
398  }
399 
400  // Anything added?
401  if (datasetNo)
402  {
403  ++blockNo;
404  }
405 
406  if (debug)
407  {
408  printMemory();
409  }
410 }
411 
412 
413 void Foam::vtkPVFoam::convertMeshFaceZones
414 (
415  vtkMultiBlockDataSet* output,
416  int& blockNo
417 )
418 {
419  arrayRange& range = arrayRangeFaceZones_;
420  range.block(blockNo); // Set output block
421  label datasetNo = 0; // Restart at dataset 0
422  const fvMesh& mesh = *meshPtr_;
423 
424  if (range.empty())
425  {
426  return;
427  }
428 
429  if (debug)
430  {
431  InfoInFunction << endl;
432  printMemory();
433  }
434 
435  const meshFaceZones& zMesh = mesh.faceZones();
436  for (int partId = range.start(); partId < range.end(); ++partId)
437  {
438  const word zoneName = getPartName(partId);
439  const label zoneId = zMesh.findZoneID(zoneName);
440 
441  if (!partStatus_[partId] || zoneId < 0)
442  {
443  continue;
444  }
445 
446  if (debug)
447  {
449  << "Creating VTKmesh for faceZone[" << zoneId << "] "
450  << zoneName << endl;
451  }
452 
453  vtkPolyData* vtkmesh = patchVTKMesh(zoneName, zMesh[zoneId]());
454 
455  if (vtkmesh)
456  {
457  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
458  vtkmesh->Delete();
459 
460  partDataset_[partId] = datasetNo++;
461  }
462  }
463 
464  // Anything added?
465  if (datasetNo)
466  {
467  ++blockNo;
468  }
469 
470  if (debug)
471  {
472  printMemory();
473  }
474 }
475 
476 
477 void Foam::vtkPVFoam::convertMeshFaceSets
478 (
479  vtkMultiBlockDataSet* output,
480  int& blockNo
481 )
482 {
483  arrayRange& range = arrayRangeFaceSets_;
484  range.block(blockNo); // Set output block
485  label datasetNo = 0; // Restart at dataset 0
486  const fvMesh& mesh = *meshPtr_;
487 
488  if (debug)
489  {
490  InfoInFunction << endl;
491  printMemory();
492  }
493 
494  for (int partId = range.start(); partId < range.end(); ++partId)
495  {
496  const word partName = getPartName(partId);
497 
498  if (!partStatus_[partId])
499  {
500  continue;
501  }
502 
503  if (debug)
504  {
506  << "Creating VTK mesh for faceSet=" << partName << endl;
507  }
508 
509  const faceSet fSet(mesh, partName);
510 
511  vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSet);
512  if (vtkmesh)
513  {
514  AddToBlock(output, vtkmesh, range, datasetNo, partName);
515  vtkmesh->Delete();
516 
517  partDataset_[partId] = datasetNo++;
518  }
519  }
520 
521  // Anything added?
522  if (datasetNo)
523  {
524  ++blockNo;
525  }
526 
527  if (debug)
528  {
529  printMemory();
530  }
531 }
532 
533 
534 void Foam::vtkPVFoam::convertMeshPointZones
535 (
536  vtkMultiBlockDataSet* output,
537  int& blockNo
538 )
539 {
540  arrayRange& range = arrayRangePointZones_;
541  range.block(blockNo); // Set output block
542  label datasetNo = 0; // Restart at dataset 0
543  const fvMesh& mesh = *meshPtr_;
544 
545  if (debug)
546  {
547  InfoInFunction << endl;
548  printMemory();
549  }
550 
551  if (range.size())
552  {
553  const meshPointZones& zMesh = mesh.pointZones();
554  for (int partId = range.start(); partId < range.end(); ++partId)
555  {
556  word zoneName = getPartName(partId);
557  label zoneId = zMesh.findZoneID(zoneName);
558 
559  if (!partStatus_[partId] || zoneId < 0)
560  {
561  continue;
562  }
563 
564  vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
565  if (vtkmesh)
566  {
567  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
568  vtkmesh->Delete();
569 
570  partDataset_[partId] = datasetNo++;
571  }
572  }
573  }
574 
575  // Anything added?
576  if (datasetNo)
577  {
578  ++blockNo;
579  }
580 
581  if (debug)
582  {
583  printMemory();
584  }
585 }
586 
587 
588 
589 void Foam::vtkPVFoam::convertMeshPointSets
590 (
591  vtkMultiBlockDataSet* output,
592  int& blockNo
593 )
594 {
595  arrayRange& range = arrayRangePointSets_;
596  range.block(blockNo); // Set output block
597  label datasetNo = 0; // Restart at dataset 0
598  const fvMesh& mesh = *meshPtr_;
599 
600  if (debug)
601  {
602  InfoInFunction << endl;
603  printMemory();
604  }
605 
606  for (int partId = range.start(); partId < range.end(); ++partId)
607  {
608  word partName = getPartName(partId);
609 
610  if (!partStatus_[partId])
611  {
612  continue;
613  }
614 
615  if (debug)
616  {
618  << "Creating VTK mesh for pointSet=" << partName << endl;
619  }
620 
621  const pointSet pSet(mesh, partName);
622 
623  vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSet);
624  if (vtkmesh)
625  {
626  AddToBlock(output, vtkmesh, range, datasetNo, partName);
627  vtkmesh->Delete();
628 
629  partDataset_[partId] = datasetNo++;
630  }
631  }
632 
633  // Anything added?
634  if (datasetNo)
635  {
636  ++blockNo;
637  }
638 
639  if (debug)
640  {
641  printMemory();
642  }
643 }
644 
645 
646 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
void setSize(const label)
Reset size of List.
Definition: List.C:281
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static void printMemory()
Simple memory used debugging information.
const fvPatchList & patches
#define InfoInFunction
Report an information message using Foam::Info.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
MeshZones< cellZone, polyMesh > meshCellZones
A MeshZones with the type cellZone.
MeshZones< pointZone, polyMesh > meshPointZones
A MeshZones with the type pointZone.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
Foam::uindirectPrimitivePatch.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
MeshZones< faceZone, polyMesh > meshFaceZones
A MeshZones with the type faceZone.
const word cloudName(propsDict.lookup("cloudName"))