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-2025 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 #include "vtkPVFoamReader.h"
28 
29 // OpenFOAM includes
30 #include "domainDecomposition.H"
31 #include "cellSet.H"
32 #include "faceSet.H"
33 #include "pointSet.H"
34 #include "fvMeshSubset.H"
36 #include "LagrangianMesh.H"
39 
40 // VTK includes
41 #include "vtkDataArraySelection.h"
42 #include "vtkMultiBlockDataSet.h"
43 #include "vtkPolyData.h"
44 #include "vtkUnstructuredGrid.h"
45 
46 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
47 
48 void Foam::vtkPVFoam::convertMeshVolume
49 (
50  vtkMultiBlockDataSet* output,
51  int& blockNo
52 )
53 {
55 
56  arrayRange& range = arrayRangeVolume_;
57  range.block(blockNo); // Set output block
58  label datasetNo = 0; // Restart at dataset 0
59 
60  const fvMesh& mesh = procMeshesPtr_->completeMesh();
61 
62  // Resize for decomposed polyhedra
63  regionPolyDecomp_.setSize(range.size());
64 
65  // Convert the internalMesh
66  // This looks like more than one part, but it isn't
67  for (int partId = range.start(); partId < range.end(); ++partId)
68  {
69  const word partName = "internalMesh";
70 
71  if (!partStatus_[partId]) continue;
72 
73  vtkUnstructuredGrid* vtkmesh =
74  volumeVTKMesh(mesh, regionPolyDecomp_[datasetNo]);
75 
76  if (vtkmesh)
77  {
78  AddToBlock(output, vtkmesh, range, datasetNo, partName);
79  vtkmesh->Delete();
80 
81  partDataset_[partId] = datasetNo++;
82  }
83  }
84 
85  // Anything added?
86  if (datasetNo) ++blockNo;
87 }
88 
89 
90 void Foam::vtkPVFoam::convertMeshlagrangian
91 (
92  vtkMultiBlockDataSet* output,
93  int& blockNo
94 )
95 {
97 
98  arrayRange& range = arrayRangelagrangian_;
99  range.block(blockNo); // Set output block
100  label datasetNo = 0; // Restart at dataset 0
101 
102  lagrangianReconstructors_.clear();
103 
104  for (int partId = range.start(); partId < range.end(); ++partId)
105  {
106  const word lagrangianName = getPartName(partId);
107 
108  if (!partStatus_[partId]) continue;
109 
110  autoPtr<lagrangianFieldReconstructor> lreconstructorPtr;
111 
112  vtkPolyData* vtkmesh =
113  lagrangianVTKMesh(lagrangianName, lreconstructorPtr);
114  if (vtkmesh)
115  {
116  AddToBlock(output, vtkmesh, range, datasetNo, lagrangianName);
117  vtkmesh->Delete();
118 
119  lagrangianReconstructors_.append(lreconstructorPtr.ptr());
120 
121  partDataset_[partId] = datasetNo++;
122  }
123  }
124 
125  // Anything added?
126  if (datasetNo) ++blockNo;
127 }
128 
129 
130 void Foam::vtkPVFoam::convertMeshLagrangian
131 (
132  vtkMultiBlockDataSet* output,
133  int& blockNo
134 )
135 {
137 
138  arrayRange& range = arrayRangeLagrangian_;
139  range.block(blockNo); // Set output block
140  label datasetNo = 0; // Restart at dataset 0
141 
142  LagrangianMeshes_.clear();
143  LagrangianReconstructors_.clear();
144 
145  for (int partId = range.start(); partId < range.end(); ++partId)
146  {
147  const word LagrangianName = getPartName(partId);
148 
149  if (!partStatus_[partId]) continue;
150 
151  autoPtr<LagrangianMesh> LmeshPtr;
152  autoPtr<LagrangianFieldReconstructor> LreconstructorPtr;
153 
154  vtkPolyData* vtkmesh =
155  LagrangianVTKMesh(LagrangianName, LmeshPtr, LreconstructorPtr);
156  if (vtkmesh)
157  {
158  AddToBlock(output, vtkmesh, range, datasetNo, LagrangianName);
159  vtkmesh->Delete();
160 
161  LagrangianMeshes_.append(LmeshPtr.ptr());
162  LagrangianReconstructors_.append(LreconstructorPtr.ptr());
163 
164  partDataset_[partId] = datasetNo++;
165  }
166  }
167 
168  // Anything added?
169  if (datasetNo) ++blockNo;
170 }
171 
172 
173 void Foam::vtkPVFoam::convertMeshPatches
174 (
175  vtkMultiBlockDataSet* output,
176  int& blockNo
177 )
178 {
179  arrayRange& range = arrayRangePatches_;
180  range.block(blockNo); // Set output block
181  label datasetNo = 0; // Restart at dataset 0
182 
183  const fvMesh& mesh = procMeshesPtr_->completeMesh();
184  const polyBoundaryMesh& patches = mesh.boundaryMesh();
185 
186  for (int partId = range.start(); partId < range.end(); ++partId)
187  {
188  if (!partStatus_[partId]) continue;
189 
190  const word patchName = getPartName(partId);
191  const labelHashSet patchIds =
192  patches.patchSet(List<wordRe>(1, wordRe(patchName)));
193 
194  DebugInFunction << "Creating VTK mesh for patch(es)[";
195  forAllConstIter(labelHashSet, patchIds, iter)
196  {
197  if (iter != patchIds.begin()) DebugInfo<< ',';
198  DebugInfo<< iter.key();
199  }
200  DebugInfo<< "]: " << patchName << endl;
201 
202  vtkPolyData* vtkmesh = nullptr;
203  if (patchIds.size() == 1)
204  {
205  vtkmesh = patchVTKMesh(patches[patchIds.begin().key()]);
206  }
207  else
208  {
209  // Patch group. Collect patch faces.
210  label sz = 0;
211  forAllConstIter(labelHashSet, patchIds, iter)
212  {
213  sz += patches[iter.key()].size();
214  }
215  labelList meshFaceLabels(sz);
216  sz = 0;
217  forAllConstIter(labelHashSet, patchIds, iter)
218  {
219  const polyPatch& pp = patches[iter.key()];
220  forAll(pp, i)
221  {
222  meshFaceLabels[sz++] = pp.start()+i;
223  }
224  }
225  UIndirectList<face> fcs(mesh.faces(), meshFaceLabels);
227 
228  vtkmesh = patchVTKMesh(pp);
229  }
230 
231  if (vtkmesh)
232  {
233  AddToBlock(output, vtkmesh, range, datasetNo, patchName);
234  vtkmesh->Delete();
235 
236  partDataset_[partId] = datasetNo++;
237  }
238  }
239 
240  // Anything added?
241  if (datasetNo)
242  {
243  ++blockNo;
244  }
245 }
246 
247 
248 void Foam::vtkPVFoam::convertMeshCellZones
249 (
250  vtkMultiBlockDataSet* output,
251  int& blockNo
252 )
253 {
254  arrayRange& range = arrayRangeCellZones_;
255  range.block(blockNo); // Set output block
256  label datasetNo = 0; // Restart at dataset 0
257 
258  const fvMesh& mesh = procMeshesPtr_->completeMesh();
259 
260  // Resize for decomposed polyhedra
261  zonePolyDecomp_.setSize(range.size());
262 
263  if (range.empty()) return;
264 
265  const cellZoneList& zMesh = mesh.cellZones();
266  for (int partId = range.start(); partId < range.end(); ++partId)
267  {
268  const word zoneName = getPartName(partId);
269  const label zoneId = zMesh.findIndex(zoneName);
270 
271  if (!partStatus_[partId] || zoneId < 0) continue;
272 
274  << "Creating VTK mesh for cellZone[" << zoneId << "]: "
275  << zoneName << endl;
276 
277  fvMeshSubset subsetter(mesh);
278  subsetter.setLargeCellSubset(zMesh[zoneId]);
279 
280  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
281  (
282  subsetter.subMesh(),
283  zonePolyDecomp_[datasetNo]
284  );
285 
286  if (vtkmesh)
287  {
288  // superCells + addPointCellLabels must contain global cell ids
290  (
291  subsetter.cellMap(),
292  zonePolyDecomp_[datasetNo].superCells()
293  );
295  (
296  subsetter.cellMap(),
297  zonePolyDecomp_[datasetNo].addPointCellLabels()
298  );
299 
300  // Copy pointMap as well, otherwise pointFields fail
301  zonePolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
302 
303  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
304  vtkmesh->Delete();
305 
306  partDataset_[partId] = datasetNo++;
307  }
308  }
309 
310  // Anything added?
311  if (datasetNo) ++blockNo;
312 }
313 
314 
315 void Foam::vtkPVFoam::convertMeshCellSets
316 (
317  vtkMultiBlockDataSet* output,
318  int& blockNo
319 )
320 {
321  arrayRange& range = arrayRangeCellSets_;
322  range.block(blockNo); // Set output block
323  label datasetNo = 0; // Restart at dataset 0
324 
325  const fvMesh& mesh = procMeshesPtr_->completeMesh();
326 
327  // Resize for decomposed polyhedra
328  setPolyDecomp_.setSize(range.size());
329 
330  for (int partId = range.start(); partId < range.end(); ++partId)
331  {
332  const word partName = getPartName(partId);
333 
334  if (!partStatus_[partId]) continue;
335 
337  << "Creating VTK mesh for cellSet: [" << partName << endl;
338 
339  const autoPtr<cellSet> cSetPtr =
340  reader_->GetDecomposedCase()
341  ? procMeshesPtr_->reconstructSet<cellSet>(partName)
342  : autoPtr<cellSet>(new cellSet(mesh, partName));
343 
344  fvMeshSubset subsetter(mesh);
345  subsetter.setLargeCellSubset(cSetPtr());
346 
347  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
348  (
349  subsetter.subMesh(),
350  setPolyDecomp_[datasetNo]
351  );
352 
353  if (vtkmesh)
354  {
355  // superCells + addPointCellLabels must contain global cell ids
357  (
358  subsetter.cellMap(),
359  setPolyDecomp_[datasetNo].superCells()
360  );
362  (
363  subsetter.cellMap(),
364  setPolyDecomp_[datasetNo].addPointCellLabels()
365  );
366 
367  // Copy pointMap as well, otherwise pointFields fail
368  setPolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
369 
370  AddToBlock(output, vtkmesh, range, datasetNo, partName);
371  vtkmesh->Delete();
372 
373  partDataset_[partId] = datasetNo++;
374  }
375  }
376 
377  // Anything added?
378  if (datasetNo) ++blockNo;
379 }
380 
381 
382 void Foam::vtkPVFoam::convertMeshFaceZones
383 (
384  vtkMultiBlockDataSet* output,
385  int& blockNo
386 )
387 {
388  arrayRange& range = arrayRangeFaceZones_;
389  range.block(blockNo); // Set output block
390  label datasetNo = 0; // Restart at dataset 0
391 
392  const fvMesh& mesh = procMeshesPtr_->completeMesh();
393 
394  if (range.empty()) return;
395 
396  const faceZoneList& zMesh = mesh.faceZones();
397  for (int partId = range.start(); partId < range.end(); ++partId)
398  {
399  const word zoneName = getPartName(partId);
400  const label zoneId = zMesh.findIndex(zoneName);
401 
402  if (!partStatus_[partId] || zoneId < 0) continue;
403 
405  << "Creating VTK mesh for faceZone[" << zoneId << "]: "
406  << zoneName << endl;
407 
408  vtkPolyData* vtkmesh = patchVTKMesh(zMesh[zoneId].patch());
409  if (vtkmesh)
410  {
411  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
412  vtkmesh->Delete();
413 
414  partDataset_[partId] = datasetNo++;
415  }
416  }
417 
418  // Anything added?
419  if (datasetNo) ++blockNo;
420 }
421 
422 
423 void Foam::vtkPVFoam::convertMeshFaceSets
424 (
425  vtkMultiBlockDataSet* output,
426  int& blockNo
427 )
428 {
429  arrayRange& range = arrayRangeFaceSets_;
430  range.block(blockNo); // Set output block
431  label datasetNo = 0; // Restart at dataset 0
432 
433  const fvMesh& mesh = procMeshesPtr_->completeMesh();
434 
435  for (int partId = range.start(); partId < range.end(); ++partId)
436  {
437  const word partName = getPartName(partId);
438 
439  if (!partStatus_[partId]) continue;
440 
442  << "Creating VTK mesh for faceSet: " << partName << endl;
443 
444  const autoPtr<faceSet> fSetPtr =
445  reader_->GetDecomposedCase()
446  ? procMeshesPtr_->reconstructSet<faceSet>(partName)
447  : autoPtr<faceSet>(new faceSet(mesh, partName));
448 
449  vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSetPtr());
450  if (vtkmesh)
451  {
452  AddToBlock(output, vtkmesh, range, datasetNo, partName);
453  vtkmesh->Delete();
454 
455  partDataset_[partId] = datasetNo++;
456  }
457  }
458 
459  // Anything added?
460  if (datasetNo) ++blockNo;
461 }
462 
463 
464 void Foam::vtkPVFoam::convertMeshPointZones
465 (
466  vtkMultiBlockDataSet* output,
467  int& blockNo
468 )
469 {
470  arrayRange& range = arrayRangePointZones_;
471  range.block(blockNo); // Set output block
472  label datasetNo = 0; // Restart at dataset 0
473 
474  const fvMesh& mesh = procMeshesPtr_->completeMesh();
475 
476  if (range.empty()) return;
477 
478  const pointZoneList& zMesh = mesh.pointZones();
479  for (int partId = range.start(); partId < range.end(); ++partId)
480  {
481  const word zoneName = getPartName(partId);
482  const label zoneId = zMesh.findIndex(zoneName);
483 
484  if (!partStatus_[partId] || zoneId < 0) continue;
485 
487  << "Creating VTK mesh for pointZone[" << zoneId << "]: "
488  << zoneName << endl;
489 
490  vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
491  if (vtkmesh)
492  {
493  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
494  vtkmesh->Delete();
495 
496  partDataset_[partId] = datasetNo++;
497  }
498  }
499 
500  // Anything added?
501  if (datasetNo) ++blockNo;
502 }
503 
504 
505 
506 void Foam::vtkPVFoam::convertMeshPointSets
507 (
508  vtkMultiBlockDataSet* output,
509  int& blockNo
510 )
511 {
512  arrayRange& range = arrayRangePointSets_;
513  range.block(blockNo); // Set output block
514  label datasetNo = 0; // Restart at dataset 0
515 
516  const fvMesh& mesh = procMeshesPtr_->completeMesh();
517 
518  for (int partId = range.start(); partId < range.end(); ++partId)
519  {
520  word partName = getPartName(partId);
521 
522  if (!partStatus_[partId]) continue;
523 
525  << "Creating VTK mesh for pointSet: " << partName << endl;
526 
527  const autoPtr<pointSet> pSetPtr =
528  reader_->GetDecomposedCase()
529  ? procMeshesPtr_->reconstructSet<pointSet>(partName)
530  : autoPtr<pointSet>(new pointSet(mesh, partName));
531 
532  vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSetPtr());
533  if (vtkmesh)
534  {
535  AddToBlock(output, vtkmesh, range, datasetNo, partName);
536  vtkmesh->Delete();
537 
538  partDataset_[partId] = datasetNo++;
539  }
540  }
541 
542  // Anything added?
543  if (datasetNo) ++blockNo;
544 }
545 
546 
547 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
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
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:437
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const fvPatchList & patches
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
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:258
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:213