VTK  9.4.20250310
vtkNek5000Reader.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
15#ifndef vtkNek5000Reader_h
16#define vtkNek5000Reader_h
17
18#include "vtkIOParallelModule.h" // For export macro
20
21VTK_ABI_NAMESPACE_BEGIN
22
23class vtkPoints;
25
26class VTKIOPARALLEL_EXPORT vtkNek5000Reader : public vtkUnstructuredGridAlgorithm
27{
28public:
31 void PrintSelf(ostream& os, vtkIndent indent) override;
32
34
35 vtkSetStringMacro(FileName);
36 vtkGetStringMacro(FileName);
37
38 vtkSetStringMacro(DataFileName);
39 vtkGetStringMacro(DataFileName);
40
41 vtkGetMacro(NumberOfTimeSteps, int);
43
46 vtkGetVector2Macro(TimeStepRange, int);
47 vtkSetVector2Macro(TimeStepRange, int);
49
54
59 const char* GetPointArrayName(int index);
63 vtkSetMacro(CleanGrid, int);
64 vtkGetMacro(CleanGrid, int);
65 vtkBooleanMacro(CleanGrid, int);
66
68
71 vtkSetMacro(SpectralElementIds, int);
72 vtkGetMacro(SpectralElementIds, int);
73 vtkBooleanMacro(SpectralElementIds, int);
75
77
80 bool GetPointArrayStatus(const char* name);
81 bool GetPointArrayStatus(int index);
82 void SetPointArrayStatus(const char* name, int status);
84
86
92
96 size_t GetVariableNamesFromData(char* varTags);
97
98 int CanReadFile(const char* fname);
99
100protected:
103
104 char* FileName;
106 // int ElementResolution;
107 // int BoundaryResolution;
109 // int my_patch_id;
110
111 int num_vars; // all vars including Pressure, Velocity, Velocity Magnitude and Temperature
112 char** var_names;
113 float** dataArray;
115
117
120
121 // Tri* T;
122 class nek5KList;
123 class nek5KObject;
124 nek5KList* myList;
125 nek5KObject* curObj;
129
131
132 std::string datafile_format;
136
137 // void setActive(); // set my_patch_id as the active one
138 // static int getNextPatchID(){return(next_patch_id++);}
139
141
142 // update which fields from the data should be used, based on GUI
145 void readData(char* dfName);
146 // copy the data from nek5000 to pv
147 void updateVtuData(vtkUnstructuredGrid* pv_ugrid); //, vtkUnstructuredGrid* pv_boundary_ugrid);
149 void addSpectralElementId(int nelements);
151 // void interpolateAndCopyContinuumData(vtkUnstructuredGrid* pv_ugrid, double **data_array, int
152 // interp_res, int num_verts);
154 // void interpolateAndCopyBoundaryPoints(int alloc_res, int interp_res, vtkPoints*
155 // boundary_points); void interpolateAndCopyBoundaryData(int alloc_res, int num_verts, int
156 // interp_res); void addCellsToBoundaryMesh(int * boundary_index, int qa); void
157 // generateBoundaryConnectivity(int * boundary_index, int res);
158 // see if the current object is missing data that was requested
160 // see if the current object matches the request
162 // see if the current object has extra data than was requested
164
166 // vtkUnstructuredGrid* Boundary_UGrid;
167 bool CALC_GEOM_FLAG; // true = need to calculate continuum geometry; false = geom is up to date
168 // bool CALC_BOUNDARY_GEOM_FLAG; // true = need to calculate boundary geometry; false = boundary
169 // geom is up to date bool HAVE_BOUNDARY_GEOM_FLAG; // true = we have boundary geometry; false =
170 // geom has not been read yet
171
172 bool READ_GEOM_FLAG; // true = need continuum geom from disk
173 // bool READ_BOUNDARY_GEOM_FLAG; // true = need boundary geom from disk
176
178 // int TimeStep;
180 int blockDims[3];
190 double TimeValue;
191 int TimeStepRange[2];
193
194 std::vector<double> TimeSteps;
195 // int UseProjection;
196 // int ExtractBoundary;
197 // int DynamicMesh;
198 // double DynamicMeshScale;
199
200 // Time query function. Called by ExecuteInformation().
201 // Fills the TimestepValues array.
203
205
206 // Description:
207 // This is called by the superclass.
208 // This is the method you should override.
210
211private:
212 vtkNek5000Reader(const vtkNek5000Reader&) = delete; // Not implemented.
213 void operator=(const vtkNek5000Reader&) = delete; // Not implemented.
214
215 int SpectralElementIds;
216 int CleanGrid;
217};
218
219VTK_ABI_NAMESPACE_END
220#endif
Store on/off settings for data arrays, etc.
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Reads Nek5000 format data files.
std::vector< double > TimeSteps
void copyContinuumPoints(vtkPoints *points)
void updateVariableStatus()
nek5KObject * curObj
static vtkNek5000Reader * New()
void updateVtuData(vtkUnstructuredGrid *pv_ugrid)
void SetPointArrayStatus(const char *name, int status)
Get/Set whether the point array with the given name or index is to be read.
vtkUnstructuredGrid * UGrid
void DisableAllPointArrays()
Turn on/off all point arrays.
bool GetAllTimesAndVariableNames(vtkInformationVector *)
int CanReadFile(const char *fname)
void addCellsToContinuumMesh()
bool objectMatchesRequest()
bool isObjectMissingData()
vtkMTimeType GetMTime() override
Return this object's modified time.
vtkDataArraySelection * PointDataArraySelection
const char * GetPointArrayName(int index)
Get the name of the point array with the given index in the input.
void addSpectralElementId(int nelements)
size_t GetVariableNamesFromData(char *varTags)
Get the names of variables stored in the data.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
bool objectHasExtraData()
std::string datafile_format
bool GetPointArrayStatus(const char *name)
Get/Set whether the point array with the given name or index is to be read.
void partitionAndReadMesh()
void copyContinuumData(vtkUnstructuredGrid *pv_ugrid)
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
int GetNumberOfPointArrays()
Get the number of point arrays available in the input.
void EnableAllPointArrays()
Turn on/off all point arrays.
~vtkNek5000Reader() override
void readData(char *dfName)
bool GetPointArrayStatus(int index)
Get/Set whether the point array with the given name or index is to be read.
represent and manipulate 3D points
Definition vtkPoints.h:139
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287