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