VTK  9.6.20260211
vtkFLUENTCFFReader.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
37
38#ifndef vtkFLUENTCFFReader_h
39#define vtkFLUENTCFFReader_h
40
41#include <memory> // std::unique_ptr
42
43#include "vtkIOFLUENTCFFModule.h" // For export macro
44
46#include "vtkNew.h" // For vtkNew
47
48VTK_ABI_NAMESPACE_BEGIN
50class vtkPoints;
51class vtkTriangle;
52class vtkTetra;
53class vtkQuad;
54class vtkHexahedron;
55class vtkPyramid;
57class vtkWedge;
58
59class VTKIOFLUENTCFF_EXPORT vtkFLUENTCFFReader : public vtkMultiBlockDataSetAlgorithm
60{
61public:
64 void PrintSelf(ostream& os, vtkIndent indent) override;
65
67
70 vtkSetMacro(FileName, std::string);
71 vtkGetMacro(FileName, std::string);
73
75
80 vtkSetMacro(RenameArrays, bool);
81 vtkGetMacro(RenameArrays, bool);
83
85
89 vtkGetMacro(NumberOfCells, vtkIdType);
91
96
101 const char* GetCellArrayName(int index);
102
104
108 int GetCellArrayStatus(const char* name);
109 void SetCellArrayStatus(const char* name, int status);
111
113
119
121 //
122 // Structures
123 //
124 struct Cell
125 {
126 int type;
127 int zone;
128 std::vector<int> faces;
130 int child;
131 std::vector<int> nodes;
132 std::vector<int> nodesOffset;
133 std::vector<int> childId;
134 };
135 struct Face
136 {
137 int type;
138 unsigned int zone;
139 std::vector<int> nodes;
140 int c0;
141 int c1;
144 int child;
149 };
150
151
152protected:
157
162 {
167 };
168
170
173 virtual bool OpenCaseFile(const std::string& filename);
174 virtual DataState OpenDataFile(const std::string& filename);
176
180 virtual void GetNumberOfCellZones();
181
185 virtual int ParseCaseFile();
186
190 virtual int GetDimension();
191
193
196 virtual void GetNodesGlobal();
197 virtual void GetCellsGlobal();
198 virtual void GetFacesGlobal();
200
204 virtual void GetNodes();
205
209 virtual void GetCells();
210
214 virtual void GetFaces();
215
221
226 virtual void GetCellOverset();
227
231 virtual void GetCellTree();
232
236 virtual void GetFaceTree();
237
242
248
252 virtual void CleanCells();
253
255
259 virtual void PopulateCellNodes();
260 virtual void PopulateCellTree();
262
264
267 virtual void PopulateTriangleCell(int i);
268 virtual void PopulateTetraCell(int i);
269 virtual void PopulateQuadCell(int i);
270 virtual void PopulateHexahedronCell(int i);
271 virtual void PopulatePyramidCell(int i);
272 virtual void PopulateWedgeCell(int i);
273 virtual void PopulatePolyhedronCell(int i);
275
279 virtual int GetData();
280
284 virtual int GetMetaData();
285
286 //
287 // Variables
288 //
290 std::string FileName;
291 bool RenameArrays = false;
294
295 struct vtkInternals;
296 std::unique_ptr<vtkInternals> HDFImpl;
297
305
306 std::vector<Cell> Cells;
307 std::vector<Face> Faces;
308 std::vector<int> CellZones;
309
313
314private:
315 vtkFLUENTCFFReader(const vtkFLUENTCFFReader&) = delete;
316 void operator=(const vtkFLUENTCFFReader&) = delete;
317
318 struct DataChunk
319 {
320 std::string variableName;
321 vtkIdType zoneId;
322 size_t dim;
323 std::vector<double> dataVector;
324 };
325
326 std::vector<DataChunk> DataChunks;
327 std::vector<std::string> PreReadData;
328 int NumberOfArrays = 0;
329
333 void ParseUDMData(
334 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& grid, const DataChunk& vectorDataChunk);
335};
336VTK_ABI_NAMESPACE_END
337#endif
Store on/off settings for data arrays, etc.
vtkNew< vtkTetra > Tetra
virtual void GetCellTree()
Get the tree (AMR) cell topology.
virtual void GetPeriodicShadowFaces()
Get the periodic shadown faces information !
virtual DataState OpenDataFile(const std::string &filename)
void DisableAllCellArrays()
Turn on/off all cell arrays.
virtual void PopulatePolyhedronCell(int i)
std::unique_ptr< vtkInternals > HDFImpl
std::vector< Face > Faces
const char * GetCellArrayName(int index)
Get the name of the cell array with the given index in the input.
virtual void CleanCells()
Removes unnecessary faces from the cells.
virtual void PopulatePyramidCell(int i)
vtkNew< vtkQuad > Quad
vtkNew< vtkWedge > Wedge
virtual void PopulateHexahedronCell(int i)
virtual void PopulateCellTree()
virtual int GetData()
Read and reconstruct data from .dat.h5 file.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual void GetNodes()
Get the size and index of node per zone.
int GetCellArrayStatus(const char *name)
Get/Set whether the cell array with the given name is to be read.
virtual void GetFacesGlobal()
~vtkFLUENTCFFReader() override
virtual int GetMetaData()
Pre-read variable name data available for selection.
virtual void GetNumberOfCellZones()
Retrieve the number of cell zones.
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual void PopulateTetraCell(int i)
vtkNew< vtkTriangle > Triangle
virtual void PopulateCellNodes()
Reconstruct and convert the Fluent data format to the VTK format.
virtual void GetCells()
Get the topology of cell per zone.
virtual void GetNodesGlobal()
Get the total number of nodes/cells/faces.
virtual void PopulateQuadCell(int i)
virtual void GetFaces()
Get the topology of face per zone.
virtual void GetFaceTree()
Get the tree (AMR) face topology.
virtual void GetNonconformalGridInterfaceFaceInformation()
Get the non conformal grid interface information !
vtkNew< vtkPyramid > Pyramid
vtkNew< vtkPoints > Points
vtkNew< vtkHexahedron > Hexahedron
int GetNumberOfCellArrays()
Get the number of cell arrays available in the input.
vtkNew< vtkDataArraySelection > CellDataArraySelection
virtual void GetCellsGlobal()
void SetCellArrayStatus(const char *name, int status)
void EnableAllCellArrays()
virtual int ParseCaseFile()
Reads necessary information from the .cas file.
virtual int GetDimension()
Get the dimension of the file (2D/3D)
virtual void PopulateWedgeCell(int i)
virtual void PopulateTriangleCell(int i)
Reconstruct VTK cell topology from Fluent format.
std::vector< Cell > Cells
static vtkFLUENTCFFReader * New()
virtual void GetInterfaceFaceParents()
Get the interface id of parent faces.
virtual bool OpenCaseFile(const std::string &filename)
Open the HDF5 file structure.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual void GetCellOverset()
Get the overset cells information !
std::vector< int > CellZones
a cell that represents a linear 3D hexahedron
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Allocate and hold a VTK object.
Definition vtkNew.h:167
represent and manipulate 3D points
Definition vtkPoints.h:139
a 3D cell that represents a linear pyramid
Definition vtkPyramid.h:95
a cell that represents a 2D quadrilateral
Definition vtkQuad.h:87
Hold a reference to a vtkObjectBase instance.
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:113
a cell that represents a triangle
dataset represents arbitrary combinations of all possible cell types
a 3D cell that represents a linear wedge
Definition vtkWedge.h:85
std::vector< int > nodesOffset
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:354