vtkPVReaders.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 Namespace
25  Foam::vtkPVReaders
26 
27 Description
28  A collection of helper functions when building a reader interface in
29  ParaView3.
30 
31 SourceFiles
32  vtkPVReaders.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef vtkPVReaders_H
37 #define vtkPVReaders_H
38 
39 // do not include legacy strstream headers
40 #ifndef VTK_EXCLUDE_STRSTREAM_HEADERS
41 # define VTK_EXCLUDE_STRSTREAM_HEADERS
42 #endif
43 
44 #include "className.H"
45 #include "fileName.H"
46 #include "stringList.H"
47 #include "wordList.H"
48 #include "HashSet.H"
49 
50 
51 // * * * * * * * * * * * * * Forward Declarations * * * * * * * * * * * * * //
52 
53 class vtkDataArraySelection;
54 class vtkDataSet;
55 class vtkPoints;
56 class vtkPVFoamReader;
57 class vtkRenderer;
58 class vtkTextActor;
59 class vtkMultiBlockDataSet;
60 class vtkPolyData;
61 class vtkUnstructuredGrid;
62 class vtkIndent;
63 
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 namespace vtkPVReaders
70 {
71  //- Declare name of the class and its debug switch
72  NamespaceName("vtkPVReaders");
73 
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_;
81 
82  public:
83 
84  partInfo(const char *name, const int blockNo=0)
85  :
86  name_(name),
87  block_(blockNo),
88  start_(-1),
89  size_(0)
90  {}
91 
92  //- Return the block holding these datasets
93  int block() const
94  {
95  return block_;
96  }
97 
98  //- Assign block number, return previous value
99  int block(int blockNo)
100  {
101  int prev = block_;
102  block_ = blockNo;
103  return prev;
104  }
105 
106  const char* name() const
107  {
108  return name_;
109  }
110 
111  int start() const
112  {
113  return start_;
114  }
115 
116  int end() const
117  {
118  return start_ + size_;
119  }
120 
121  int size() const
122  {
123  return size_;
124  }
125 
126  bool empty() const
127  {
128  return !size_;
129  }
130 
131  void reset()
132  {
133  start_ = -1;
134  size_ = 0;
135  }
136 
137  //- Assign new start and reset the size
138  void operator=(const int i)
139  {
140  start_ = i;
141  size_ = 0;
142  }
143 
144  //- Increment the size
145  void operator+=(const int n)
146  {
147  size_ += n;
148  }
149  };
150 
151 
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  );
162 
163 
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  );
172 
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
177  (
178  vtkMultiBlockDataSet* output,
179  const partInfo& selector
180  );
181 
182 
183  //- Retrieve the current selections as a wordHashSet
185  (
186  vtkDataArraySelection* select
187  );
188 
189 
190  //- Retrieve a sub-list of the current selections
192  (
193  vtkDataArraySelection*,
194  const partInfo&
195  );
196 
197 
198  //- Retrieve the current selections
199  stringList getSelectedArrayEntries(vtkDataArraySelection*);
200 
201  //- Retrieve a sub-list of the current selections
203  (
204  vtkDataArraySelection* select,
205  const partInfo& selector
206  );
207 
208 
209  //- Set selection(s)
211  (
212  vtkDataArraySelection*,
213  const stringList&
214  );
215 
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 } // End namespace vtkPV
221 
222 } // End namespace Foam
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 #endif
227 
228 // ************************************************************************* //
A HashTable with keys but without contents.
Definition: HashSet.H:59
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
NamespaceName("vtkPVReaders")
Declare name of the class and its debug switch.
label GetNumberOfDataSets(vtkMultiBlockDataSet *output, const partInfo &selector)
Convenience method use to convert the readers from VTK 5.
Definition: vtkPVReaders.C:162
partInfo(const char *name, const int blockNo=0)
Definition: vtkPVReaders.H:84
void operator=(const int i)
Assign new start and reset the size.
Definition: vtkPVReaders.H:138
vtkDataSet * GetDataSetFromBlock(vtkMultiBlockDataSet *output, const partInfo &selector, const label datasetNo)
Convenience method use to convert the readers from VTK 5.
Definition: vtkPVReaders.C:140
const char * name() const
Definition: vtkPVReaders.H:106
int block() const
Return the block holding these datasets.
Definition: vtkPVReaders.H:93
wordHashSet getSelected(vtkDataArraySelection *select)
Retrieve the current selections as a wordHashSet.
Definition: vtkPVReaders.C:187
int block(int blockNo)
Assign block number, return previous value.
Definition: vtkPVReaders.H:99
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Bookkeeping for GUI checklists and the multi-block organization.
Definition: vtkPVReaders.H:75
label n
void setSelectedArrayEntries(vtkDataArraySelection *, const stringList &)
Set selection(s)
Definition: vtkPVReaders.C:306
void AddToBlock(vtkMultiBlockDataSet *output, vtkDataSet *dataset, const partInfo &selector, const label datasetNo, const std::string &datasetName)
Convenience method use to convert the readers from VTK 5.
Definition: vtkPVReaders.C:80
Namespace for OpenFOAM.
void operator+=(const int n)
Increment the size.
Definition: vtkPVReaders.H:145
stringList getSelectedArrayEntries(vtkDataArraySelection *)
Retrieve the current selections.
Definition: vtkPVReaders.C:228