vtkPV3FoamMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 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 "vtkPV3Foam.H"
27 
28 // OpenFOAM includes
29 #include "cellSet.H"
30 #include "faceSet.H"
31 #include "pointSet.H"
32 #include "fvMeshSubset.H"
33 #include "vtkPV3FoamReader.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::vtkPV3Foam::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::vtkPV3Foam::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::vtkPV3Foam::convertMeshVolume" << endl;
99  printMemory();
100  }
101 }
102 
103 
104 void Foam::vtkPV3Foam::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::vtkPV3Foam::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::vtkPV3Foam::convertMeshLagrangian" << endl;
150  printMemory();
151  }
152 }
153 
154 
155 void Foam::vtkPV3Foam::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::vtkPV3Foam::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 
182  labelHashSet patchIds
183  (
184  patches.patchSet(List<wordRe>(1, wordRe(patchName)))
185  );
186 
187  if (debug)
188  {
189  Info<< "Creating VTK mesh for patches [" << patchIds <<"] "
190  << patchName << endl;
191  }
192 
193  vtkPolyData* vtkmesh = NULL;
194  if (patchIds.size() == 1)
195  {
196  vtkmesh = patchVTKMesh(patchName, patches[patchIds.begin().key()]);
197  }
198  else
199  {
200  // Patch group. Collect patch faces.
201  label sz = 0;
202  forAllConstIter(labelHashSet, patchIds, iter)
203  {
204  sz += patches[iter.key()].size();
205  }
206  labelList meshFaceLabels(sz);
207  sz = 0;
208  forAllConstIter(labelHashSet, patchIds, iter)
209  {
210  const polyPatch& pp = patches[iter.key()];
211  forAll(pp, i)
212  {
213  meshFaceLabels[sz++] = pp.start()+i;
214  }
215  }
216  UIndirectList<face> fcs(mesh.faces(), meshFaceLabels);
217  uindirectPrimitivePatch pp(fcs, mesh.points());
218 
219  vtkmesh = patchVTKMesh(patchName, pp);
220  }
221 
222 
223  if (vtkmesh)
224  {
225  AddToBlock(output, vtkmesh, range, datasetNo, patchName);
226  vtkmesh->Delete();
227 
228  partDataset_[partId] = datasetNo++;
229  }
230  }
231 
232  // anything added?
233  if (datasetNo)
234  {
235  ++blockNo;
236  }
237 
238  if (debug)
239  {
240  Info<< "<end> Foam::vtkPV3Foam::convertMeshPatches" << endl;
241  printMemory();
242  }
243 }
244 
245 
246 void Foam::vtkPV3Foam::convertMeshCellZones
247 (
248  vtkMultiBlockDataSet* output,
249  int& blockNo
250 )
251 {
252  arrayRange& range = arrayRangeCellZones_;
253  range.block(blockNo); // set output block
254  label datasetNo = 0; // restart at dataset 0
255  const fvMesh& mesh = *meshPtr_;
256 
257  // resize for decomposed polyhedra
258  zonePolyDecomp_.setSize(range.size());
259 
260  if (range.empty())
261  {
262  return;
263  }
264 
265  if (debug)
266  {
267  Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
268  printMemory();
269  }
270 
271  const cellZoneMesh& zMesh = mesh.cellZones();
272  for (int partId = range.start(); partId < range.end(); ++partId)
273  {
274  const word zoneName = getPartName(partId);
275  const label zoneId = zMesh.findZoneID(zoneName);
276 
277  if (!partStatus_[partId] || zoneId < 0)
278  {
279  continue;
280  }
281 
282  if (debug)
283  {
284  Info<< "Creating VTK mesh for cellZone[" << zoneId << "] "
285  << zoneName << endl;
286  }
287 
288  fvMeshSubset subsetter(mesh);
289  subsetter.setLargeCellSubset(zMesh[zoneId]);
290 
291  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
292  (
293  subsetter.subMesh(),
294  zonePolyDecomp_[datasetNo]
295  );
296 
297  if (vtkmesh)
298  {
299  // superCells + addPointCellLabels must contain global cell ids
301  (
302  subsetter.cellMap(),
303  zonePolyDecomp_[datasetNo].superCells()
304  );
306  (
307  subsetter.cellMap(),
308  zonePolyDecomp_[datasetNo].addPointCellLabels()
309  );
310 
311  // copy pointMap as well, otherwise pointFields fail
312  zonePolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
313 
314  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
315  vtkmesh->Delete();
316 
317  partDataset_[partId] = datasetNo++;
318  }
319  }
320 
321  // anything added?
322  if (datasetNo)
323  {
324  ++blockNo;
325  }
326 
327  if (debug)
328  {
329  Info<< "<end> Foam::vtkPV3Foam::convertMeshCellZones" << endl;
330  printMemory();
331  }
332 }
333 
334 
335 void Foam::vtkPV3Foam::convertMeshCellSets
336 (
337  vtkMultiBlockDataSet* output,
338  int& blockNo
339 )
340 {
341  arrayRange& range = arrayRangeCellSets_;
342  range.block(blockNo); // set output block
343  label datasetNo = 0; // restart at dataset 0
344  const fvMesh& mesh = *meshPtr_;
345 
346  // resize for decomposed polyhedra
347  csetPolyDecomp_.setSize(range.size());
348 
349  if (debug)
350  {
351  Info<< "<beg> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
352  printMemory();
353  }
354 
355  for (int partId = range.start(); partId < range.end(); ++partId)
356  {
357  const word partName = getPartName(partId);
358 
359  if (!partStatus_[partId])
360  {
361  continue;
362  }
363 
364  if (debug)
365  {
366  Info<< "Creating VTK mesh for cellSet=" << partName << endl;
367  }
368 
369  const cellSet cSet(mesh, partName);
370  fvMeshSubset subsetter(mesh);
371  subsetter.setLargeCellSubset(cSet);
372 
373  vtkUnstructuredGrid* vtkmesh = volumeVTKMesh
374  (
375  subsetter.subMesh(),
376  csetPolyDecomp_[datasetNo]
377  );
378 
379  if (vtkmesh)
380  {
381  // superCells + addPointCellLabels must contain global cell ids
383  (
384  subsetter.cellMap(),
385  csetPolyDecomp_[datasetNo].superCells()
386  );
388  (
389  subsetter.cellMap(),
390  csetPolyDecomp_[datasetNo].addPointCellLabels()
391  );
392 
393  // copy pointMap as well, otherwise pointFields fail
394  csetPolyDecomp_[datasetNo].pointMap() = subsetter.pointMap();
395 
396  AddToBlock(output, vtkmesh, range, datasetNo, partName);
397  vtkmesh->Delete();
398 
399  partDataset_[partId] = datasetNo++;
400  }
401  }
402 
403  // anything added?
404  if (datasetNo)
405  {
406  ++blockNo;
407  }
408 
409  if (debug)
410  {
411  Info<< "<end> Foam::vtkPV3Foam::convertMeshCellSets" << endl;
412  printMemory();
413  }
414 }
415 
416 
417 void Foam::vtkPV3Foam::convertMeshFaceZones
418 (
419  vtkMultiBlockDataSet* output,
420  int& blockNo
421 )
422 {
423  arrayRange& range = arrayRangeFaceZones_;
424  range.block(blockNo); // set output block
425  label datasetNo = 0; // restart at dataset 0
426  const fvMesh& mesh = *meshPtr_;
427 
428  if (range.empty())
429  {
430  return;
431  }
432 
433  if (debug)
434  {
435  Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
436  printMemory();
437  }
438 
439  const faceZoneMesh& zMesh = mesh.faceZones();
440  for (int partId = range.start(); partId < range.end(); ++partId)
441  {
442  const word zoneName = getPartName(partId);
443  const label zoneId = zMesh.findZoneID(zoneName);
444 
445  if (!partStatus_[partId] || zoneId < 0)
446  {
447  continue;
448  }
449 
450  if (debug)
451  {
452  Info<< "Creating VTKmesh for faceZone[" << zoneId << "] "
453  << zoneName << endl;
454  }
455 
456  vtkPolyData* vtkmesh = patchVTKMesh(zoneName, zMesh[zoneId]());
457 
458  if (vtkmesh)
459  {
460  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
461  vtkmesh->Delete();
462 
463  partDataset_[partId] = datasetNo++;
464  }
465  }
466 
467  // anything added?
468  if (datasetNo)
469  {
470  ++blockNo;
471  }
472 
473  if (debug)
474  {
475  Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceZones" << endl;
476  printMemory();
477  }
478 }
479 
480 
481 void Foam::vtkPV3Foam::convertMeshFaceSets
482 (
483  vtkMultiBlockDataSet* output,
484  int& blockNo
485 )
486 {
487  arrayRange& range = arrayRangeFaceSets_;
488  range.block(blockNo); // set output block
489  label datasetNo = 0; // restart at dataset 0
490  const fvMesh& mesh = *meshPtr_;
491 
492  if (debug)
493  {
494  Info<< "<beg> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
495  printMemory();
496  }
497 
498  for (int partId = range.start(); partId < range.end(); ++partId)
499  {
500  const word partName = getPartName(partId);
501 
502  if (!partStatus_[partId])
503  {
504  continue;
505  }
506 
507  if (debug)
508  {
509  Info<< "Creating VTK mesh for faceSet=" << partName << endl;
510  }
511 
512  const faceSet fSet(mesh, partName);
513 
514  vtkPolyData* vtkmesh = faceSetVTKMesh(mesh, fSet);
515  if (vtkmesh)
516  {
517  AddToBlock(output, vtkmesh, range, datasetNo, partName);
518  vtkmesh->Delete();
519 
520  partDataset_[partId] = datasetNo++;
521  }
522  }
523 
524  // anything added?
525  if (datasetNo)
526  {
527  ++blockNo;
528  }
529 
530  if (debug)
531  {
532  Info<< "<end> Foam::vtkPV3Foam::convertMeshFaceSets" << endl;
533  printMemory();
534  }
535 }
536 
537 
538 void Foam::vtkPV3Foam::convertMeshPointZones
539 (
540  vtkMultiBlockDataSet* output,
541  int& blockNo
542 )
543 {
544  arrayRange& range = arrayRangePointZones_;
545  range.block(blockNo); // set output block
546  label datasetNo = 0; // restart at dataset 0
547  const fvMesh& mesh = *meshPtr_;
548 
549  if (debug)
550  {
551  Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
552  printMemory();
553  }
554 
555  if (range.size())
556  {
557  const pointZoneMesh& zMesh = mesh.pointZones();
558  for (int partId = range.start(); partId < range.end(); ++partId)
559  {
560  word zoneName = getPartName(partId);
561  label zoneId = zMesh.findZoneID(zoneName);
562 
563  if (!partStatus_[partId] || zoneId < 0)
564  {
565  continue;
566  }
567 
568  vtkPolyData* vtkmesh = pointZoneVTKMesh(mesh, zMesh[zoneId]);
569  if (vtkmesh)
570  {
571  AddToBlock(output, vtkmesh, range, datasetNo, zoneName);
572  vtkmesh->Delete();
573 
574  partDataset_[partId] = datasetNo++;
575  }
576  }
577  }
578 
579  // anything added?
580  if (datasetNo)
581  {
582  ++blockNo;
583  }
584 
585  if (debug)
586  {
587  Info<< "<end> Foam::vtkPV3Foam::convertMeshPointZones" << endl;
588  printMemory();
589  }
590 }
591 
592 
593 
594 void Foam::vtkPV3Foam::convertMeshPointSets
595 (
596  vtkMultiBlockDataSet* output,
597  int& blockNo
598 )
599 {
600  arrayRange& range = arrayRangePointSets_;
601  range.block(blockNo); // set output block
602  label datasetNo = 0; // restart at dataset 0
603  const fvMesh& mesh = *meshPtr_;
604 
605  if (debug)
606  {
607  Info<< "<beg> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
608  printMemory();
609  }
610 
611  for (int partId = range.start(); partId < range.end(); ++partId)
612  {
613  word partName = getPartName(partId);
614 
615  if (!partStatus_[partId])
616  {
617  continue;
618  }
619 
620  if (debug)
621  {
622  Info<< "Creating VTK mesh for pointSet=" << partName << endl;
623  }
624 
625  const pointSet pSet(mesh, partName);
626 
627  vtkPolyData* vtkmesh = pointSetVTKMesh(mesh, pSet);
628  if (vtkmesh)
629  {
630  AddToBlock(output, vtkmesh, range, datasetNo, partName);
631  vtkmesh->Delete();
632 
633  partDataset_[partId] = datasetNo++;
634  }
635  }
636 
637  // anything added?
638  if (datasetNo)
639  {
640  ++blockNo;
641  }
642 
643  if (debug)
644  {
645  Info<< "<end> Foam::vtkPV3Foam::convertMeshPointSets" << endl;
646  printMemory();
647  }
648 }
649 
650 
651 // ************************************************************************* //
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
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:253
static void printMemory()
Simple memory used debugging information.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
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:295
PrimitivePatch< face, UIndirectList, const pointField & > uindirectPrimitivePatch
Foam::uindirectPrimitivePatch.
messageStream Info
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.