vtkPVblockMesh.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-2021 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::vtkPVblockMesh
26 
27 Description
28  Provides a reader interface for OpenFOAM blockMesh to VTK interaction
29 
30 SourceFiles
31  vtkPVblockMesh.C
32  vtkPVblockMeshConvert.C
33  vtkPVblockMeshUpdate.C
34  vtkPVblockMeshUtils.C
35 
36  // Needed by VTK:
37  vtkDataArrayTemplateImplicit.txx
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef vtkPVblockMesh_H
42 #define vtkPVblockMesh_H
43 
44 // do not include legacy strstream headers
45 #ifndef VTK_EXCLUDE_STRSTREAM_HEADERS
46  #define VTK_EXCLUDE_STRSTREAM_HEADERS
47 #endif
48 
49 #include "className.H"
50 #include "fileName.H"
51 #include "stringList.H"
52 #include "wordList.H"
53 
54 #include "primitivePatch.H"
55 
56 // * * * * * * * * * * * * * Forward Declarations * * * * * * * * * * * * * //
57 
58 class vtkDataArraySelection;
59 class vtkDataSet;
60 class vtkPoints;
61 class vtkPVblockMeshReader;
62 class vtkRenderer;
63 class vtkTextActor;
64 class vtkMultiBlockDataSet;
65 class vtkPolyData;
66 class vtkUnstructuredGrid;
67 class vtkIndent;
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 
74 // OpenFOAM class forward declarations
75 class argList;
76 class Time;
77 class blockMesh;
78 
79 template<class Type> class List;
80 
81 /*---------------------------------------------------------------------------*\
82  Class vtkPVblockMesh Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 class vtkPVblockMesh
86 {
87  // Private classes
88 
89  //- Bookkeeping for GUI checklists and the multi-block organisation
90  class arrayRange
91  {
92  const char *name_;
93  int block_;
94  int start_;
95  int size_;
96 
97  public:
98 
99  arrayRange(const char *name, const int blockNo=0)
100  :
101  name_(name),
102  block_(blockNo),
103  start_(0),
104  size_(0)
105  {}
106 
107  //- Return the block holding these datasets
108  int block() const
109  {
110  return block_;
111  }
112 
113  //- Assign block number, return previous value
114  int block(int blockNo)
115  {
116  int prev = block_;
117  block_ = blockNo;
118  return prev;
119  }
120 
121  //- Return block name
122  const char* name() const
123  {
124  return name_;
125  }
126 
127  //- Return array start index
128  int start() const
129  {
130  return start_;
131  }
132 
133  //- Return array end index
134  int end() const
135  {
136  return start_ + size_;
137  }
138 
139  //- Return sublist size
140  int size() const
141  {
142  return size_;
143  }
144 
145  bool empty() const
146  {
147  return !size_;
148  }
149 
150  //- Reset the size to zero and optionally assign a new start
151  void reset(const int startAt = 0)
152  {
153  start_ = startAt;
154  size_ = 0;
155  }
156 
157  //- Increment the size
158  void operator+=(const int n)
159  {
160  size_ += n;
161  }
162  };
163 
164 
165  // Private Data
166 
167  //- Access to the controlling vtkPVblockMeshReader
168  vtkPVblockMeshReader* reader_;
169 
170  //- OpenFOAM time control
171  autoPtr<Time> dbPtr_;
172 
173  //- OpenFOAM mesh
174  blockMesh* meshPtr_;
175 
176  //- The mesh region
177  word meshRegion_;
178 
179  //- The mesh directory for the region
180  fileName meshDir_;
181 
182  //- Selected geometrical parts
183  boolList blockStatus_;
184 
185  //- Selected curved edges
186  boolList edgeStatus_;
187 
188  //- First instance and size of blockMesh blocks
189  // used to index into blockStatus_
190  arrayRange arrayRangeBlocks_;
191 
192  //- First instance and size of CurvedEdges (only partially used)
193  arrayRange arrayRangeEdges_;
194 
195  //- First instance and size of block corners (only partially used)
196  arrayRange arrayRangeCorners_;
197 
198  //- List of point numbers for rendering to window
199  List<vtkTextActor*> pointNumberTextActorsPtrs_;
200 
201  // Private Member Functions
202 
203  // Convenience method use to convert the readers from VTK 5
204  // multiblock API to the current composite data infrastructure
205  static void AddToBlock
206  (
207  vtkMultiBlockDataSet* output,
208  vtkDataSet* dataset,
209  const arrayRange&,
210  const label datasetNo,
211  const std::string& datasetName
212  );
213 
214  // Convenience method use to convert the readers from VTK 5
215  // multiblock API to the current composite data infrastructure
216  static vtkDataSet* GetDataSetFromBlock
217  (
218  vtkMultiBlockDataSet* output,
219  const arrayRange&,
220  const label datasetNo
221  );
222 
223  // Convenience method use to convert the readers from VTK 5
224  // multiblock API to the current composite data infrastructure
225  static label GetNumberOfDataSets
226  (
227  vtkMultiBlockDataSet* output,
228  const arrayRange&
229  );
230 
231  //- Update boolList from GUI selection
232  static void updateBoolListStatus
233  (
234  boolList&,
235  vtkDataArraySelection*
236  );
237 
238  //- Reset data counters
239  void resetCounters();
240 
241  // Update information helper functions
242 
243  //- Internal block info
244  void updateInfoBlocks(vtkDataArraySelection*);
245 
246  //- Block curved edges info
247  void updateInfoEdges(vtkDataArraySelection*);
248 
249  // Update helper functions
250 
251  //- OpenFOAM mesh
252  void updateFoamMesh();
253 
254  // Mesh conversion functions
255 
256  //- Mesh blocks
257  void convertMeshBlocks(vtkMultiBlockDataSet*, int& blockNo);
258 
259  //- Mesh curved edges
260  void convertMeshEdges(vtkMultiBlockDataSet*, int& blockNo);
261 
262  //- Mesh corners
263  void convertMeshCorners(vtkMultiBlockDataSet*, int& blockNo);
264 
265 
266  // GUI selection helper functions
267 
268  //- Retrieve the current selections
269  static wordHashSet getSelected(vtkDataArraySelection*);
270 
271  //- Retrieve a sub-list of the current selections
272  static wordHashSet getSelected
273  (
274  vtkDataArraySelection*,
275  const arrayRange&
276  );
277 
278  //- Retrieve the current selections
279  static stringList getSelectedArrayEntries(vtkDataArraySelection*);
280 
281  //- Retrieve a sub-list of the current selections
282  static stringList getSelectedArrayEntries
283  (
284  vtkDataArraySelection*,
285  const arrayRange&
286  );
287 
288  //- Set selection(s)
289  static void setSelectedArrayEntries
290  (
291  vtkDataArraySelection*,
292  const stringList&
293  );
294 
295 
296 public:
297 
298  //- Static data members
299 
300  ClassName("vtkPVblockMesh");
301 
302 
303  // Constructors
304 
305  //- Construct from components
307  (
308  const char* const FileName,
309  vtkPVblockMeshReader* reader
310  );
311 
312  //- Disallow default bitwise copy construction
313  vtkPVblockMesh(const vtkPVblockMesh&) = delete;
314 
315 
316  //- Destructor
317  ~vtkPVblockMesh();
318 
319 
320  // Member Functions
321 
322  //- Update
323  void updateInfo();
324 
325  void Update(vtkMultiBlockDataSet* output);
326 
327  //- Clean any storage
328  void CleanUp();
329 
330  //- Add/remove point numbers to/from the view
331  void renderPointNumbers(vtkRenderer*, const bool show);
332 
333  // Access
334 
335  //- Debug information
336  void PrintSelf(ostream&, vtkIndent) const;
337 
338 
339  // Member Operators
340 
341  //- Disallow default bitwise assignment
342  void operator=(const vtkPVblockMesh&) = delete;
343 };
344 
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 } // End namespace Foam
349 
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351 
352 #endif
353 
354 // ************************************************************************* //
label n
A HashTable with keys but without contents.
Definition: HashSet.H:62
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A multi-block mesh generator.
Definition: blockMesh.H:65
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:66
A class for handling file names.
Definition: fileName.H:82
Provides a reader interface for OpenFOAM blockMesh to VTK interaction.
~vtkPVblockMesh()
Destructor.
ClassName("vtkPVblockMesh")
Static data members.
void Update(vtkMultiBlockDataSet *output)
void PrintSelf(ostream &, vtkIndent) const
Debug information.
void CleanUp()
Clean any storage.
void renderPointNumbers(vtkRenderer *, const bool show)
Add/remove point numbers to/from the view.
vtkPVblockMesh(const char *const FileName, vtkPVblockMeshReader *reader)
Construct from components.
void operator=(const vtkPVblockMesh &)=delete
Disallow default bitwise assignment.
void updateInfo()
Update.
A class for handling words, derived from string.
Definition: word.H:62
Macro definitions for declaring ClassName(), NamespaceName(), etc.
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