ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / graphics / PV3Readers / vtkPV3Readers / vtkPV3Readers.H
blobe955be79fc86aacd6e39486c94a3dae6c0477748
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Namespace
25     Foam::vtkPV3Readers
27 Description
28     A collection of helper functions when building a reader interface in
29     ParaView3.
31 SourceFiles
32     vtkPV3Readers.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef vtkPV3Readers_H
37 #define vtkPV3Readers_H
39 // do not include legacy strstream headers
40 #ifndef  VTK_EXCLUDE_STRSTREAM_HEADERS
41 # define VTK_EXCLUDE_STRSTREAM_HEADERS
42 #endif
44 #include "className.H"
45 #include "fileName.H"
46 #include "stringList.H"
47 #include "wordList.H"
48 #include "HashSet.H"
51 // * * * * * * * * * * * * * Forward Declarations  * * * * * * * * * * * * * //
53 class vtkDataArraySelection;
54 class vtkDataSet;
55 class vtkPoints;
56 class vtkPV3FoamReader;
57 class vtkRenderer;
58 class vtkTextActor;
59 class vtkMultiBlockDataSet;
60 class vtkPolyData;
61 class vtkUnstructuredGrid;
62 class vtkIndent;
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 namespace Foam
69 namespace vtkPV3Readers
71     //- Declare name of the class and its debug switch
72     NamespaceName("vtkPV3Readers");
74     //- Bookkeeping for GUI checklists and the multi-block organization
75     class partInfo
76     {
77         const char *name_;
78         int block_;
79         int start_;
80         int size_;
82     public:
84         partInfo(const char *name, const int blockNo=0)
85         :
86             name_(name),
87             block_(blockNo),
88             start_(-1),
89             size_(0)
90         {}
92         //- Return the block holding these datasets
93         int block() const
94         {
95             return block_;
96         }
98         //- Assign block number, return previous value
99         int block(int blockNo)
100         {
101             int prev = block_;
102             block_ = blockNo;
103             return prev;
104         }
106         const char* name() const
107         {
108             return name_;
109         }
111         int start() const
112         {
113             return start_;
114         }
116         int end() const
117         {
118             return start_ + size_;
119         }
121         int size() const
122         {
123             return size_;
124         }
126         bool empty() const
127         {
128             return !size_;
129         }
131         void reset()
132         {
133             start_ = -1;
134             size_ = 0;
135         }
137         //- Assign new start and reset the size
138         void operator=(const int i)
139         {
140             start_ = i;
141             size_ = 0;
142         }
144         //- Increment the size
145         void operator+=(const int n)
146         {
147             size_ += n;
148         }
149     };
152     //- Convenience method use to convert the readers from VTK 5
153     //  multiblock API to the current composite data infrastructure
154     void AddToBlock
155     (
156         vtkMultiBlockDataSet* output,
157         vtkDataSet* dataset,
158         const partInfo& selector,
159         const label datasetNo,
160         const std::string& datasetName
161     );
164     //- Convenience method use to convert the readers from VTK 5
165     // multiblock API to the current composite data infrastructure
166     vtkDataSet* GetDataSetFromBlock
167     (
168         vtkMultiBlockDataSet* output,
169         const partInfo& selector,
170         const label datasetNo
171     );
173     //- Convenience method use to convert the readers from VTK 5
174     // multiblock API to the current composite data infrastructure
175     // ununsed at the moment
176     label GetNumberOfDataSets
177     (
178         vtkMultiBlockDataSet* output,
179         const partInfo& selector
180     );
183     //- Retrieve the current selections as a wordHashSet
184     wordHashSet getSelected
185     (
186         vtkDataArraySelection* select
187     );
190     //- Retrieve a sub-list of the current selections
191     wordHashSet getSelected
192     (
193         vtkDataArraySelection*,
194         const partInfo&
195     );
198     //- Retrieve the current selections
199     stringList getSelectedArrayEntries(vtkDataArraySelection*);
201     //- Retrieve a sub-list of the current selections
202     stringList getSelectedArrayEntries
203     (
204         vtkDataArraySelection* select,
205         const partInfo& selector
206     );
209     //- Set selection(s)
210     void setSelectedArrayEntries
211     (
212         vtkDataArraySelection*,
213         const stringList&
214     );
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 } // End namespace vtkPV3
222 } // End namespace Foam
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 #endif
228 // ************************************************************************* //