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