VTK  9.5.20250803
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
32#ifndef vtkFLUENTCFFReader_h
33#define vtkFLUENTCFFReader_h
34
35#include <memory> // std::unique_ptr
36
37#include "vtkIOFLUENTCFFModule.h" // For export macro
38
40#include "vtkNew.h" // For vtkNew
41
42VTK_ABI_NAMESPACE_BEGIN
44class vtkPoints;
45class vtkTriangle;
46class vtkTetra;
47class vtkQuad;
48class vtkHexahedron;
49class vtkPyramid;
50class vtkWedge;
51
52class VTKIOFLUENTCFF_EXPORT vtkFLUENTCFFReader : public vtkMultiBlockDataSetAlgorithm
53{
54public:
57 void PrintSelf(ostream& os, vtkIndent indent) override;
58
60
63 vtkSetMacro(FileName, std::string);
64 vtkGetMacro(FileName, std::string);
66
68
72 vtkGetMacro(NumberOfCells, vtkIdType);
74
79
84 const char* GetCellArrayName(int index);
85
87
91 int GetCellArrayStatus(const char* name);
92 void SetCellArrayStatus(const char* name, int status);
94
96
102
104 //
105 // Structures
106 //
107 struct Cell
108 {
109 int type;
110 int zone;
111 std::vector<int> faces;
113 int child;
114 std::vector<int> nodes;
115 std::vector<int> nodesOffset;
116 std::vector<int> childId;
117 };
118 struct Face
119 {
120 int type;
121 unsigned int zone;
122 std::vector<int> nodes;
123 int c0;
124 int c1;
127 int child;
132 };
134 {
135 std::string variableName;
137 std::vector<double> scalarData;
138 };
140 {
141 std::string variableName;
143 size_t dim;
144 std::vector<double> vectorData;
145 };
147
148protected:
153
158 {
159 NOT_LOADED = 0,
160 AVAILABLE = 1,
161 LOADED = 2,
162 ERROR = 3
163 };
164
166
169 virtual bool OpenCaseFile(const std::string& filename);
170 virtual DataState OpenDataFile(const std::string& filename);
172
176 virtual void GetNumberOfCellZones();
177
181 virtual int ParseCaseFile();
182
186 virtual int GetDimension();
187
189
192 virtual void GetNodesGlobal();
193 virtual void GetCellsGlobal();
194 virtual void GetFacesGlobal();
196
200 virtual void GetNodes();
201
205 virtual void GetCells();
206
210 virtual void GetFaces();
211
217
222 virtual void GetCellOverset();
223
227 virtual void GetCellTree();
228
232 virtual void GetFaceTree();
233
238
244
248 virtual void CleanCells();
249
251
255 virtual void PopulateCellNodes();
256 virtual void PopulateCellTree();
258
260
263 virtual void PopulateTriangleCell(int i);
264 virtual void PopulateTetraCell(int i);
265 virtual void PopulateQuadCell(int i);
266 virtual void PopulateHexahedronCell(int i);
267 virtual void PopulatePyramidCell(int i);
268 virtual void PopulateWedgeCell(int i);
269 virtual void PopulatePolyhedronCell(int i);
271
275 virtual int GetData();
276
280 virtual int GetMetaData();
281
282 //
283 // Variables
284 //
286 std::string FileName;
287 vtkIdType NumberOfCells = 0;
288 int NumberOfCellArrays = 0;
289
290 struct vtkInternals;
291 std::unique_ptr<vtkInternals> HDFImpl;
292
300
301 std::vector<Cell> Cells;
302 std::vector<Face> Faces;
303 std::vector<int> CellZones;
304 std::vector<ScalarDataChunk> ScalarDataChunks;
305 std::vector<VectorDataChunk> VectorDataChunks;
306 std::vector<std::string> PreReadScalarData;
307 std::vector<std::string> PreReadVectorData;
308
309 vtkTypeBool SwapBytes = 0;
310 int GridDimension = 0;
311 DataState FileState = DataState::NOT_LOADED;
312 int NumberOfScalars = 0;
313 int NumberOfVectors = 0;
314
315private:
316 vtkFLUENTCFFReader(const vtkFLUENTCFFReader&) = delete;
317 void operator=(const vtkFLUENTCFFReader&) = delete;
318};
319VTK_ABI_NAMESPACE_END
320#endif
Store on/off settings for data arrays, etc.
reads a dataset in Fluent CFF file format
vtkNew< vtkTetra > Tetra
virtual void GetCellTree()
Get the tree (AMR) cell topology.
virtual void GetPeriodicShadowFaces()
Get the periodic shadown faces information !!! NOT IMPLEMENTED YET !!!
virtual DataState OpenDataFile(const std::string &filename)
std::vector< std::string > PreReadScalarData
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.
std::vector< VectorDataChunk > VectorDataChunks
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
std::vector< std::string > PreReadVectorData
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.
std::vector< ScalarDataChunk > ScalarDataChunks
virtual void GetFaceTree()
Get the tree (AMR) face topology.
virtual void GetNonconformalGridInterfaceFaceInformation()
Get the non conformal grid interface information !!! NOT IMPLEMENTED YET !!!
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 !!! NOT IMPLEMENTED YET !!!
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.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
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
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:113
a cell that represents a triangle
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:332