VTK  9.5.20251217
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
33#ifndef vtkFLUENTCFFReader_h
34#define vtkFLUENTCFFReader_h
35
36#include <memory> // std::unique_ptr
37
38#include "vtkIOFLUENTCFFModule.h" // For export macro
39
41#include "vtkNew.h" // For vtkNew
42
43VTK_ABI_NAMESPACE_BEGIN
45class vtkPoints;
46class vtkTriangle;
47class vtkTetra;
48class vtkQuad;
49class vtkHexahedron;
50class vtkPyramid;
52class vtkWedge;
53
54class VTKIOFLUENTCFF_EXPORT vtkFLUENTCFFReader : public vtkMultiBlockDataSetAlgorithm
55{
56public:
59 void PrintSelf(ostream& os, vtkIndent indent) override;
60
62
65 vtkSetMacro(FileName, std::string);
66 vtkGetMacro(FileName, std::string);
68
70
74 vtkGetMacro(NumberOfCells, vtkIdType);
76
81
86 const char* GetCellArrayName(int index);
87
89
93 int GetCellArrayStatus(const char* name);
94 void SetCellArrayStatus(const char* name, int status);
96
98
104
106 //
107 // Structures
108 //
109 struct Cell
110 {
111 int type;
112 int zone;
113 std::vector<int> faces;
115 int child;
116 std::vector<int> nodes;
117 std::vector<int> nodesOffset;
118 std::vector<int> childId;
119 };
120 struct Face
121 {
122 int type;
123 unsigned int zone;
124 std::vector<int> nodes;
125 int c0;
126 int c1;
129 int child;
134 };
135
136
137protected:
142
147 {
152 };
153
155
158 virtual bool OpenCaseFile(const std::string& filename);
159 virtual DataState OpenDataFile(const std::string& filename);
161
165 virtual void GetNumberOfCellZones();
166
170 virtual int ParseCaseFile();
171
175 virtual int GetDimension();
176
178
181 virtual void GetNodesGlobal();
182 virtual void GetCellsGlobal();
183 virtual void GetFacesGlobal();
185
189 virtual void GetNodes();
190
194 virtual void GetCells();
195
199 virtual void GetFaces();
200
206
211 virtual void GetCellOverset();
212
216 virtual void GetCellTree();
217
221 virtual void GetFaceTree();
222
227
233
237 virtual void CleanCells();
238
240
244 virtual void PopulateCellNodes();
245 virtual void PopulateCellTree();
247
249
252 virtual void PopulateTriangleCell(int i);
253 virtual void PopulateTetraCell(int i);
254 virtual void PopulateQuadCell(int i);
255 virtual void PopulateHexahedronCell(int i);
256 virtual void PopulatePyramidCell(int i);
257 virtual void PopulateWedgeCell(int i);
258 virtual void PopulatePolyhedronCell(int i);
260
264 virtual int GetData();
265
269 virtual int GetMetaData();
270
271 //
272 // Variables
273 //
275 std::string FileName;
278
279 struct vtkInternals;
280 std::unique_ptr<vtkInternals> HDFImpl;
281
289
290 std::vector<Cell> Cells;
291 std::vector<Face> Faces;
292 std::vector<int> CellZones;
293
297
298private:
299 vtkFLUENTCFFReader(const vtkFLUENTCFFReader&) = delete;
300 void operator=(const vtkFLUENTCFFReader&) = delete;
301
302 struct DataChunk
303 {
304 std::string variableName;
305 vtkIdType zoneId;
306 size_t dim;
307 std::vector<double> dataVector;
308 };
309
310 std::vector<DataChunk> DataChunks;
311 std::vector<std::string> PreReadData;
312 int NumberOfArrays = 0;
313
317 void ParseUDMData(
318 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& grid, const DataChunk& vectorDataChunk);
319};
320VTK_ABI_NAMESPACE_END
321#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:368