vtkPVFoam.H
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 Class
25  Foam::vtkPVFoam
26 
27 Description
28  Provides a reader interface for OpenFOAM to VTK interaction.
29 
30 SourceFiles
31  vtkPVFoam.C
32  vtkPVFoam.H
33  vtkPVFoamFields.C
34  vtkPVFoamMesh.C
35  vtkPVFoamMeshLagrangian.C
36  vtkPVFoamTemplates.C
37  vtkPVFoamMeshSet.C
38  vtkPVFoamMeshVolume.C
39  vtkPVFoamMeshZone.C
40  vtkPVFoamFaceField.H
41  vtkPVFoamLagrangianFields.H
42  vtkPVFoamPatchField.H
43  vtkPVFoamPointFields.H
44  vtkPVFoamPoints.H
45  vtkPVFoamUpdateInfo.C
46  vtkPVFoamUpdateInfoFields.H
47  vtkPVFoamUtils.C
48  vtkPVFoamVolFields.H
49  vtkPVFoamAddToSelection.H
50 
51  // Needed by VTK:
52  vtkDataArrayTemplateImplicit.txx
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef vtkPVFoam_H
57 #define vtkPVFoam_H
58 
59 // do not include legacy strstream headers
60 #ifndef VTK_EXCLUDE_STRSTREAM_HEADERS
61 # define VTK_EXCLUDE_STRSTREAM_HEADERS
62 #endif
63 
64 #include "className.H"
65 #include "fileName.H"
66 #include "stringList.H"
67 #include "wordList.H"
68 #include "primitivePatch.H"
70 #include "volPointInterpolation.H"
71 
72 // * * * * * * * * * * * * * Forward Declarations * * * * * * * * * * * * * //
73 
74 class vtkDataArraySelection;
75 class vtkDataSet;
76 class vtkPoints;
77 class vtkPVFoamReader;
78 class vtkRenderer;
79 class vtkTextActor;
80 class vtkMultiBlockDataSet;
81 class vtkPolyData;
82 class vtkUnstructuredGrid;
83 class vtkIndent;
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 namespace Foam
88 {
89 
90 // OpenFOAM class forward declarations
91 class argList;
92 class Time;
93 class fvMesh;
94 class IOobjectList;
95 class polyPatch;
96 class faceSet;
97 class pointSet;
98 
99 template<class Type> class IOField;
100 template<class Type> class List;
101 
102 /*---------------------------------------------------------------------------*\
103  Class vtkPVFoam Declaration
104 \*---------------------------------------------------------------------------*/
106 class vtkPVFoam
107 {
108  // Private classes
109 
110  //- Bookkeeping for GUI checklists and the multi-block organization
111  class arrayRange
112  {
113  const char *name_;
114  int block_;
115  int start_;
116  int size_;
117 
118  public:
119 
120  arrayRange(const char *name, const int blockNo=0)
121  :
122  name_(name),
123  block_(blockNo),
124  start_(0),
125  size_(0)
126  {}
127 
128  //- Return the block holding these datasets
129  int block() const
130  {
131  return block_;
132  }
133 
134  //- Assign block number, return previous value
135  int block(int blockNo)
136  {
137  int prev = block_;
138  block_ = blockNo;
139  return prev;
140  }
141 
142  //- Return block name
143  const char* name() const
144  {
145  return name_;
146  }
147 
148  //- Return array start index
149  int start() const
150  {
151  return start_;
152  }
153 
154  //- Return array end index
155  int end() const
156  {
157  return start_ + size_;
158  }
159 
160  //- Return sublist size
161  int size() const
162  {
163  return size_;
164  }
165 
166  bool empty() const
167  {
168  return !size_;
169  }
170 
171  //- Reset the size to zero and optionally assign a new start
172  void reset(const int startAt = 0)
173  {
174  start_ = startAt;
175  size_ = 0;
176  }
177 
178  //- Increment the size
179  void operator+=(const int n)
180  {
181  size_ += n;
182  }
183  };
184 
185  //- Bookkeeping for polyhedral cell decomposition
186  // hide in extra pointMap (cellSet/cellZone) for now
187  class polyDecomp
188  {
189  labelList superCells_;
190  labelList addPointCellLabels_;
191  labelList pointMap_;
192 
193  public:
194 
195  polyDecomp()
196  {}
197 
198  //- Label of original cell for decomposed cells
199  labelList& superCells()
200  {
201  return superCells_;
202  }
203 
204  //- Label of original cell for decomposed cells
205  const labelList& superCells() const
206  {
207  return superCells_;
208  }
209 
210  //- Cell-centre labels for additional points of decomposed cells
211  labelList& addPointCellLabels()
212  {
213  return addPointCellLabels_;
214  }
215 
216  //- Cell-centre labels for additional points of decomposed cells
217  const labelList& addPointCellLabels() const
218  {
219  return addPointCellLabels_;
220  }
221 
222  //- Point labels for subsetted meshes
223  labelList& pointMap()
224  {
225  return pointMap_;
226  }
227 
228  //- Point labels for subsetted meshes
229  const labelList& pointMap() const
230  {
231  return pointMap_;
232  }
233 
234 
235  //- Clear
236  void clear()
237  {
238  superCells_.clear();
239  addPointCellLabels_.clear();
240  pointMap_.clear();
241  }
242  };
243 
244 
245  // Private Data
246 
247  //- Access to the controlling vtkPVFoamReader
248  vtkPVFoamReader* reader_;
249 
250  //- OpenFOAM time control
251  autoPtr<Time> dbPtr_;
252 
253  //- Configuration dictionary
254  dictionary configDict_;
255 
256  //- OpenFOAM mesh
257  fvMesh* meshPtr_;
258 
259  //- The mesh region
260  word meshRegion_;
261 
262  //- The mesh directory for the region
263  fileName meshDir_;
264 
265  //- The time index
266  int timeIndex_;
267 
268  //- Track changes in mesh geometry
269  bool meshChanged_;
270 
271  //- Track changes in fields
272  bool fieldsChanged_;
273 
274  //- Selected geometrical parts (internalMesh, patches, ...)
275  boolList partStatus_;
276 
277  //- Datasets corresponding to selected geometrical pieces
278  // a negative number indicates that no vtkmesh exists for this piece
279  labelList partDataset_;
280 
281  //- First instance and size of various mesh parts
282  // used to index into partStatus_ and partDataset_
283  arrayRange arrayRangeVolume_;
284  arrayRange arrayRangePatches_;
285  arrayRange arrayRangeLagrangian_;
286  arrayRange arrayRangeCellZones_;
287  arrayRange arrayRangeFaceZones_;
288  arrayRange arrayRangePointZones_;
289  arrayRange arrayRangeCellSets_;
290  arrayRange arrayRangeFaceSets_;
291  arrayRange arrayRangePointSets_;
292 
293  //- Decomposed cells information (mesh regions)
294  // TODO: regions
295  List<polyDecomp> regionPolyDecomp_;
296 
297  //- Decomposed cells information (cellZone meshes)
298  List<polyDecomp> zonePolyDecomp_;
299 
300  //- Decomposed cells information (cellSet meshes)
301  List<polyDecomp> csetPolyDecomp_;
302 
303  //- List of patch names for rendering to window
304  List<vtkTextActor*> patchTextActorsPtrs_;
305 
306  // Private Member Functions
307 
308  // Convenience method use to convert the readers from VTK 5
309  // multiblock API to the current composite data infrastructure
310  static void AddToBlock
311  (
312  vtkMultiBlockDataSet* output,
313  vtkDataSet* dataset,
314  const arrayRange&,
315  const label datasetNo,
316  const std::string& datasetName
317  );
318 
319  // Convenience method use to convert the readers from VTK 5
320  // multiblock API to the current composite data infrastructure
321  static vtkDataSet* GetDataSetFromBlock
322  (
323  vtkMultiBlockDataSet* output,
324  const arrayRange&,
325  const label datasetNo
326  );
327 
328  // Convenience method use to convert the readers from VTK 5
329  // multiblock API to the current composite data infrastructure
330  static label GetNumberOfDataSets
331  (
332  vtkMultiBlockDataSet* output,
333  const arrayRange&
334  );
335 
336  //- Reset data counters
337  void resetCounters();
338 
339  // Update information helper functions
340 
341  //- Update the mesh parts selected in the GUI
342  void updateMeshPartsStatus();
343 
344  //- Internal mesh info
345  void updateInfoInternalMesh(vtkDataArraySelection*);
346 
347  //- Lagrangian info
348  void updateInfoLagrangian(vtkDataArraySelection*);
349 
350  //- Patch info
351  void updateInfoPatches
352  (
353  vtkDataArraySelection*,
354  stringList&,
355  const bool
356  );
357 
358  //- Set info
359  void updateInfoSets(vtkDataArraySelection*);
360 
361  //- Zone info
362  void updateInfoZones(vtkDataArraySelection*);
363 
364  //- Get non-empty zone names for zoneType from file
365  wordList getZoneNames(const word& zoneType) const;
366 
367  //- Get non-empty zone names from mesh info
368  template<class ZoneType>
369  wordList getZoneNames
370  (
372  ) const;
373 
374  //- Add objects of Type to paraview array selection
375  template<class Type>
376  label addToSelection
377  (
378  vtkDataArraySelection*,
379  const IOobjectList&,
380  const string& suffix=string::null
381  );
382 
383  //- Field info
384  template<template<class> class patchType, class meshType>
385  void updateInfoFields(vtkDataArraySelection*);
386 
387  //- Lagrangian field info
388  void updateInfoLagrangianFields();
389 
390 
391  // Update helper functions
392 
393  //- OpenFOAM mesh
394  void updateFoamMesh();
395 
396  //- Reduce memory footprint after conversion
397  void reduceMemory();
398 
399  //- Volume fields
400  void updateVolFields(vtkMultiBlockDataSet*);
401 
402  //- Point fields
403  void updatePointFields(vtkMultiBlockDataSet*);
404 
405  //- Lagrangian fields
406  void updateLagrangianFields(vtkMultiBlockDataSet*);
407 
408 
409  // Mesh conversion functions
410 
411  //- Volume mesh
412  void convertMeshVolume(vtkMultiBlockDataSet*, int& blockNo);
413 
414  //- Lagrangian mesh
415  void convertMeshLagrangian(vtkMultiBlockDataSet*, int& blockNo);
416 
417  //- Patch meshes
418  void convertMeshPatches(vtkMultiBlockDataSet*, int& blockNo);
419 
420  //- Cell zone meshes
421  void convertMeshCellZones(vtkMultiBlockDataSet*, int& blockNo);
422 
423  //- Face zone meshes
424  void convertMeshFaceZones(vtkMultiBlockDataSet*, int& blockNo);
425 
426  //- Point zone meshes
427  void convertMeshPointZones(vtkMultiBlockDataSet*, int& blockNo);
428 
429  //- Cell set meshes
430  void convertMeshCellSets(vtkMultiBlockDataSet*, int& blockNo);
431 
432  //- Face set meshes
433  void convertMeshFaceSets(vtkMultiBlockDataSet*, int& blockNo);
434 
435  //- Point set meshes
436  void convertMeshPointSets(vtkMultiBlockDataSet*, int& blockNo);
437 
438 
439  // Add mesh functions
440 
441  //- Add internal mesh/cell set meshes
442  vtkUnstructuredGrid* volumeVTKMesh(const fvMesh&, polyDecomp&);
443 
444  //- Add Lagrangian mesh
445  vtkPolyData* lagrangianVTKMesh
446  (
447  const fvMesh&,
448  const word& cloudName
449  );
450 
451  //- Add patch mesh
452  template<class PatchType>
453  vtkPolyData* patchVTKMesh(const word& name, const PatchType&);
454 
455  //- Add point zone
456  vtkPolyData* pointZoneVTKMesh
457  (
458  const fvMesh&,
459  const labelList& pointLabels
460  );
461 
462  //- Add face set mesh
463  vtkPolyData* faceSetVTKMesh
464  (
465  const fvMesh&,
466  const faceSet&
467  );
468 
469  //- Add point mesh
470  vtkPolyData* pointSetVTKMesh
471  (
472  const fvMesh&,
473  const pointSet&
474  );
475 
476  // Field conversion functions
477 
478  //- Convert volume fields
479  void convertVolFields(vtkMultiBlockDataSet*);
480 
481  //- Convert point fields
482  void convertPointFields(vtkMultiBlockDataSet*);
483 
484  //- Convert Lagrangian fields
485  void convertLagrangianFields(vtkMultiBlockDataSet*);
486 
487 
488  //- Add the fields in the selected time directory to the selection
489  // lists
490  template<class GeoField>
491  label addObjectsToSelection
492  (
493  vtkDataArraySelection*,
494  const IOobjectList&,
495  const string& suffix=string::null
496  );
497 
498 
499  // Convert OpenFOAM fields
500 
501  //- Volume fields - all types
502  template<class Type>
503  void convertVolFields
504  (
505  const fvMesh&,
507  const IOobjectList&,
508  const bool interpFields,
509  vtkMultiBlockDataSet* output
510  );
511 
512  //- Volume field - all selected parts
513  template<class Type>
514  void convertVolFieldBlock
515  (
518  vtkMultiBlockDataSet* output,
519  const arrayRange&,
520  const List<polyDecomp>& decompLst
521  );
522 
523  //- Volume field
524  template<class Type>
525  void convertVolField
526  (
528  vtkMultiBlockDataSet* output,
529  const arrayRange&,
530  const label datasetNo,
531  const polyDecomp&
532  );
533 
534  //- Patch field
535  template<class Type>
536  void convertPatchField
537  (
538  const word& name,
539  const Field<Type>&,
540  vtkMultiBlockDataSet* output,
541  const arrayRange&,
542  const label datasetNo
543  );
544 
545  //- Face set/zone field
546  template<class Type>
547  void convertFaceField
548  (
550  vtkMultiBlockDataSet* output,
551  const arrayRange&,
552  const label datasetNo,
553  const fvMesh&,
554  const labelList& faceLabels
555  );
556 
557  //- Lagrangian fields - all types
558  template<class Type>
559  void convertLagrangianFields
560  (
561  const IOobjectList&,
562  vtkMultiBlockDataSet* output,
563  const label datasetNo
564  );
565 
566  //- Lagrangian field
567  template<class Type>
568  void convertLagrangianField
569  (
570  const IOField<Type>&,
571  vtkMultiBlockDataSet* output,
572  const arrayRange&,
573  const label datasetNo
574  );
575 
576  //- Point fields - all types
577  template<class Type>
578  void convertPointFields
579  (
580  const fvMesh&,
581  const pointMesh&,
582  const IOobjectList&,
583  vtkMultiBlockDataSet* output
584  );
585 
586  //- Point field - all selected parts
587  template<class Type>
588  void convertPointFieldBlock
589  (
591  vtkMultiBlockDataSet* output,
592  const arrayRange&,
593  const List<polyDecomp>&
594  );
595 
596  //- Point fields
597  template<class Type>
598  void convertPointField
599  (
602  vtkMultiBlockDataSet* output,
603  const arrayRange&,
604  const label datasetNo,
605  const polyDecomp&
606  );
607 
608  //- Patch point field
609  template<class Type>
610  void convertPatchPointField
611  (
612  const word& name,
613  const Field<Type>&,
614  vtkMultiBlockDataSet* output,
615  const arrayRange&,
616  const label datasetNo
617  );
618 
619 
620  // GUI selection helper functions
621 
622  //- Only keep what is listed in hashSet
623  static void pruneObjectList
624  (
625  IOobjectList&,
626  const wordHashSet&
627  );
628 
629  //- Get the list of selected objects
630  IOobjectList getObjects
631  (
632  const wordHashSet& selected,
633  const fvMesh& mesh,
634  const fileName& instance,
635  const fileName& local = ""
636  );
637 
638  //- Retrieve the current selections
639  static wordHashSet getSelected(vtkDataArraySelection*);
640 
641  //- Retrieve a sub-list of the current selections
642  static wordHashSet getSelected
643  (
644  vtkDataArraySelection*,
645  const arrayRange&
646  );
647 
648  //- Retrieve the current selections
649  static stringList getSelectedArrayEntries(vtkDataArraySelection*);
650 
651  //- Retrieve a sub-list of the current selections
652  static stringList getSelectedArrayEntries
653  (
654  vtkDataArraySelection*,
655  const arrayRange&
656  );
657 
658  //- Set selection(s)
659  static void setSelectedArrayEntries
660  (
661  vtkDataArraySelection*,
662  const stringList&
663  );
664 
665  //- Get the first word from the mesh parts selection
666  word getPartName(const int);
667 
668 
669  //- Disallow default bitwise copy construct
670  vtkPVFoam(const vtkPVFoam&);
671 
672  //- Disallow default bitwise assignment
673  void operator=(const vtkPVFoam&);
674 
675 
676 public:
677 
678  //- Static data members
679 
680  ClassName("vtkPVFoam");
681 
682 
683  // Constructors
684 
685  //- Construct from components
686  vtkPVFoam
687  (
688  const char* const FileName,
689  vtkPVFoamReader* reader
690  );
691 
692 
693  //- Destructor
694  ~vtkPVFoam();
695 
696 
697  // Member Functions
698 
699  //- Update
700  void updateInfo();
701 
702  void Update
703  (
704  vtkMultiBlockDataSet* output,
705  vtkMultiBlockDataSet* lagrangianOutput
706  );
707 
708  //- Clean any storage
709  void CleanUp();
710 
711  //- Allocate and return a list of selected times
712  // returns the count via the parameter
713  double* findTimes(int& nTimeSteps);
714 
715  //- Add/remove patch names to/from the view
716  void renderPatchNames(vtkRenderer*, const bool show);
717 
718  //- Set the runTime to the first plausible request time,
719  // returns the timeIndex
720  // sets to "constant" on error
721  int setTime(int count, const double requestTimes[]);
722 
723 
724  //- The current time index
725  int timeIndex() const
726  {
727  return timeIndex_;
728  }
729 
730 
731  // Access
732 
733  //- Debug information
734  void PrintSelf(ostream&, vtkIndent) const;
735 
736  //- Simple memory used debugging information
737  static void printMemory();
738 
739 };
740 
741 
742 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
743 
744 } // End namespace Foam
745 
746 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
747 
748 #ifdef NoRepository
749  #include "vtkPVFoamTemplates.C"
750 #endif
751 
752 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
753 
754 #endif
755 
756 // ************************************************************************* //
A HashTable with keys but without contents.
Definition: HashSet.H:59
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
int timeIndex() const
The current time index.
Definition: vtkPVFoam.H:724
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.
A class for handling file names.
Definition: fileName.H:69
A list of face labels.
Definition: faceSet.H:48
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
A set of point labels.
Definition: pointSet.H:48
labelList pointLabels(nPoints, -1)
Provides a reader interface for OpenFOAM to VTK interaction.
Definition: vtkPVFoam.H:105
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
tUEqn clear()
double * findTimes(int &nTimeSteps)
Allocate and return a list of selected times.
void updateInfo()
Update.
void renderPatchNames(vtkRenderer *, const bool show)
Add/remove patch names to/from the view.
A list of mesh zones.
Generic GeometricField class.
void PrintSelf(ostream &, vtkIndent) const
Debug information.
void Update(vtkMultiBlockDataSet *output, vtkMultiBlockDataSet *lagrangianOutput)
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
ClassName("vtkPVFoam")
Static data members.
void CleanUp()
Clean any storage.
dynamicFvMesh & mesh
Pre-declare SubField and related Field type.
Definition: Field.H:57
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:124
A class for handling words, derived from string.
Definition: word.H:59
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:63
static const string null
An empty string.
Definition: string.H:86
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
int setTime(int count, const double requestTimes[])
Set the runTime to the first plausible request time,.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label n
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const word cloudName(propsDict.lookup("cloudName"))
~vtkPVFoam()
Destructor.
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Namespace for OpenFOAM.