vtkPVFoam.H
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-2017 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  //- OpenFOAM mesh
254  fvMesh* meshPtr_;
255 
256  //- The mesh region
257  word meshRegion_;
258 
259  //- The mesh directory for the region
260  fileName meshDir_;
261 
262  //- The time index
263  int timeIndex_;
264 
265  //- Track changes in mesh geometry
266  bool meshChanged_;
267 
268  //- Track changes in fields
269  bool fieldsChanged_;
270 
271  //- Selected geometrical parts (internalMesh, patches, ...)
272  boolList partStatus_;
273 
274  //- Datasets corresponding to selected geometrical pieces
275  // a negative number indicates that no vtkmesh exists for this piece
276  labelList partDataset_;
277 
278  //- First instance and size of various mesh parts
279  // used to index into partStatus_ and partDataset_
280  arrayRange arrayRangeVolume_;
281  arrayRange arrayRangePatches_;
282  arrayRange arrayRangeLagrangian_;
283  arrayRange arrayRangeCellZones_;
284  arrayRange arrayRangeFaceZones_;
285  arrayRange arrayRangePointZones_;
286  arrayRange arrayRangeCellSets_;
287  arrayRange arrayRangeFaceSets_;
288  arrayRange arrayRangePointSets_;
289 
290  //- Decomposed cells information (mesh regions)
291  // TODO: regions
292  List<polyDecomp> regionPolyDecomp_;
293 
294  //- Decomposed cells information (cellZone meshes)
295  List<polyDecomp> zonePolyDecomp_;
296 
297  //- Decomposed cells information (cellSet meshes)
298  List<polyDecomp> csetPolyDecomp_;
299 
300  //- List of patch names for rendering to window
301  List<vtkTextActor*> patchTextActorsPtrs_;
302 
303  // Private Member Functions
304 
305  // Convenience method use to convert the readers from VTK 5
306  // multiblock API to the current composite data infrastructure
307  static void AddToBlock
308  (
309  vtkMultiBlockDataSet* output,
310  vtkDataSet* dataset,
311  const arrayRange&,
312  const label datasetNo,
313  const std::string& datasetName
314  );
315 
316  // Convenience method use to convert the readers from VTK 5
317  // multiblock API to the current composite data infrastructure
318  static vtkDataSet* GetDataSetFromBlock
319  (
320  vtkMultiBlockDataSet* output,
321  const arrayRange&,
322  const label datasetNo
323  );
324 
325  // Convenience method use to convert the readers from VTK 5
326  // multiblock API to the current composite data infrastructure
327  static label GetNumberOfDataSets
328  (
329  vtkMultiBlockDataSet* output,
330  const arrayRange&
331  );
332 
333  //- Reset data counters
334  void resetCounters();
335 
336  // Update information helper functions
337 
338  //- Update the mesh parts selected in the GUI
339  void updateMeshPartsStatus();
340 
341  //- Internal mesh info
342  void updateInfoInternalMesh(vtkDataArraySelection*);
343 
344  //- Lagrangian info
345  void updateInfoLagrangian(vtkDataArraySelection*);
346 
347  //- Patch info
348  void updateInfoPatches(vtkDataArraySelection*, stringList&);
349 
350  //- Set info
351  void updateInfoSets(vtkDataArraySelection*);
352 
353  //- Zone info
354  void updateInfoZones(vtkDataArraySelection*);
355 
356  //- Get non-empty zone names for zoneType from file
357  wordList getZoneNames(const word& zoneType) const;
358 
359  //- Get non-empty zone names from mesh info
360  template<class ZoneType>
361  wordList getZoneNames
362  (
364  ) const;
365 
366  //- Add objects of Type to paraview array selection
367  template<class Type>
368  label addToSelection
369  (
370  vtkDataArraySelection*,
371  const IOobjectList&,
372  const string& suffix=string::null
373  );
374 
375  //- Field info
376  template<template<class> class patchType, class meshType>
377  void updateInfoFields(vtkDataArraySelection*);
378 
379  //- Lagrangian field info
380  void updateInfoLagrangianFields();
381 
382 
383  // Update helper functions
384 
385  //- OpenFOAM mesh
386  void updateFoamMesh();
387 
388  //- Reduce memory footprint after conversion
389  void reduceMemory();
390 
391  //- Volume fields
392  void updateVolFields(vtkMultiBlockDataSet*);
393 
394  //- Point fields
395  void updatePointFields(vtkMultiBlockDataSet*);
396 
397  //- Lagrangian fields
398  void updateLagrangianFields(vtkMultiBlockDataSet*);
399 
400 
401  // Mesh conversion functions
402 
403  //- Volume mesh
404  void convertMeshVolume(vtkMultiBlockDataSet*, int& blockNo);
405 
406  //- Lagrangian mesh
407  void convertMeshLagrangian(vtkMultiBlockDataSet*, int& blockNo);
408 
409  //- Patch meshes
410  void convertMeshPatches(vtkMultiBlockDataSet*, int& blockNo);
411 
412  //- Cell zone meshes
413  void convertMeshCellZones(vtkMultiBlockDataSet*, int& blockNo);
414 
415  //- Face zone meshes
416  void convertMeshFaceZones(vtkMultiBlockDataSet*, int& blockNo);
417 
418  //- Point zone meshes
419  void convertMeshPointZones(vtkMultiBlockDataSet*, int& blockNo);
420 
421  //- Cell set meshes
422  void convertMeshCellSets(vtkMultiBlockDataSet*, int& blockNo);
423 
424  //- Face set meshes
425  void convertMeshFaceSets(vtkMultiBlockDataSet*, int& blockNo);
426 
427  //- Point set meshes
428  void convertMeshPointSets(vtkMultiBlockDataSet*, int& blockNo);
429 
430 
431  // Add mesh functions
432 
433  //- Add internal mesh/cell set meshes
434  vtkUnstructuredGrid* volumeVTKMesh(const fvMesh&, polyDecomp&);
435 
436  //- Add Lagrangian mesh
437  vtkPolyData* lagrangianVTKMesh
438  (
439  const fvMesh&,
440  const word& cloudName
441  );
442 
443  //- Add patch mesh
444  template<class PatchType>
445  vtkPolyData* patchVTKMesh(const word& name, const PatchType&);
446 
447  //- Add point zone
448  vtkPolyData* pointZoneVTKMesh
449  (
450  const fvMesh&,
451  const labelList& pointLabels
452  );
453 
454  //- Add face set mesh
455  vtkPolyData* faceSetVTKMesh
456  (
457  const fvMesh&,
458  const faceSet&
459  );
460 
461  //- Add point mesh
462  vtkPolyData* pointSetVTKMesh
463  (
464  const fvMesh&,
465  const pointSet&
466  );
467 
468  // Field conversion functions
469 
470  //- Convert volume fields
471  void convertVolFields(vtkMultiBlockDataSet*);
472 
473  //- Convert point fields
474  void convertPointFields(vtkMultiBlockDataSet*);
475 
476  //- Convert Lagrangian fields
477  void convertLagrangianFields(vtkMultiBlockDataSet*);
478 
479 
480  //- Add the fields in the selected time directory to the selection
481  // lists
482  template<class GeoField>
483  label addObjectsToSelection
484  (
485  vtkDataArraySelection*,
486  const IOobjectList&,
487  const string& suffix=string::null
488  );
489 
490 
491  // Convert OpenFOAM fields
492 
493  //- Volume fields - all types
494  template<class Type>
495  void convertVolFields
496  (
497  const fvMesh&,
499  const IOobjectList&,
500  const bool interpFields,
501  vtkMultiBlockDataSet* output
502  );
503 
504  //- Volume field - all selected parts
505  template<class Type>
506  void convertVolFieldBlock
507  (
510  vtkMultiBlockDataSet* output,
511  const arrayRange&,
512  const List<polyDecomp>& decompLst
513  );
514 
515  //- Volume field
516  template<class Type>
517  void convertVolField
518  (
520  vtkMultiBlockDataSet* output,
521  const arrayRange&,
522  const label datasetNo,
523  const polyDecomp&
524  );
525 
526  //- Patch field
527  template<class Type>
528  void convertPatchField
529  (
530  const word& name,
531  const Field<Type>&,
532  vtkMultiBlockDataSet* output,
533  const arrayRange&,
534  const label datasetNo
535  );
536 
537  //- Face set/zone field
538  template<class Type>
539  void convertFaceField
540  (
542  vtkMultiBlockDataSet* output,
543  const arrayRange&,
544  const label datasetNo,
545  const fvMesh&,
546  const labelList& faceLabels
547  );
548 
549  //- Lagrangian fields - all types
550  template<class Type>
551  void convertLagrangianFields
552  (
553  const IOobjectList&,
554  vtkMultiBlockDataSet* output,
555  const label datasetNo
556  );
557 
558  //- Lagrangian field
559  template<class Type>
560  void convertLagrangianField
561  (
562  const IOField<Type>&,
563  vtkMultiBlockDataSet* output,
564  const arrayRange&,
565  const label datasetNo
566  );
567 
568  //- Point fields - all types
569  template<class Type>
570  void convertPointFields
571  (
572  const fvMesh&,
573  const pointMesh&,
574  const IOobjectList&,
575  vtkMultiBlockDataSet* output
576  );
577 
578  //- Point field - all selected parts
579  template<class Type>
580  void convertPointFieldBlock
581  (
583  vtkMultiBlockDataSet* output,
584  const arrayRange&,
585  const List<polyDecomp>&
586  );
587 
588  //- Point fields
589  template<class Type>
590  void convertPointField
591  (
594  vtkMultiBlockDataSet* output,
595  const arrayRange&,
596  const label datasetNo,
597  const polyDecomp&
598  );
599 
600  //- Patch point field
601  template<class Type>
602  void convertPatchPointField
603  (
604  const word& name,
605  const Field<Type>&,
606  vtkMultiBlockDataSet* output,
607  const arrayRange&,
608  const label datasetNo
609  );
610 
611 
612  // GUI selection helper functions
613 
614  //- Only keep what is listed in hashSet
615  static void pruneObjectList
616  (
617  IOobjectList&,
618  const wordHashSet&
619  );
620 
621  //- Retrieve the current selections
622  static wordHashSet getSelected(vtkDataArraySelection*);
623 
624  //- Retrieve a sub-list of the current selections
625  static wordHashSet getSelected
626  (
627  vtkDataArraySelection*,
628  const arrayRange&
629  );
630 
631  //- Retrieve the current selections
632  static stringList getSelectedArrayEntries(vtkDataArraySelection*);
633 
634  //- Retrieve a sub-list of the current selections
635  static stringList getSelectedArrayEntries
636  (
637  vtkDataArraySelection*,
638  const arrayRange&
639  );
640 
641  //- Set selection(s)
642  static void setSelectedArrayEntries
643  (
644  vtkDataArraySelection*,
645  const stringList&
646  );
647 
648  //- Get the first word from the mesh parts selection
649  word getPartName(const int);
650 
651 
652  //- Disallow default bitwise copy construct
653  vtkPVFoam(const vtkPVFoam&);
654 
655  //- Disallow default bitwise assignment
656  void operator=(const vtkPVFoam&);
657 
658 
659 public:
660 
661  //- Static data members
662 
663  ClassName("vtkPVFoam");
664 
665 
666  // Constructors
667 
668  //- Construct from components
669  vtkPVFoam
670  (
671  const char* const FileName,
672  vtkPVFoamReader* reader
673  );
674 
675 
676  //- Destructor
677  ~vtkPVFoam();
678 
679 
680  // Member Functions
681 
682  //- Update
683  void updateInfo();
684 
685  void Update
686  (
687  vtkMultiBlockDataSet* output,
688  vtkMultiBlockDataSet* lagrangianOutput
689  );
690 
691  //- Clean any storage
692  void CleanUp();
693 
694  //- Allocate and return a list of selected times
695  // returns the count via the parameter
696  double* findTimes(int& nTimeSteps);
697 
698  //- Add/remove patch names to/from the view
699  void renderPatchNames(vtkRenderer*, const bool show);
700 
701  //- Set the runTime to the first plausible request time,
702  // returns the timeIndex
703  // sets to "constant" on error
704  int setTime(int count, const double requestTimes[]);
705 
706 
707  //- The current time index
708  int timeIndex() const
709  {
710  return timeIndex_;
711  }
712 
713 
714  // Access
715 
716  //- Debug information
717  void PrintSelf(ostream&, vtkIndent) const;
718 
719  //- Simple memory used debugging information
720  static void printMemory();
721 
722 };
723 
724 
725 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
726 
727 } // End namespace Foam
728 
729 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
730 
731 #ifdef NoRepository
732  #include "vtkPVFoamTemplates.C"
733 #endif
734 
735 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
736 
737 #endif
738 
739 // ************************************************************************* //
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:707
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
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.
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.