BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3Foam.H
blob7efb47990bd73b360d9c54b949b2963896e590a9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::vtkPV3Foam
28 Description
29     Provides a reader interface for OpenFOAM to VTK interaction.
31 SourceFiles
32     vtkPV3Foam.C
33     vtkPV3Foam.H
34     vtkPV3FoamI.H
35     vtkPV3FoamFields.C
36     vtkPV3FoamMesh.C
37     vtkPV3FoamMeshLagrangian.C
38     vtkPV3FoamMeshPatch.C
39     vtkPV3FoamMeshSet.C
40     vtkPV3FoamMeshVolume.C
41     vtkPV3FoamMeshZone.C
42     vtkPV3FoamFaceField.H
43     vtkPV3FoamLagrangianFields.H
44     vtkPV3FoamPatchField.H
45     vtkPV3FoamPointFields.H
46     vtkPV3FoamPoints.H
47     vtkPV3FoamUpdateInfo.C
48     vtkPV3FoamUpdateInfoFields.H
49     vtkPV3FoamUtilities.C
50     vtkPV3FoamVolFields.H
51     vtkPV3FoamAddToSelection.H
53     // Needed by VTK:
54     vtkDataArrayTemplateImplicit.txx
56 \*---------------------------------------------------------------------------*/
58 #ifndef vtkPV3Foam_H
59 #define vtkPV3Foam_H
61 // do not include legacy strstream headers
62 #ifndef  VTK_EXCLUDE_STRSTREAM_HEADERS
63 # define VTK_EXCLUDE_STRSTREAM_HEADERS
64 #endif
66 #include "className.H"
67 #include "fileName.H"
68 #include "stringList.H"
69 #include "wordList.H"
70 #include "primitivePatch.H"
71 #include "PrimitivePatchInterpolation.H"
72 #include "volPointInterpolation.H"
74 #undef VTKPV3FOAM_DUALPORT
76 // * * * * * * * * * * * * * Forward Declarations  * * * * * * * * * * * * * //
78 class vtkDataArraySelection;
79 class vtkDataSet;
80 class vtkPoints;
81 class vtkPV3FoamReader;
82 class vtkRenderer;
83 class vtkTextActor;
84 class vtkMultiBlockDataSet;
85 class vtkPolyData;
86 class vtkUnstructuredGrid;
87 class vtkIndent;
89 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 namespace Foam
94 // Foam class forward declarations
95 class argList;
96 class Time;
97 class fvMesh;
98 class IOobjectList;
99 class polyPatch;
100 class faceSet;
101 class pointSet;
103 template<class Type> class IOField;
104 template<class Type> class List;
106 /*---------------------------------------------------------------------------*\
107                         Class vtkPV3Foam Declaration
108 \*---------------------------------------------------------------------------*/
110 class vtkPV3Foam
112     // Private classes
114         //- Bookkeeping for GUI checklists and the multi-block organization
115         class partInfo
116         {
117             const char *name_;
118             int block_;
119             int start_;
120             int size_;
122         public:
124             partInfo(const char *name, const int blockNo=0)
125             :
126                 name_(name),
127                 block_(blockNo),
128                 start_(-1),
129                 size_(0)
130             {}
132             //- Return the block holding these datasets
133             int block() const
134             {
135                 return block_;
136             }
138             //- Assign block number, return previous value
139             int block(int blockNo)
140             {
141                 int prev = block_;
142                 block_ = blockNo;
143                 return prev;
144             }
146             const char* name() const
147             {
148                 return name_;
149             }
151             int start() const
152             {
153                 return start_;
154             }
156             int end() const
157             {
158                 return start_ + size_;
159             }
161             int size() const
162             {
163                 return size_;
164             }
166             bool empty() const
167             {
168                 return !size_;
169             }
171             void reset()
172             {
173                 start_ = -1;
174                 size_ = 0;
175             }
177             //- Assign new start and reset the size
178             void operator=(const int i)
179             {
180                 start_ = i;
181                 size_ = 0;
182             }
184             //- Increment the size
185             void operator+=(const int n)
186             {
187                 size_ += n;
188             }
189         };
191         //- bookkeeping for polyhedral cell decomposition
192         //  hide in extra pointMap (cellSet/cellZone) for now
193         class polyDecomp
194         {
195             labelList superCells_;
196             labelList addPointCellLabels_;
197             labelList pointMap_;
199         public:
201             polyDecomp()
202             {}
204             //- Label of original cell for decomposed cells
205             labelList& superCells()
206             {
207                 return superCells_;
208             }
210             //- Label of original cell for decomposed cells
211             const labelList& superCells() const
212             {
213                 return superCells_;
214             }
216             //- Cell-centre labels for additional points of decomposed cells
217             labelList& addPointCellLabels()
218             {
219                 return addPointCellLabels_;
220             }
222             //- Cell-centre labels for additional points of decomposed cells
223             const labelList& addPointCellLabels() const
224             {
225                 return addPointCellLabels_;
226             }
228             //- Point labels for subsetted meshes
229             labelList& pointMap()
230             {
231                 return pointMap_;
232             }
234             //- Point labels for subsetted meshes
235             const labelList& pointMap() const
236             {
237                 return pointMap_;
238             }
241             //- Clear
242             void clear()
243             {
244                 superCells_.clear();
245                 addPointCellLabels_.clear();
246                 pointMap_.clear();
247             }
248         };
251     // Private Data
253         //- Access to the controlling vtkPV3FoamReader
254         vtkPV3FoamReader* reader_;
256         //- Foam time control
257         autoPtr<Time> dbPtr_;
259         //- Foam mesh
260         fvMesh* meshPtr_;
262         //- The mesh region
263         word meshRegion_;
265         //- The mesh directory for the region
266         fileName meshDir_;
268         //- The time index
269         int timeIndex_;
271         //- Track changes in mesh geometry
272         bool meshChanged_;
274         //- Track changes in fields
275         bool fieldsChanged_;
277         //- Selected geometrical parts (internalMesh, patches, ...)
278         boolList partStatus_;
280         //- Datasets corresponding to selected geometrical pieces
281         //  a negative number indicates that no vtkmesh exists for this piece
282         labelList partDataset_;
284         //- First instance and size of various mesh parts
285         //  used to index into partStatus_ and partDataset_
286         partInfo partInfoVolume_;
287         partInfo partInfoPatches_;
288         partInfo partInfoLagrangian_;
289         partInfo partInfoCellZones_;
290         partInfo partInfoFaceZones_;
291         partInfo partInfoPointZones_;
292         partInfo partInfoCellSets_;
293         partInfo partInfoFaceSets_;
294         partInfo partInfoPointSets_;
296         //- Decomposed cells information (mesh regions)
297         //  TODO: regions
298         List<polyDecomp> regionPolyDecomp_;
300         //- Decomposed cells information (cellZone meshes)
301         List<polyDecomp> zonePolyDecomp_;
303         //- Decomposed cells information (cellSet meshes)
304         List<polyDecomp> csetPolyDecomp_;
306         //- List of patch names for rendering to window
307         List<vtkTextActor*> patchTextActorsPtrs_;
309     // Private Member Functions
311         // Convenience method use to convert the readers from VTK 5
312         // multiblock API to the current composite data infrastructure
313         static void AddToBlock
314         (
315             vtkMultiBlockDataSet* output,
316             vtkDataSet* dataset,
317             const partInfo&,
318             const label datasetNo,
319             const string& datasetName
320         );
322         // Convenience method use to convert the readers from VTK 5
323         // multiblock API to the current composite data infrastructure
324         static vtkDataSet* GetDataSetFromBlock
325         (
326             vtkMultiBlockDataSet* output,
327             const partInfo&,
328             const label datasetNo
329         );
331         // Convenience method use to convert the readers from VTK 5
332         // multiblock API to the current composite data infrastructure
333         static label GetNumberOfDataSets
334         (
335             vtkMultiBlockDataSet* output,
336             const partInfo&
337         );
339         //- Reset data counters
340         void resetCounters();
342         // Update information helper functions
344             //- Update the mesh parts selected in the GUI
345             void updateMeshPartsStatus();
347             //- Internal mesh info
348             void updateInfoInternalMesh();
350             //- Lagrangian info
351             void updateInfoLagrangian();
353             //- Patch info
354             void updateInfoPatches();
356             //- Set info
357             void updateInfoSets();
359             //- Zone info
360             void updateInfoZones();
362             //- Read zone names for zoneType from file
363             wordList readZoneNames(const word& zoneType);
365             //- Add objects of Type to paraview array selection
366             template<class Type>
367             label addToSelection
368             (
369                 vtkDataArraySelection*,
370                 const IOobjectList&,
371                 const string& suffix=string::null
372             );
374             //- Field info
375             template<template<class> class patchType, class meshType>
376             void updateInfoFields(vtkDataArraySelection*);
378             //- Lagrangian field info
379             void updateInfoLagrangianFields();
382         // Update helper functions
384             //- Foam mesh
385             void updateFoamMesh();
387             //- Reduce memory footprint after conversion
388             void reduceMemory();
390             //- Volume fields
391             void updateVolFields(vtkMultiBlockDataSet*);
393             //- Point fields
394             void updatePointFields(vtkMultiBlockDataSet*);
396             //- Lagrangian fields
397             void updateLagrangianFields(vtkMultiBlockDataSet*);
400         // Mesh conversion functions
402             //- Volume mesh
403             void convertMeshVolume(vtkMultiBlockDataSet*, int& blockNo);
405             //- Lagrangian mesh
406             void convertMeshLagrangian(vtkMultiBlockDataSet*, int& blockNo);
408             //- Patch meshes
409             void convertMeshPatches(vtkMultiBlockDataSet*, int& blockNo);
411             //- Cell zone meshes
412             void convertMeshCellZones(vtkMultiBlockDataSet*, int& blockNo);
414             //- Face zone meshes
415             void convertMeshFaceZones(vtkMultiBlockDataSet*, int& blockNo);
417             //- Point zone meshes
418             void convertMeshPointZones(vtkMultiBlockDataSet*, int& blockNo);
420             //- Cell set meshes
421             void convertMeshCellSets(vtkMultiBlockDataSet*, int& blockNo);
423             //- Face set meshes
424             void convertMeshFaceSets(vtkMultiBlockDataSet*, int& blockNo);
426             //- Point set meshes
427             void convertMeshPointSets(vtkMultiBlockDataSet*, int& blockNo);
430         // Add mesh functions
432             //- Add internal mesh/cell set meshes
433             vtkUnstructuredGrid* volumeVTKMesh(const fvMesh&, polyDecomp&);
435             //- Add Lagrangian mesh
436             vtkPolyData* lagrangianVTKMesh
437             (
438                 const fvMesh&,
439                 const word& cloudName
440             );
442             //- Add patch mesh
443             vtkPolyData* patchVTKMesh(const polyPatch&);
445             //- Add face zone mesh
446             vtkPolyData* faceZoneVTKMesh
447             (
448                 const fvMesh&,
449                 const labelList& faceLabels
450             );
452             //- Add point zone
453             vtkPolyData* pointZoneVTKMesh
454             (
455                 const fvMesh&,
456                 const labelList& pointLabels
457             );
459             //- Add face set mesh
460             vtkPolyData* faceSetVTKMesh
461             (
462                 const fvMesh&,
463                 const faceSet&
464             );
466             //- Add point mesh
467             vtkPolyData* pointSetVTKMesh
468             (
469                 const fvMesh&,
470                 const pointSet&
471             );
473         // Field conversion functions
475             //- Convert volume fields
476             void convertVolFields(vtkMultiBlockDataSet*);
478             //- Convert point fields
479             void convertPointFields(vtkMultiBlockDataSet*);
481             //- Convert Lagrangian fields
482             void convertLagrangianFields(vtkMultiBlockDataSet*);
485         //- Add the fields in the selected time directory to the selection
486         //  lists
487         template<class GeoField>
488         label addObjectsToSelection
489         (
490             vtkDataArraySelection*,
491             const IOobjectList&,
492             const string& suffix=string::null
493         );
496         // Convert Foam fields
498             //- Volume fields - all types
499             template<class Type>
500             void convertVolFields
501             (
502                 const fvMesh&,
503                 const PtrList<PrimitivePatchInterpolation<primitivePatch> >&,
504                 const IOobjectList&,
505                 vtkMultiBlockDataSet* output
506             );
508             //- Volume field - all selected parts
509             template<class Type>
510             void convertVolFieldBlock
511             (
512                 const GeometricField<Type, fvPatchField, volMesh>&,
513                 autoPtr<GeometricField<Type, pointPatchField, pointMesh> >&,
514                 vtkMultiBlockDataSet* output,
515                 const partInfo& selector,
516                 const List<polyDecomp>& decompLst
517             );
519             //- Volume field
520             template<class Type>
521             void convertVolField
522             (
523                 const GeometricField<Type, fvPatchField, volMesh>&,
524                 vtkMultiBlockDataSet* output,
525                 const partInfo&,
526                 const label datasetNo,
527                 const polyDecomp&
528             );
530             //- Patch field
531             template<class Type>
532             void convertPatchField
533             (
534                 const word& name,
535                 const Field<Type>&,
536                 vtkMultiBlockDataSet* output,
537                 const partInfo&,
538                 const label datasetNo
539             );
541             //- face set/zone field
542             template<class Type>
543             void convertFaceField
544             (
545                 const GeometricField<Type, fvPatchField, volMesh>&,
546                 vtkMultiBlockDataSet* output,
547                 const partInfo&,
548                 const label datasetNo,
549                 const fvMesh&,
550                 const labelList& faceLabels
551             );
553             //- face set/zone field
554             template<class Type>
555             void convertFaceField
556             (
557                 const GeometricField<Type, fvPatchField, volMesh>&,
558                 vtkMultiBlockDataSet* output,
559                 const partInfo&,
560                 const label datasetNo,
561                 const fvMesh&,
562                 const faceSet&
563             );
565             //- Lagrangian fields - all types
566             template<class Type>
567             void convertLagrangianFields
568             (
569                 const IOobjectList&,
570                 vtkMultiBlockDataSet* output,
571                 const label datasetNo
572             );
574             //- Lagrangian field
575             template<class Type>
576             void convertLagrangianField
577             (
578                 const IOField<Type>&,
579                 vtkMultiBlockDataSet* output,
580                 const partInfo&,
581                 const label datasetNo
582             );
584             //- Point fields - all types
585             template<class Type>
586             void convertPointFields
587             (
588                 const fvMesh&,
589                 const pointMesh&,
590                 const IOobjectList&,
591                 vtkMultiBlockDataSet* output
592             );
594             //- Point field - all selected parts
595             template<class Type>
596             void convertPointFieldBlock
597             (
598                 const GeometricField<Type, pointPatchField, pointMesh>&,
599                 vtkMultiBlockDataSet* output,
600                 const partInfo& selector,
601                 const List<polyDecomp>&
602             );
604             //- Point fields
605             template<class Type>
606             void convertPointField
607             (
608                 const GeometricField<Type, pointPatchField, pointMesh>&,
609                 const GeometricField<Type, fvPatchField, volMesh>&,
610                 vtkMultiBlockDataSet* output,
611                 const partInfo&,
612                 const label datasetNo,
613                 const polyDecomp&
614             );
616             //- Patch point field
617             template<class Type>
618             void convertPatchPointField
619             (
620                 const word& name,
621                 const Field<Type>&,
622                 vtkMultiBlockDataSet* output,
623                 const partInfo&,
624                 const label datasetNo
625             );
628        // GUI selection helper functions
630             //- Extract up to the first non-word characters
631             inline static word getFirstWord(const char*);
633             //- Only keep what is listed in hashSet
634             static void pruneObjectList
635             (
636                 IOobjectList&,
637                 const wordHashSet&
638             );
640             //- Retrieve the current selections
641             static wordHashSet getSelected(vtkDataArraySelection*);
643             //- Retrieve a sub-list of the current selections
644             static wordHashSet getSelected
645             (
646                 vtkDataArraySelection*,
647                 const partInfo&
648             );
650             //- Retrieve the current selections
651             static stringList getSelectedArrayEntries(vtkDataArraySelection*);
653             //- Retrieve a sub-list of the current selections
654             static stringList getSelectedArrayEntries
655             (
656                 vtkDataArraySelection*,
657                 const partInfo&
658             );
660             //- Set selection(s)
661             static void setSelectedArrayEntries
662             (
663                 vtkDataArraySelection*,
664                 const stringList&
665             );
667             //- Get the first word from the mesh parts selection
668             word getPartName(int);
671         //- Disallow default bitwise copy construct
672         vtkPV3Foam(const vtkPV3Foam&);
674         //- Disallow default bitwise assignment
675         void operator=(const vtkPV3Foam&);
678 public:
680     //- Static data members
682         ClassName("vtkPV3Foam");
685     // Constructors
687         //- Construct from components
688         vtkPV3Foam
689         (
690             const char* const FileName,
691             vtkPV3FoamReader* reader
692         );
695     //- Destructor
697         ~vtkPV3Foam();
700     // Member Functions
702         //- Update
703         void updateInfo();
705         void Update
706         (
707             vtkMultiBlockDataSet* output,
708             vtkMultiBlockDataSet* lagrangianOutput
709         );
711         //- Clean any storage
712         void CleanUp();
714         //- Allocate and return a list of selected times
715         //  returns the count via the parameter
716         double* findTimes(int& nTimeSteps);
718         //- Add patch names to the display
719         void addPatchNames(vtkRenderer* renderer);
721         //- Remove patch names from the display
722         void removePatchNames(vtkRenderer* renderer);
724         //- set the runTime to the first plausible request time,
725         //  returns the timeIndex
726         //  sets to "constant" on error
727         int setTime(int count, const double requestTimes[]);
730         //- The current time index
731         int timeIndex() const
732         {
733            return timeIndex_;
734         }
737      // Access
739         //- Debug information
740         void PrintSelf(ostream&, vtkIndent) const;
742         //- Simple memory used debugging information
743         static void printMemory();
748 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
750 } // End namespace Foam
752 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
754 #    include "vtkPV3FoamI.H"
756 #endif
758 // ************************************************************************* //