1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200
|
// ************************************************************************* //
// avtFluentFileFormat.h //
// ************************************************************************* //
#ifndef AVT_Fluent_FILE_FORMAT_H
#define AVT_Fluent_FILE_FORMAT_H
#include <avtSTMDFileFormat.h>
#include <vector>
#include <string>
#include <map>
// Fluent plugin
#include <visitstream.h>
#include <stdlib.h>
#include <sstream>
#include "vtkPoints.h"
#include "vtkTriangle.h"
#include "vtkTetra.h"
#include "vtkQuad.h"
#include "vtkHexahedron.h"
#include "vtkPyramid.h"
#include "vtkWedge.h"
#include "vtkConvexPointSet.h"
#include "vtkUnstructuredGrid.h"
#include "vtkDoubleArray.h"
using namespace std;
// ****************************************************************************
// Class: avtFluentFileFormat
//
// Purpose:
// Reads in Fluent files as a plugin to VisIt.
//
// Programmer: bdotson -- generated by xml2avt
// Creation: Fri Jun 30 15:02:33 PST 2006
//
// Modifications:
//
// Hank Childs, Fri Sep 8 14:30:04 PDT 2006
// Added methods HasInvariantMetaData and HasInvariantSIL on the advice
// of Terry Jordan.
//
// ****************************************************************************
class avtFluentFileFormat : public avtSTMDFileFormat
{
public:
avtFluentFileFormat(const char *);
virtual ~avtFluentFileFormat() {;};
virtual const char *GetType(void) { return "Fluent"; };
virtual void FreeUpResources(void);
virtual vtkDataSet *GetMesh(int, const char *);
virtual vtkDataArray *GetVar(int, const char *);
virtual vtkDataArray *GetVectorVar(int, const char *);
//Number of variables is dynamic
virtual bool HasInvariantMetaData(void) const
{ return false; };
//Number of Domains is dynamic
virtual bool HasInvariantSIL(void) const
{ return false; };
protected:
virtual void PopulateDatabaseMetaData(avtDatabaseMetaData *);
int OpenCaseFile(const char *filename);
int OpenDataFile(const char *filename);
int GetCaseChunk ();
void GetNumberOfCellZones();
int GetCaseIndex();
void LoadVariableNames();
int GetDataIndex();
int GetDataChunk();
void ParseCaseFile();
int GetDimension();
void GetLittleEndianFlag();
void GetNodesAscii();
void GetNodesSinglePrecision();
void GetNodesDoublePrecision();
void GetCellsAscii();
void GetCellsBinary();
void GetFacesAscii();
void GetFacesBinary();
void GetPeriodicShadowFacesAscii();
void GetPeriodicShadowFacesBinary();
void GetCellTreeAscii();
void GetCellTreeBinary();
void GetFaceTreeAscii();
void GetFaceTreeBinary();
void GetInterfaceFaceParentsAscii();
void GetInterfaceFaceParentsBinary();
void GetNonconformalGridInterfaceFaceInformationAscii();
void GetNonconformalGridInterfaceFaceInformationBinary();
void CleanCells();
void PopulateCellNodes();
int GetCaseBufferInt(int ptr);
float GetCaseBufferFloat(int ptr);
double GetCaseBufferDouble(int ptr);
void PopulateTriangleCell(int i);
void PopulateTetraCell(int i);
void PopulateQuadCell(int i);
void PopulateHexahedronCell(int i);
void PopulatePyramidCell(int i);
void PopulateWedgeCell(int i);
void PopulatePolyhedronCell(int i);
void ParseDataFile();
int GetDataBufferInt(int ptr);
float GetDataBufferFloat(int ptr);
double GetDataBufferDouble(int ptr);
void GetData(int dataType);
//
// Variables
//
ifstream FluentCaseFile;
ifstream FluentDataFile;
string CaseBuffer;
string DataBuffer;
vtkPoints *Points;
vtkTriangle *Triangle;
vtkTetra *Tetra;
vtkQuad *Quad;
vtkHexahedron *Hexahedron;
vtkPyramid *Pyramid;
vtkWedge *Wedge;
vtkConvexPointSet *ConvexPointSet;
struct Cell {
int type;
int zone;
vector<int> faces;
int parent;
int child;
vector<int> nodes;
};
struct Face {
int type;
int zone;
vector<int> nodes;
int c0;
int c1;
int periodicShadow;
int parent;
int child;
int interfaceFaceParent;
int interfaceFaceChild;
int ncgParent;
int ncgChild;
};
struct ScalarDataChunk {
int subsectionId;
int zoneId;
vector<double> scalarData;
};
struct VectorDataChunk {
int subsectionId;
int zoneId;
vector<double> iComponentData;
vector<double> jComponentData;
vector<double> kComponentData;
};
vector< Cell > Cells;
vector< Face > Faces;
map< int, string > VariableNames;
vector< int > CellZones;
vector< ScalarDataChunk > ScalarDataChunks;
vector< VectorDataChunk > VectorDataChunks;
vector< vector<int> > SubSectionZones;
vector< int > SubSectionIds;
vector< int > SubSectionSize;
vector< string > ScalarVariableNames;
vector< int > ScalarSubSectionIds;
vector< string > VectorVariableNames;
vector< int > VectorSubSectionIds;
int LittleEndianFlag;
int GridDimension;
int DataPass;
int NumberOfScalars;
int NumberOfVectors;
};
#endif
|