VTK  9.3.0
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 
42 VTK_ABI_NAMESPACE_BEGIN
44 class vtkPoints;
45 class vtkTriangle;
46 class vtkTetra;
47 class vtkQuad;
48 class vtkHexahedron;
49 class vtkPyramid;
50 class vtkWedge;
51 
52 class VTKIOFLUENTCFF_EXPORT vtkFLUENTCFFReader : public vtkMultiBlockDataSetAlgorithm
53 {
54 public:
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;
112  int parent;
113  int child;
114  std::vector<int> nodes;
115  std::vector<int> childId;
116  };
117  struct Face
118  {
119  int type;
120  unsigned int zone;
121  std::vector<int> nodes;
122  int c0;
123  int c1;
125  int parent;
126  int child;
130  int ncgChild;
131  };
133  {
136  std::vector<double> scalarData;
137  };
139  {
142  size_t dim;
143  std::vector<double> vectorData;
144  };
146 
147 protected:
152 
157  {
158  NOT_LOADED = 0,
159  AVAILABLE = 1,
160  LOADED = 2,
161  ERROR = 3
162  };
163 
165 
168  virtual bool OpenCaseFile(const std::string& filename);
169  virtual DataState OpenDataFile(const std::string& filename);
171 
175  virtual void GetNumberOfCellZones();
176 
180  virtual int ParseCaseFile();
181 
185  virtual int GetDimension();
186 
188 
191  virtual void GetNodesGlobal();
192  virtual void GetCellsGlobal();
193  virtual void GetFacesGlobal();
195 
199  virtual void GetNodes();
200 
204  virtual void GetCells();
205 
209  virtual void GetFaces();
210 
215  virtual void GetPeriodicShadowFaces();
216 
221  virtual void GetCellOverset();
222 
226  virtual void GetCellTree();
227 
231  virtual void GetFaceTree();
232 
236  virtual void GetInterfaceFaceParents();
237 
243 
247  virtual void CleanCells();
248 
250 
254  virtual void PopulateCellNodes();
255  virtual void PopulateCellTree();
257 
259 
262  virtual void PopulateTriangleCell(int i);
263  virtual void PopulateTetraCell(int i);
264  virtual void PopulateQuadCell(int i);
265  virtual void PopulateHexahedronCell(int i);
266  virtual void PopulatePyramidCell(int i);
267  virtual void PopulateWedgeCell(int i);
268  virtual void PopulatePolyhedronCell(int i);
270 
274  virtual int GetData();
275 
279  virtual int GetMetaData();
280 
281  //
282  // Variables
283  //
286  vtkIdType NumberOfCells = 0;
287  int NumberOfCellArrays = 0;
288 
289  struct vtkInternals;
290  std::unique_ptr<vtkInternals> HDFImpl;
291 
299 
300  std::vector<Cell> Cells;
301  std::vector<Face> Faces;
302  std::vector<int> CellZones;
303  std::vector<ScalarDataChunk> ScalarDataChunks;
304  std::vector<VectorDataChunk> VectorDataChunks;
305  std::vector<std::string> PreReadScalarData;
306  std::vector<std::string> PreReadVectorData;
307 
308  vtkTypeBool SwapBytes = 0;
309  int GridDimension = 0;
310  DataState FileState = DataState::NOT_LOADED;
311  int NumberOfScalars = 0;
312  int NumberOfVectors = 0;
313 
314 private:
315  vtkFLUENTCFFReader(const vtkFLUENTCFFReader&) = delete;
316  void operator=(const vtkFLUENTCFFReader&) = delete;
317 };
318 VTK_ABI_NAMESPACE_END
319 #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
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
const char * GetCellArrayName(int index)
Get the name of the cell array with the given index in the input.
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)
static vtkFLUENTCFFReader * New()
virtual void PopulateWedgeCell(int i)
virtual void PopulateTriangleCell(int i)
Reconstruct VTK cell topology from Fluent format.
std::vector< Cell > Cells
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
Definition: vtkHexahedron.h:34
a simple class to control print indentation
Definition: vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
represent and manipulate 3D points
Definition: vtkPoints.h:29
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:36
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:28
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:34
a cell that represents a triangle
Definition: vtkTriangle.h:28
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:36
@ name
Definition: vtkX3D.h:219
@ index
Definition: vtkX3D.h:246
@ string
Definition: vtkX3D.h:490
std::vector< int > childId
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315