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-2023 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  wordList getZoneNames(const word& zoneType) const;
364 
365  //- Get non-empty zone names from mesh info
366  template<class ZoneType>
367  wordList getZoneNames
368  (
370  ) const;
371 
372  //- Add objects of Type to paraview array selection
373  template<class Type>
374  label addToSelection
375  (
376  vtkDataArraySelection*,
377  const IOobjectList&,
378  const string& suffix=string::null
379  );
380 
381  //- Add objects of Type to paraview array selection
382  template<template<class> class patchType, class meshType>
383  void addFieldsToSelection
384  (
385  vtkDataArraySelection *select,
386  const IOobjectList& objects,
387  const string& suffix=string::null
388  );
389 
390  //- Add objects of Type to paraview array selection
391  template<template<class> class patchType, class meshType>
392  void addInternalFieldsToSelection
393  (
394  vtkDataArraySelection *select,
395  const IOobjectList& objects,
396  const string& suffix=string::null
397  );
398 
399  //- Field info
400  void updateInfoFields();
401 
402  //- Lagrangian field info
403  void updateInfoLagrangianFields();
404 
405 
406  // Update helper functions
407 
408  //- OpenFOAM mesh
409  void updateFoamMesh();
410 
411  //- Reduce memory footprint after conversion
412  void reduceMemory();
413 
414  //- Volume fields
415  void updateVolFields(vtkMultiBlockDataSet*);
416 
417  //- Point fields
418  void updatePointFields(vtkMultiBlockDataSet*);
419 
420  //- Lagrangian fields
421  void updateLagrangianFields(vtkMultiBlockDataSet*);
422 
423 
424  // Mesh conversion functions
425 
426  //- Volume mesh
427  void convertMeshVolume(vtkMultiBlockDataSet*, int& blockNo);
428 
429  //- Lagrangian mesh
430  void convertMeshLagrangian(vtkMultiBlockDataSet*, int& blockNo);
431 
432  //- Patches
433  void convertMeshPatches(vtkMultiBlockDataSet*, int& blockNo);
434 
435  //- Cell zones
436  void convertMeshCellZones(vtkMultiBlockDataSet*, int& blockNo);
437 
438  //- Face zones
439  void convertMeshFaceZones(vtkMultiBlockDataSet*, int& blockNo);
440 
441  //- Point zones
442  void convertMeshPointZones(vtkMultiBlockDataSet*, int& blockNo);
443 
444  //- Cell sets
445  void convertMeshCellSets(vtkMultiBlockDataSet*, int& blockNo);
446 
447  //- Face sets
448  void convertMeshFaceSets(vtkMultiBlockDataSet*, int& blockNo);
449 
450  //- Point sets
451  void convertMeshPointSets(vtkMultiBlockDataSet*, int& blockNo);
452 
453 
454  // Add mesh functions
455 
456  //- Add internal mesh/cell set
457  vtkUnstructuredGrid* volumeVTKMesh(const fvMesh&, polyDecomp&);
458 
459  //- Add Lagrangian mesh
460  vtkPolyData* lagrangianVTKMesh
461  (
462  const fvMesh&,
463  const word& cloudName
464  );
465 
466  //- Add patch mesh
467  template<class PatchType>
468  vtkPolyData* patchVTKMesh(const word& name, const PatchType&);
469 
470  //- Add point zone
471  vtkPolyData* pointZoneVTKMesh
472  (
473  const fvMesh&,
474  const labelList& pointLabels
475  );
476 
477  //- Add face set
478  vtkPolyData* faceSetVTKMesh
479  (
480  const fvMesh&,
481  const faceSet&
482  );
483 
484  //- Add point mesh
485  vtkPolyData* pointSetVTKMesh
486  (
487  const fvMesh&,
488  const pointSet&
489  );
490 
491  // Field conversion functions
492 
493  //- Convert fields
494  void convertFields(vtkMultiBlockDataSet*);
495 
496  //- Convert Lagrangian fields
497  void convertLagrangianFields(vtkMultiBlockDataSet*);
498 
499 
500  //- Add the fields in the selected time directory to the selection
501  // lists
502  template<class GeoField>
503  label addObjectsToSelection
504  (
505  vtkDataArraySelection*,
506  const IOobjectList&,
507  const string& suffix=string::null
508  );
509 
510 
511  // Convert OpenFOAM fields
512 
513  //- Volume fields - all types
514  template<class Type>
515  void convertVolFields
516  (
517  const fvMesh&,
519  const IOobjectList&,
520  const bool interpFields,
521  vtkMultiBlockDataSet* output
522  );
523 
524  //- Volume internal fields - all types
525  template<class Type>
526  void convertVolInternalFields
527  (
528  const fvMesh&,
529  const IOobjectList&,
530  vtkMultiBlockDataSet* output
531  );
532 
533  //- Volume field - all selected parts
534  template<class Type>
535  void convertVolFieldBlock
536  (
537  const VolField<Type>&,
539  vtkMultiBlockDataSet* output,
540  const arrayRange&,
541  const List<polyDecomp>& decompLst
542  );
543 
544  //- Volume internal field - all selected parts
545  template<class Type>
546  void convertVolInternalFieldBlock
547  (
548  const typename VolField<Type>
549  ::Internal&,
550  vtkMultiBlockDataSet* output,
551  const arrayRange&,
552  const List<polyDecomp>& decompLst
553  );
554 
555  //- Volume field
556  template<class Type>
557  void convertVolInternalField
558  (
559  const typename VolField<Type>
560  ::Internal&,
561  vtkMultiBlockDataSet* output,
562  const arrayRange&,
563  const label datasetNo,
564  const polyDecomp&
565  );
566 
567  //- Patch field
568  template<class Type>
569  void convertPatchField
570  (
571  const word& name,
572  const Field<Type>&,
573  vtkMultiBlockDataSet* output,
574  const arrayRange&,
575  const label datasetNo
576  );
577 
578  //- Face set/zone field from volField
579  template<class Type>
580  void convertSurfaceField
581  (
582  const VolField<Type>&,
583  vtkMultiBlockDataSet* output,
584  const arrayRange&,
585  const label datasetNo,
586  const fvMesh&,
587  const labelList& faceLabels
588  );
589 
590  //- Face set/zone field from surfaceField
591  template<class Type>
592  void convertSurfaceField
593  (
594  const SurfaceField<Type>& tf,
595  vtkMultiBlockDataSet* output,
596  const arrayRange& range,
597  const label datasetNo,
598  const fvMesh& mesh,
599  const labelList& faceLabels
600  );
601 
602  //- Surface fields - all types
603  template<class Type>
604  void convertSurfaceFields
605  (
606  const fvMesh& mesh,
607  const IOobjectList& objects,
608  vtkMultiBlockDataSet* output
609  );
610 
611  //- Lagrangian fields - all types
612  template<class Type>
613  void convertLagrangianFields
614  (
615  const IOobjectList&,
616  vtkMultiBlockDataSet* output,
617  const label datasetNo
618  );
619 
620  //- Lagrangian field
621  template<class Type>
622  void convertLagrangianField
623  (
624  const IOField<Type>&,
625  vtkMultiBlockDataSet* output,
626  const arrayRange&,
627  const label datasetNo
628  );
629 
630  //- Point fields - all types
631  template<class Type>
632  void convertPointFields
633  (
634  const fvMesh&,
635  const pointMesh&,
636  const IOobjectList&,
637  vtkMultiBlockDataSet* output
638  );
639 
640  //- Point field - all selected parts
641  template<class Type>
642  void convertPointFieldBlock
643  (
644  const PointField<Type>&,
645  vtkMultiBlockDataSet* output,
646  const arrayRange&,
647  const List<polyDecomp>&
648  );
649 
650  //- Point fields
651  template<class Type>
652  void convertPointField
653  (
654  const PointField<Type>&,
655  const VolField<Type>&,
656  vtkMultiBlockDataSet* output,
657  const arrayRange&,
658  const label datasetNo,
659  const polyDecomp&
660  );
661 
662  //- Patch point field
663  template<class Type>
664  void convertPatchPointField
665  (
666  const word& name,
667  const Field<Type>&,
668  vtkMultiBlockDataSet* output,
669  const arrayRange&,
670  const label datasetNo
671  );
672 
673 
674  // GUI selection helper functions
675 
676  //- Only keep what is listed in hashSet
677  static void pruneObjectList
678  (
679  IOobjectList&,
680  const wordHashSet&
681  );
682 
683  //- Get the list of selected objects
684  IOobjectList getObjects
685  (
686  const wordHashSet& selected,
687  const fvMesh& mesh,
688  const fileName& instance,
689  const fileName& local = ""
690  );
691 
692  //- Retrieve the current selections
693  static wordHashSet getSelected(vtkDataArraySelection*);
694 
695  //- Retrieve a sub-list of the current selections
696  static wordHashSet getSelected
697  (
698  vtkDataArraySelection*,
699  const arrayRange&
700  );
701 
702  //- Retrieve the current selections
703  static stringList getSelectedArrayEntries(vtkDataArraySelection*);
704 
705  //- Retrieve a sub-list of the current selections
706  static stringList getSelectedArrayEntries
707  (
708  vtkDataArraySelection*,
709  const arrayRange&
710  );
711 
712  //- Set selection(s)
713  static void setSelectedArrayEntries
714  (
715  vtkDataArraySelection*,
716  const stringList&
717  );
718 
719  //- Get the first word from the mesh parts selection
720  word getPartName(const int);
721 
722 
723 public:
724 
725  //- Static data members
726 
727  ClassName("vtkPVFoam");
728 
729 
730  // Constructors
731 
732  //- Construct from components
733  vtkPVFoam
734  (
735  const char* const FileName,
736  vtkPVFoamReader* reader
737  );
738 
739  //- Disallow default bitwise copy construction
740  vtkPVFoam(const vtkPVFoam&) = delete;
741 
742 
743  //- Destructor
744  ~vtkPVFoam();
745 
746 
747  // Member Functions
748 
749  //- Update information for selection dialogs
750  void updateInfo();
751 
752  //- Update
753  void Update
754  (
755  vtkMultiBlockDataSet* output,
756  vtkMultiBlockDataSet* lagrangianOutput
757  );
758 
759  //- Clean any storage
760  void CleanUp();
761 
762  //- Allocate and return a list of selected times
763  // returns the count via the parameter
764  double* findTimes(int& nTimeSteps);
765 
766  //- Add/remove patch names to/from the view
767  void renderPatchNames(vtkRenderer*, const bool show);
768 
769  //- Set the runTime to the first plausible request time,
770  // returns the timeIndex
771  // sets to "constant" on error
772  int setTime(int count, const double requestTimes[]);
773 
774  //- The current time index
775  int timeIndex() const
776  {
777  return timeIndex_;
778  }
779 
780 
781  // Access
782 
783  //- Debug information
784  void PrintSelf(ostream&, vtkIndent) const;
785 
786  //- Simple memory used debugging information
787  static void printMemory();
788 
789 
790  // Member Operators
791 
792 
793  //- Disallow default bitwise assignment
794  void operator=(const vtkPVFoam&) = delete;
795 };
796 
797 
798 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
799 
800 } // End namespace Foam
801 
802 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
803 
804 #ifdef NoRepository
805  #include "vtkPVFoamTemplates.C"
806 #endif
807 
808 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
809 
810 #endif
811 
812 // ************************************************************************* //
scalar range
label n
Pre-declare SubField and related Field type.
Definition: Field.H:82
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
A list of mesh zones.
Definition: MeshZones.H:66
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:160
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:101
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
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:774
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"))