VTK  9.6.20260417
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 <map> // std::map
42#include <memory> // std::unique_ptr
43
44#include "vtkIOFLUENTCFFModule.h" // For export macro
45
47#include "vtkNew.h" // For vtkNew
48
49VTK_ABI_NAMESPACE_BEGIN
52class vtkPoints;
53class vtkTriangle;
54class vtkTetra;
55class vtkQuad;
56class vtkHexahedron;
57class vtkPyramid;
59class vtkWedge;
60
61class VTKIOFLUENTCFF_EXPORT vtkFLUENTCFFReader : public vtkMultiBlockDataSetAlgorithm
62{
63public:
66 void PrintSelf(ostream& os, vtkIndent indent) override;
67
69
72 vtkSetMacro(FileName, std::string);
73 vtkGetMacro(FileName, std::string);
75
77
82 vtkSetMacro(RenameArrays, bool);
83 vtkGetMacro(RenameArrays, bool);
85
87
91 vtkGetMacro(NumberOfCells, vtkIdType);
93
98
103 const char* GetCellArrayName(int index);
104
106
110 int GetCellArrayStatus(const char* name);
111 void SetCellArrayStatus(const char* name, int status);
113
115
121
126
131 const char* GetFaceArrayName(int index);
132
134
138 int GetFaceArrayStatus(const char* name);
139 void SetFaceArrayStatus(const char* name, int status);
141
143
149
151
155 vtkSetMacro(ReadFaces, bool);
156 vtkGetMacro(ReadFaces, bool);
157 vtkBooleanMacro(ReadFaces, bool);
159
164
166 //
167 // Structures
168 //
169 struct Cell
170 {
171 int type;
172 int zone;
173 std::vector<int> faces;
175 int child;
176 std::vector<int> nodes;
177 std::vector<int> nodesOffset;
178 std::vector<int> childId;
179 };
180 struct Face
181 {
182 int type;
183 unsigned int zone;
184 std::vector<int> nodes;
185 int c0;
186 int c1;
189 int child;
194 };
195
196
197protected:
202
207 {
212 };
213
215
218 virtual bool OpenCaseFile(const std::string& filename);
219 virtual DataState OpenDataFile(const std::string& filename);
221
225 virtual void GetNumberOfCellZones();
226
230 virtual int ParseCaseFile();
231
235 virtual int GetDimension();
236
238
241 virtual void GetNodesGlobal();
242 virtual void GetCellsGlobal();
243 virtual void GetFacesGlobal();
245
249 virtual void GetNodes();
250
254 virtual void GetCells();
255
259 virtual void GetFaces();
260
266
271 virtual void GetCellOverset();
272
276 virtual void GetCellTree();
277
281 virtual void GetFaceTree();
282
287
293
297 virtual void CleanCells();
298
300
304 virtual void PopulateCellNodes();
305 virtual void PopulateCellTree();
307
309
312 virtual void PopulateTriangleCell(int i);
313 virtual void PopulateTetraCell(int i);
314 virtual void PopulateQuadCell(int i);
315 virtual void PopulateHexahedronCell(int i);
316 virtual void PopulatePyramidCell(int i);
317 virtual void PopulateWedgeCell(int i);
318 virtual void PopulatePolyhedronCell(int i);
320
324 virtual int GetData();
325
329 virtual int GetMetaData();
330
331 //
332 // Variables
333 //
335 std::string FileName;
336 bool RenameArrays = false;
339
340 struct vtkInternals;
341 std::unique_ptr<vtkInternals> HDFImpl;
342
350
351 std::vector<Cell> Cells;
352 std::vector<Face> Faces;
353 std::vector<int> CellZones;
354
358
359private:
360 vtkFLUENTCFFReader(const vtkFLUENTCFFReader&) = delete;
361 void operator=(const vtkFLUENTCFFReader&) = delete;
362
363 struct FaceZone
364 {
365 std::string name;
366 int minId;
367 int maxId;
368 int zoneType;
369 };
370
371 struct DataChunk
372 {
373 std::string variableName;
374 vtkIdType zoneId;
375 size_t dim;
376 std::vector<double> dataVector;
377 };
378
379 typedef long long hid_t;
380 void ReadFieldMetaData(hid_t parentGroup, int iphase, std::vector<std::string>& container);
381
385 void ParseUDMData(
386 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& grid, const DataChunk& vectorDataChunk);
387
391 void CreateFaces(vtkMultiBlockDataSet* output);
392
398 void FillDataArrayForFaceZone(const FaceZone& zone, vtkUnstructuredGrid* faceGrid,
399 const std::map<unsigned int, vtkIdType>& faceIdToLocalIndex);
400
401 std::vector<DataChunk> DataChunks;
402 std::vector<std::string> PreReadData;
403 std::vector<std::string> PreReadFaceData;
404 int NumberOfArrays = 0;
405
406 std::map<int, FaceZone> FaceZonesById;
407
408 vtkNew<vtkDataArraySelection> FaceDataArraySelection;
409 bool ReadFaces = false;
410};
411VTK_ABI_NAMESPACE_END
412#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)
int GetFaceArrayStatus(const char *name)
Get/Set whether the face array with the given name is to be read.
void DisableAllCellArrays()
Turn on/off all cell arrays.
virtual void PopulatePolyhedronCell(int i)
const char * GetFaceArrayName(int index)
Get the name of the face array with the given index in the input.
int GetNumberOfFaceArrays()
Get the number of face arrays available in the input.
std::unique_ptr< vtkInternals > HDFImpl
std::vector< Face > Faces
void DisableAllFaceArrays()
Turn on/off all face arrays.
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.
vtkMTimeType GetMTime() override
Overridden to take into account mtimes for vtkDataArraySelection instances.
virtual void GetFacesGlobal()
void SetFaceArrayStatus(const char *name, int status)
~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.
void EnableAllFaceArrays()
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.
Composite dataset that organizes datasets into blocks.
Allocate and hold a VTK object.
Definition vtkNew.h:168
represent and manipulate 3D points
Definition vtkPoints.h:140
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
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
int64_t hid_t
int vtkIdType
Definition vtkType.h:363
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318