VTK  9.3.20240419
vtkFLUENTReader.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
26 #ifndef vtkFLUENTReader_h
27 #define vtkFLUENTReader_h
28 
29 #include "vtkIOGeometryModule.h" // For export macro
31 
32 VTK_ABI_NAMESPACE_BEGIN
34 class vtkPoints;
35 class vtkTriangle;
36 class vtkTetra;
37 class vtkQuad;
38 class vtkHexahedron;
39 class vtkPyramid;
40 class vtkWedge;
41 class vtkConvexPointSet;
42 
43 class VTKIOGEOMETRY_EXPORT vtkFLUENTReader : public vtkMultiBlockDataSetAlgorithm
44 {
45 public:
46  static vtkFLUENTReader* New();
48  void PrintSelf(ostream& os, vtkIndent indent) override;
49 
51 
57 
59 
63  vtkGetMacro(NumberOfCells, vtkIdType);
65 
70 
75  const char* GetCellArrayName(int index);
76 
78 
82  int GetCellArrayStatus(const char* name);
83  void SetCellArrayStatus(const char* name, int status);
85 
87 
93 
95 
112  void SetDataByteOrder(int);
114  //
115  // Structures
116  //
117  struct Cell;
118  struct Face;
119  struct ScalarDataChunk;
120  struct VectorDataChunk;
121  struct stdString;
122  struct intVector;
123  struct doubleVector;
124  struct stringVector;
125  struct cellVector;
126  struct faceVector;
127  struct stdMap;
128  struct scalarDataVector;
129  struct vectorDataVector;
130  struct intVectorVector;
132 
133 protected:
135  ~vtkFLUENTReader() override;
138 
140 
144  vtkSetMacro(SwapBytes, vtkTypeBool);
145  vtkTypeBool GetSwapBytes() { return this->SwapBytes; }
146  vtkBooleanMacro(SwapBytes, vtkTypeBool);
148 
149  virtual bool OpenCaseFile(const char* filename);
150  virtual bool OpenDataFile(const char* filename);
151  virtual int GetCaseChunk();
152  virtual void GetNumberOfCellZones();
153  virtual int GetCaseIndex();
154  virtual void LoadVariableNames();
155  virtual int GetDataIndex();
156  virtual int GetDataChunk();
157  virtual void GetSpeciesVariableNames();
158 
159  virtual bool ParseCaseFile();
160  virtual int GetDimension();
161  virtual void GetLittleEndianFlag();
162  virtual void GetNodesAscii();
163  virtual void GetNodesSinglePrecision();
164  virtual void GetNodesDoublePrecision();
165  virtual void GetCellsAscii();
166  virtual void GetCellsBinary();
167  virtual bool GetFacesAscii();
168  virtual void GetFacesBinary();
171  virtual void GetCellTreeAscii();
172  virtual void GetCellTreeBinary();
173  virtual void GetFaceTreeAscii();
174  virtual void GetFaceTreeBinary();
179  virtual void GetPartitionInfo() {}
180  virtual void CleanCells();
181  virtual void PopulateCellNodes();
182  virtual int GetCaseBufferInt(int ptr);
183  virtual float GetCaseBufferFloat(int ptr);
184  virtual double GetCaseBufferDouble(int ptr);
185  virtual void PopulateTriangleCell(int i);
186  virtual void PopulateTetraCell(int i);
187  virtual void PopulateQuadCell(int i);
188  virtual void PopulateHexahedronCell(int i);
189  virtual void PopulatePyramidCell(int i);
190  virtual void PopulateWedgeCell(int i);
191  virtual void PopulatePolyhedronCell(int i);
192  virtual void ParseDataFile();
193  virtual int GetDataBufferInt(int ptr);
194  virtual float GetDataBufferFloat(int ptr);
195  virtual double GetDataBufferDouble(int ptr);
196  virtual void GetData(int dataType);
197  virtual bool ParallelCheckCell(int vtkNotUsed(i)) { return true; }
198 
199  //
200  // Variables
201  //
203  char* FileName;
206 
207  istream* FluentCaseFile;
208  istream* FluentDataFile;
209  stdString* CaseBuffer;
210  stdString* DataBuffer;
211 
220 
221  cellVector* Cells;
222  faceVector* Faces;
223  stdMap* VariableNames;
224  intVector* CellZones;
225  scalarDataVector* ScalarDataChunks;
226  vectorDataVector* VectorDataChunks;
227 
228  intVectorVector* SubSectionZones;
229  intVector* SubSectionIds;
230  intVector* SubSectionSize;
231 
232  stringVector* ScalarVariableNames;
234  stringVector* VectorVariableNames;
236 
239  int DataPass;
242 
243 private:
253  void FillMultiBlockFromFaces(vtkMultiBlockDataSet* output);
254 
255  vtkFLUENTReader(const vtkFLUENTReader&) = delete;
256  void operator=(const vtkFLUENTReader&) = delete;
257 
258  bool Parsed = false;
259 };
260 VTK_ABI_NAMESPACE_END
261 #endif
a 3D cell defined by a set of convex points
Store on/off settings for data arrays, etc.
reads a dataset in Fluent file format
virtual void GetPeriodicShadowFacesBinary()
intVector * SubSectionIds
virtual void GetNumberOfCellZones()
vtkConvexPointSet * ConvexPointSet
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkHexahedron * Hexahedron
virtual void GetPartitionInfo()
faceVector * Faces
virtual void ParseDataFile()
~vtkFLUENTReader() override
virtual void PopulateWedgeCell(int i)
stdMap * VariableNames
const char * GetCellArrayName(int index)
Get the name of the cell array with the given index in the input.
virtual void GetData(int dataType)
virtual void GetSpeciesVariableNames()
virtual void GetFaceTreeAscii()
virtual void GetNonconformalGridInterfaceFaceInformationBinary()
vtkPoints * Points
virtual void GetInterfaceFaceParentsAscii()
virtual void GetNodesDoublePrecision()
virtual void CleanCells()
vtkGetFilePathMacro(FileName)
Specify the file name of the Fluent case file to read.
stdString * CaseBuffer
virtual float GetDataBufferFloat(int ptr)
virtual bool OpenDataFile(const char *filename)
virtual void GetCellsBinary()
virtual void GetPeriodicShadowFacesAscii()
stdString * DataBuffer
virtual double GetDataBufferDouble(int ptr)
istream * FluentDataFile
scalarDataVector * ScalarDataChunks
virtual void PopulateCellNodes()
virtual void GetCellTreeBinary()
void SetDataByteOrder(int)
These methods should be used instead of the SwapBytes methods.
void SetDataByteOrderToLittleEndian()
These methods should be used instead of the SwapBytes methods.
virtual int GetDataBufferInt(int ptr)
virtual void PopulateTriangleCell(int i)
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
stringVector * VectorVariableNames
vtkTypeBool GetSwapBytes()
Set/Get the byte swapping to explicitly swap the bytes of a file.
vtkTypeBool SwapBytes
const char * GetDataByteOrderAsString()
These methods should be used instead of the SwapBytes methods.
intVectorVector * SubSectionZones
virtual bool ParallelCheckCell(int vtkNotUsed(i))
virtual void PopulatePyramidCell(int i)
virtual void PopulatePolyhedronCell(int i)
virtual int GetCaseBufferInt(int ptr)
virtual void GetInterfaceFaceParentsBinary()
void SetCellArrayStatus(const char *name, int status)
Get/Set whether the cell array with the given name is to be read.
virtual void GetFaceTreeBinary()
static vtkFLUENTReader * New()
vtkPyramid * Pyramid
virtual void GetCellsAscii()
virtual void PopulateHexahedronCell(int i)
virtual bool ParseCaseFile()
int GetCellArrayStatus(const char *name)
Get/Set whether the cell array with the given name is to be read.
intVector * VectorSubSectionIds
void EnableAllCellArrays()
Turn on/off all cell arrays.
virtual double GetCaseBufferDouble(int ptr)
virtual float GetCaseBufferFloat(int ptr)
virtual void GetFacesBinary()
intVector * CellZones
virtual void GetNodesAscii()
virtual int GetCaseIndex()
virtual int GetCaseChunk()
virtual void GetNonconformalGridInterfaceFaceInformationAscii()
int GetDataByteOrder()
These methods should be used instead of the SwapBytes methods.
vtkSetFilePathMacro(FileName)
Specify the file name of the Fluent case file to read.
void SetDataByteOrderToBigEndian()
These methods should be used instead of the SwapBytes methods.
istream * FluentCaseFile
intVector * ScalarSubSectionIds
intVector * SubSectionSize
vtkDataArraySelection * CellDataArraySelection
vtkIdType NumberOfCells
virtual void PopulateQuadCell(int i)
vectorDataVector * VectorDataChunks
virtual void GetLittleEndianFlag()
virtual int GetDataChunk()
virtual void GetNodesSinglePrecision()
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkTriangle * Triangle
stringVector * ScalarVariableNames
virtual bool GetFacesAscii()
virtual void GetCellTreeAscii()
virtual bool OpenCaseFile(const char *filename)
virtual int GetDimension()
virtual void LoadVariableNames()
int GetNumberOfCellArrays()
Get the number of cell arrays available in the input.
cellVector * Cells
virtual void PopulateTetraCell(int i)
virtual int GetDataIndex()
void DisableAllCellArrays()
Turn on/off all cell arrays.
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.
Composite dataset that organizes datasets into blocks.
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
Definition: vtkTriangle.h:137
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:85
@ name
Definition: vtkX3D.h:219
@ index
Definition: vtkX3D.h:246
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315