VTK  9.4.20250114
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
40#ifndef vtkFLUENTReader_h
41#define vtkFLUENTReader_h
42
43#include "vtkDeprecation.h" // For deprecation macro
44#include "vtkIOGeometryModule.h" // For export macro
46#include "vtkNew.h" // For vtkNew
47
48#include <set>
49#include <unordered_map>
50
51VTK_ABI_NAMESPACE_BEGIN
53class vtkPoints;
54class vtkTriangle;
55class vtkTetra;
56class vtkQuad;
57class vtkHexahedron;
58class vtkPyramid;
59class vtkWedge;
62
63class VTKIOGEOMETRY_EXPORT vtkFLUENTReader : public vtkMultiBlockDataSetAlgorithm
64{
65public:
68 void PrintSelf(ostream& os, vtkIndent indent) override;
69
71
77
79
83 vtkGetMacro(NumberOfCells, vtkIdType);
85
87
91 vtkGetMacro(CacheData, bool);
92 vtkSetMacro(CacheData, bool);
93 vtkBooleanMacro(CacheData, bool);
95
100
105 const char* GetCellArrayName(int index);
106
108
112 int GetCellArrayStatus(const char* name);
113 void SetCellArrayStatus(const char* name, int status);
115
117
123
125
131
133
152 //
153 // Structures
154 //
155 struct Cell;
156 struct Face;
157 struct Zone;
158 struct ZoneSection;
159 struct ScalarDataChunk;
160 struct VectorDataChunk;
161 struct SubSection;
163
170
171protected:
176
178
182 vtkSetMacro(SwapBytes, vtkTypeBool);
183 vtkTypeBool GetSwapBytes() { return this->SwapBytes; }
184 vtkBooleanMacro(SwapBytes, vtkTypeBool);
186
187 virtual bool OpenCaseFile(const char* filename);
188 virtual bool OpenDataFile(const char* filename);
189 virtual int GetCaseChunk();
190 virtual int GetCaseIndex();
191 virtual void LoadVariableNames();
192 virtual int GetDataIndex();
193 virtual int GetDataChunk();
195
196 virtual int GetDimension();
197 virtual void GetLittleEndianFlag();
198 virtual void GetNodesAscii();
201 virtual void GetCellsAscii();
202 virtual void GetCellsBinary();
203 virtual bool GetFacesAscii();
204 virtual void GetFacesBinary();
207 virtual void GetCellTreeAscii();
208 virtual void GetCellTreeBinary();
209 virtual void GetFaceTreeAscii();
210 virtual void GetFaceTreeBinary();
215 virtual void GetPartitionInfo() {}
216 virtual void CleanCells();
217 virtual void PopulateCellNodes();
218 virtual int GetCaseBufferInt(int ptr);
219 virtual float GetCaseBufferFloat(int ptr);
220 virtual double GetCaseBufferDouble(int ptr);
221 virtual void PopulateTriangleCell(size_t cellIdx);
222 virtual void PopulateTetraCell(size_t cellIdx);
223 virtual void PopulateQuadCell(size_t cellIdx);
224 virtual void PopulateHexahedronCell(size_t cellIdx);
225 virtual void PopulatePyramidCell(size_t cellIdx);
226 virtual void PopulateWedgeCell(size_t cellIdx);
227 virtual void PopulatePolyhedronCell(size_t cellIdx);
228 virtual int GetDataBufferInt(int ptr);
229 virtual float GetDataBufferFloat(int ptr);
230 virtual double GetDataBufferDouble(int ptr);
231 virtual void GetData(int dataType);
232 virtual bool ParallelCheckCell(int vtkNotUsed(i)) { return true; }
233
235 "ReadZone is deprecated. It was an internal method an should not be used.")
236 virtual void ReadZone();
238 "ParseCaseFile is deprecated. It was an internal method an should not be used.")
239 virtual bool ParseCaseFile();
241 "ParseDataFile is deprecated. It was an internal method an should not be used.")
242 virtual void ParseDataFile();
243
244private:
248 bool AreCellsEnabled();
252 void DisableCellsAndFaces(std::vector<unsigned int>& disabledZones);
256 void DisableZones(std::vector<unsigned int>& disabledZones, bool& areAllZonesDisabled);
260 bool FillMultiblock(std::vector<unsigned int>& disabledZones,
261 std::vector<size_t>& zoneIDToBlockIdx,
262 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs);
266 void FillMultiblockData(std::vector<unsigned int>& disabledZones,
267 std::vector<size_t>& zoneIDToBlockIdx,
268 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs);
272 void GetArraysFromSubSections();
276 void InitOutputBlocks(vtkMultiBlockDataSet* output, std::vector<size_t>& zoneIDToBlockIdx,
277 std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs);
281 void ParseDataZone(int index);
288 void ParseDataZones(bool areCellsEnabled);
292 void ParseZone(int index);
299 void ParseZones(bool areCellsEnabled);
304 bool PreParseDataFile();
310 bool PreParseFluentFile();
315 bool ReadDataZoneSectionId(unsigned int& zoneSectionId);
320 bool ReadZoneSectionId(unsigned int& zoneSectionId);
325 bool ReadZoneSection(int limit);
326
330 void UpdateZoneSectionSelection();
331
343 void FillMultiBlockFromFaces(std::vector<vtkSmartPointer<vtkUnstructuredGrid>>& blockUGs,
344 const std::vector<size_t>& zoneIDToBlockIdx, std::vector<unsigned int> disabledZones);
345
346 vtkFLUENTReader(const vtkFLUENTReader&) = delete;
347 void operator=(const vtkFLUENTReader&) = delete;
348
349 //
350 // Variables
351 //
352 vtkNew<vtkDataArraySelection> ZoneSectionSelection;
353 vtkNew<vtkDataArraySelection> CellDataArraySelection;
354 char* FileName = nullptr;
355 vtkIdType NumberOfCells = 0;
356 bool CacheData = true;
357
358 istream* FluentFile = nullptr;
359 istream* FluentDataFile = nullptr;
360 std::string FluentBuffer;
361 std::string DataBuffer;
362
363 // File data cache
364 vtkNew<vtkPoints> Points;
365 std::vector<Cell> Cells;
366 std::vector<Face> Faces;
367 std::vector<ZoneSection> ZoneSections;
368 std::map<size_t, std::string> VariableNames;
369 std::vector<ScalarDataChunk> ScalarDataChunks;
370 std::vector<VectorDataChunk> VectorDataChunks;
371 std::vector<SubSection> SubSections;
372
373 std::vector<std::string> ScalarVariableNames;
374 std::vector<int> ScalarSubSectionIds;
375 std::vector<std::string> VectorVariableNames;
376 std::vector<int> VectorSubSectionIds;
377
378 vtkTypeBool SwapBytes;
379 int GridDimension = 0;
380 int NumberOfScalars = 0;
381 int NumberOfVectors = 0;
382
383 std::vector<Zone> Zones;
384 std::vector<Zone> DataZones;
385 std::vector<Cell> CurrentCells;
386 std::vector<Face> CurrentFaces;
387 std::vector<ZoneSection> CurrentZoneSections;
388
389 bool IsFilePreParsed = false;
390};
391
392VTK_ABI_NAMESPACE_END
393#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 float GetDataBufferFloat(int ptr)
virtual bool OpenDataFile(const char *filename)
virtual void GetCellsBinary()
virtual void GetPeriodicShadowFacesAscii()
virtual double GetDataBufferDouble(int ptr)
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.
virtual int GetDataBufferInt(int ptr)
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 bool ParallelCheckCell(int vtkNotUsed(i))
virtual int GetCaseBufferInt(int ptr)
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 double GetCaseBufferDouble(int ptr)
virtual float GetCaseBufferFloat(int ptr)
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.
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.
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
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:139
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.
#define VTK_DEPRECATED_IN_9_5_0(reason)
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270