VTK  9.3.1
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:
55  static vtkFLUENTCFFReader* New();
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 
78  int GetNumberOfCellArrays();
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 
99  void DisableAllCellArrays();
100  void EnableAllCellArrays();
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:
149  ~vtkFLUENTCFFReader() override;
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 
242  virtual void GetNonconformalGridInterfaceFaceInformation();
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
std::vector< std::string > PreReadVectorData
vtkNew< vtkPoints > Points
Store vtkAlgorithm input/output information.
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:35
std::vector< int > childId
std::unique_ptr< vtkInternals > HDFImpl
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:27
int vtkIdType
Definition: vtkType.h:315
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
std::vector< Cell > Cells
static vtkMultiBlockDataSetAlgorithm * New()
int vtkTypeBool
Definition: vtkABI.h:64
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:33
std::vector< VectorDataChunk > VectorDataChunks
vtkNew< vtkHexahedron > Hexahedron
reads a dataset in Fluent CFF file format
a simple class to control print indentation
Definition: vtkIndent.h:28
vtkNew< vtkTriangle > Triangle
Store on/off settings for data arrays, etc.
vtkNew< vtkDataArraySelection > CellDataArraySelection
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:33
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSmartPointer< vtkDataArray > GetData(const Ioss::GroupingEntity *entity, const std::string &fieldname, Ioss::Transform *transform=nullptr, Cache *cache=nullptr, const std::string &cachekey=std::string())
Returns a VTK array for a given field (fieldname) on the chosen block (or set) entity.
vtkNew< vtkPyramid > Pyramid
std::vector< Face > Faces
std::vector< ScalarDataChunk > ScalarDataChunks
vtkNew< vtkTetra > Tetra
std::vector< int > CellZones
a cell that represents a triangle
Definition: vtkTriangle.h:27
Store zero or more vtkInformation instances.
std::vector< std::string > PreReadScalarData
vtkNew< vtkQuad > Quad
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:35
represent and manipulate 3D points
Definition: vtkPoints.h:28
vtkNew< vtkWedge > Wedge
virtual int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.