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