VTK  9.6.20260328
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
39
40#ifndef vtkFLUENTReader_h
41#define vtkFLUENTReader_h
42
43#include "vtkIOGeometryModule.h" // For export macro
45#include "vtkNew.h" // For vtkNew
46
47#include <map>
48#include <unordered_set>
49
50VTK_ABI_NAMESPACE_BEGIN
52class vtkPoints;
53class vtkTriangle;
54class vtkTetra;
55class vtkQuad;
56class vtkHexahedron;
57class vtkPyramid;
58class vtkWedge;
61
62class VTKIOGEOMETRY_EXPORT vtkFLUENTReader : public vtkMultiBlockDataSetAlgorithm
63{
64public:
67 void PrintSelf(ostream& os, vtkIndent indent) override;
68
70
76
78
82 vtkGetMacro(NumberOfCells, vtkIdType);
84
86
90 vtkGetMacro(CacheData, bool);
91 vtkSetMacro(CacheData, bool);
92 vtkBooleanMacro(CacheData, bool);
94
99
104 const char* GetCellArrayName(int index);
105
107
111 int GetCellArrayStatus(const char* name);
112 void SetCellArrayStatus(const char* name, int status);
114
116
122
124
130
132
151 //
152 // Structures
153 //
154 struct Cell;
155 struct Face;
156 struct Zone;
157 struct ZoneSection;
158 struct ScalarDataChunk;
159 struct VectorDataChunk;
160 struct SubSection;
162
169
170protected:
175
177
181 vtkSetMacro(SwapBytes, vtkTypeBool);
182 vtkTypeBool GetSwapBytes() { return this->SwapBytes; }
183 vtkBooleanMacro(SwapBytes, vtkTypeBool);
185
186 virtual bool OpenCaseFile(const char* filename);
187 virtual bool OpenDataFile(const char* filename);
188 virtual int GetCaseChunk();
189 virtual int GetCaseIndex();
190 virtual void LoadVariableNames();
191 virtual int GetDataIndex();
192 virtual int GetDataChunk();
194
195 virtual int GetDimension();
196 virtual void GetLittleEndianFlag();
197 virtual void GetNodesAscii();
200 virtual void GetCellsAscii();
201 virtual void GetCellsBinary();
202 virtual bool GetFacesAscii();
203 virtual void GetFacesBinary();
206 virtual void GetCellTreeAscii();
207 virtual void GetCellTreeBinary();
208 virtual void GetFaceTreeAscii();
209 virtual void GetFaceTreeBinary();
214 virtual void GetPartitionInfo() {}
215 virtual void CleanCells();
216 virtual void PopulateCellNodes();
217 virtual void PopulateTriangleCell(size_t cellIdx);
218 virtual void PopulateTetraCell(size_t cellIdx);
219 virtual void PopulateQuadCell(size_t cellIdx);
220 virtual void PopulateHexahedronCell(size_t cellIdx);
221 virtual void PopulatePyramidCell(size_t cellIdx);
222 virtual void PopulateWedgeCell(size_t cellIdx);
223 virtual void PopulatePolyhedronCell(size_t cellIdx);
224 virtual void GetData(int dataType);
225 virtual bool ParallelCheckCell(int vtkNotUsed(i)) { return true; }
226
227private:
228 int GetCaseBufferInt(size_t ptr);
229 float GetCaseBufferFloat(size_t ptr);
230 double GetCaseBufferDouble(size_t ptr);
231 int GetDataBufferInt(size_t ptr);
232 float GetDataBufferFloat(size_t ptr);
233 double GetDataBufferDouble(size_t ptr);
234
238 bool AreCellsEnabled();
242 void DisableZones(std::unordered_set<unsigned int>& disabledZones, bool& areAllZonesDisabled);
246 bool FillMultiblock(std::unordered_set<unsigned int>& disabledZones,
247 const std::map<unsigned int, size_t>& zoneIDToBlockIdx,
248 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs);
252 void FillMultiblockData(std::unordered_set<unsigned int>& disabledZones,
253 const std::map<unsigned int, size_t>& zoneIDToBlockIdx,
254 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs);
258 void GetArraysFromSubSections();
262 void InitOutputBlocks(vtkMultiBlockDataSet* output,
263 std::map<unsigned int, size_t>& zoneIDToBlockIdx,
264 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs);
268 void ParseDataZone(int index);
275 void ParseDataZones(bool areCellsEnabled);
279 void ParseZone(int index);
286 void ParseZones(bool areCellsEnabled);
291 bool PreParseDataFile();
297 bool PreParseFluentFile();
302 bool ReadDataZoneSectionId(unsigned int& zoneSectionId);
307 bool ReadZoneSectionId(unsigned int& zoneSectionId);
312 bool ReadZoneSection(int limit);
313
317 void UpdateZoneSectionSelection();
318
330 void FillMultiBlockFromFaces(std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs,
331 const std::map<unsigned int, size_t>& zoneIDToBlockIdx,
332 std::unordered_set<unsigned int>& disabledZones);
333
334 vtkFLUENTReader(const vtkFLUENTReader&) = delete;
335 void operator=(const vtkFLUENTReader&) = delete;
336
337 //
338 // Variables
339 //
340 vtkNew<vtkDataArraySelection> ZoneSectionSelection;
341 vtkNew<vtkDataArraySelection> CellDataArraySelection;
342 char* FileName = nullptr;
343 vtkIdType NumberOfCells = 0;
344 bool CacheData = true;
345
346 istream* FluentFile = nullptr;
347 istream* FluentDataFile = nullptr;
348 std::string FluentBuffer;
349 std::string DataBuffer;
350
351 // File data cache
352 vtkNew<vtkPoints> Points;
353 std::vector<Cell> Cells;
354 std::vector<Face> Faces;
355 std::vector<ZoneSection> ZoneSections;
356 std::map<size_t, std::string> VariableNames;
357 std::vector<ScalarDataChunk> ScalarDataChunks;
358 std::vector<VectorDataChunk> VectorDataChunks;
359 std::vector<SubSection> SubSections;
360
361 std::vector<std::string> ScalarVariableNames;
362 std::vector<int> ScalarSubSectionIds;
363 std::vector<std::string> VectorVariableNames;
364 std::vector<int> VectorSubSectionIds;
365
366 vtkTypeBool SwapBytes;
367 int GridDimension = 0;
368 int NumberOfScalars = 0;
369 int NumberOfVectors = 0;
370
371 std::vector<Zone> Zones;
372 std::vector<Zone> DataZones;
373 std::vector<ZoneSection> CurrentZoneSections;
374
375 bool IsFilePreParsed = false;
376};
377
378VTK_ABI_NAMESPACE_END
379#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 PopulateWedgeCell(size_t cellIdx)
virtual void GetPeriodicShadowFacesBinary()
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual void GetPartitionInfo()
virtual void PopulateQuadCell(size_t cellIdx)
~vtkFLUENTReader() override
virtual void GetData(int dataType)
virtual void GetSpeciesVariableNames()
virtual void GetFaceTreeAscii()
virtual void PopulatePolyhedronCell(size_t cellIdx)
virtual void GetNonconformalGridInterfaceFaceInformationBinary()
virtual void GetInterfaceFaceParentsAscii()
virtual void PopulatePyramidCell(size_t cellIdx)
virtual void GetNodesDoublePrecision()
virtual void CleanCells()
vtkGetFilePathMacro(FileName)
Specify the file name of the Fluent file to read.
virtual bool OpenDataFile(const char *filename)
virtual void GetCellsBinary()
virtual void GetPeriodicShadowFacesAscii()
virtual void PopulateCellNodes()
virtual void GetCellTreeBinary()
const char * GetCellArrayName(int index)
Get the name of the cell array with the given index in the input.
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.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkTypeBool GetSwapBytes()
Set/Get the byte swapping to explicitly swap the bytes of a file.
virtual void PopulateHexahedronCell(size_t cellIdx)
virtual void GetInterfaceFaceParentsBinary()
vtkDataArraySelection * GetZoneSectionSelection()
Zone section selection, to determine which zone sections are loaded.
void SetCellArrayStatus(const char *name, int status)
Get/Set whether the cell array with the given name is to be read.
virtual void GetFaceTreeBinary()
virtual void GetCellsAscii()
int GetCellArrayStatus(const char *name)
Get/Set whether the cell array with the given name is to be read.
virtual void PopulateTetraCell(size_t cellIdx)
void EnableAllCellArrays()
Turn on/off all cell arrays.
virtual void GetFacesBinary()
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 file to read.
void SetDataByteOrderToBigEndian()
These methods should be used instead of the SwapBytes methods.
virtual bool ParallelCheckCell(int i)
const char * GetDataByteOrderAsString()
These methods should be used instead of the SwapBytes methods.
virtual void PopulateTriangleCell(size_t cellIdx)
virtual void GetLittleEndianFlag()
virtual int GetDataChunk()
virtual void GetNodesSinglePrecision()
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual bool GetFacesAscii()
virtual void GetCellTreeAscii()
virtual bool OpenCaseFile(const char *filename)
virtual int GetDimension()
virtual void LoadVariableNames()
vtkMTimeType GetMTime() override
Get the last modified time of this filter.
int GetNumberOfCellArrays()
Get the number of cell arrays available in the input.
static vtkFLUENTReader * New()
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.
Composite dataset that organizes datasets into blocks.
Allocate and hold a VTK object.
Definition vtkNew.h:167
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
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
int vtkTypeBool
Definition vtkABI.h:64
@ Cells
A cell specified by degrees of freedom held in arrays.
int vtkIdType
Definition vtkType.h:363
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318