VTK  9.0.1
vtkWindBladeReader.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkWindBladeReader.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
30 #ifndef vtkWindBladeReader_h
31 #define vtkWindBladeReader_h
32 
33 #include "vtkIOGeometryModule.h" // For export macro
35 
37 class vtkCallbackCommand;
38 class vtkStringArray;
39 class vtkFloatArray;
40 class vtkIntArray;
41 class vtkPoints;
42 class vtkStructuredGrid;
44 class vtkMultiBlockDataSetAglorithm;
46 class WindBladeReaderInternal;
47 
48 class VTKIOGEOMETRY_EXPORT vtkWindBladeReader : public vtkStructuredGridAlgorithm
49 {
50 public:
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
55  vtkSetStringMacro(Filename);
56  vtkGetStringMacro(Filename);
57 
58  vtkSetVector6Macro(WholeExtent, int);
59  vtkGetVector6Macro(WholeExtent, int);
60 
61  vtkSetVector6Macro(SubExtent, int);
62  vtkGetVector6Macro(SubExtent, int);
63 
67  vtkStructuredGrid* GetFieldOutput(); // Output port 0
68  vtkUnstructuredGrid* GetBladeOutput(); // Output port 1
69  vtkStructuredGrid* GetGroundOutput(); // Output port 2
70 
72 
78  const char* GetPointArrayName(int index);
80 
81  int GetPointArrayStatus(const char* name);
82  void SetPointArrayStatus(const char* name, int status);
83 
86 
87 protected:
88  static float DRY_AIR_CONSTANT;
89  static int NUM_PART_SIDES; // Blade parts rhombus
90  static const int NUM_BASE_SIDES; // Base pyramid
91  static const int LINE_SIZE;
92  static int DIMENSION;
93  static int BYTES_PER_DATA;
94  static int SCALAR;
95  static int VECTOR;
96  static int FLOAT;
97  static int INTEGER;
98 
101 
102  char* Filename; // Base file name
103 
104  // Extent information
105  vtkIdType NumberOfTuples; // Number of tuples in subextent
106 
107  // Field
108  int WholeExtent[6]; // Extents of entire grid
109  int SubExtent[6]; // Processor grid extent
110  int UpdateExtent[6];
111  int Dimension[3]; // Size of entire grid
112  int SubDimension[3]; // Size of processor grid
113 
114  // Ground
115  int GExtent[6]; // Extents of ground grid
116  int GSubExtent[6]; // Processor grid extent
117  int GDimension[3]; // Size of ground grid
118 
119  float Step[3]; // Spacing delta
120  int UseTopographyFile; // Topography or flat
121  vtkStdString TopographyFile; // Name of topography data file
122  vtkPoints* Points; // Structured grid geometry
123  vtkPoints* GPoints; // Structured grid geometry for ground
124  vtkPoints* BPoints; // Unstructured grid geometry
125  float Compression; // Stretching at Z surface [0,1]
126  float Fit; // Cubic or quadratic [0,1]
127 
128  // Rectilinear coordinate spacing
133  float ZMinValue;
134 
135  // Variable information
136  int NumberOfFileVariables; // Number of variables in data file
137  int NumberOfDerivedVariables; // Number of variables derived from file
138  int NumberOfVariables; // Number of variables to display
139 
140  vtkStringArray* DivideVariables; // Divide data by density at read
141  vtkStdString* VariableName; // Names of each variable
142  int* VariableStruct; // SCALAR or VECTOR
143  int* VariableCompSize; // Number of components
144  int* VariableBasicType; // FLOAT or INTEGER
145  int* VariableByteCount; // Number of bytes in basic type
146  long int* VariableOffset; // Offset into data file
147  size_t BlockSize; // Size of every data block
148  size_t GBlockSize; // Size of every data block
149 
150  vtkFloatArray** Data; // Actual data arrays
151  vtkStdString RootDirectory; // Directory where the .wind file is.
152  vtkStdString DataDirectory; // Location of actual data
153  vtkStdString DataBaseName; // Base name of files
154 
155  // Time step information
156  int NumberOfTimeSteps; // Number of time steps
157  int TimeStepFirst; // First time step
158  int TimeStepLast; // Last time step
159  int TimeStepDelta; // Delta on time steps
160  double* TimeSteps; // Actual times available for request
161 
162  // Turbine information
163  int NumberOfBladeTowers; // Number of turbines
164  int NumberOfBladePoints; // Points for drawing parts of blades
165  int NumberOfBladeCells; // Turbines * Blades * Parts
166 
167  vtkFloatArray* XPosition; // Location of tower
168  vtkFloatArray* YPosition; // Location of tower
169  vtkFloatArray* HubHeight; // Height of tower
170  vtkFloatArray* AngularVeloc; // Angular Velocity
171  vtkFloatArray* BladeLength; // Blade length
172  vtkIntArray* BladeCount; // Number of blades per tower
173 
174  int UseTurbineFile; // Turbine data available
175  vtkStdString TurbineDirectory; // Turbine unstructured data
176  vtkStdString TurbineTowerName; // Name of tower file
177  vtkStdString TurbineBladeName; // Base name of time series blade data
178  int NumberOfLinesToSkip; // New format has lines that need to be skipped in
179  // blade files
180 
181  // Selected field of interest
183 
184  // Observer to modify this object when array selections are modified
186 
187  // Read the header file describing the dataset
188  virtual bool ReadGlobalData();
189  void ReadDataVariables(istream& inStr);
190  virtual bool FindVariableOffsets();
191 
192  // Turbine methods
193  virtual void SetupBladeData();
194  virtual void LoadBladeData(int timeStep);
195 
196  // Calculate the coordinates
200  virtual void CreateZTopography(float* zdata);
201  float GDeform(float sigma, float sigmaMax, int flag);
202  void Spline(float* x, float* y, int n, float yp1, float ypn, float* y2);
203  void Splint(float* xa, float* ya, float* y2a, int n, float x, float* y, int);
204 
205  // Load a variable from data file
206  virtual void LoadVariableData(int var);
207 
208  // Variables which must be divided by density after being read from file
209  void DivideByDensity(const char* name);
210 
211  // Calculate derived variables
212  virtual void CalculatePressure(int pres, int prespre, int tempg, int density);
213  virtual void CalculateVorticity(int vort, int uvw, int density);
214 
215  // convenience functions shared between serial and parallel version
217  vtkInformationVector* outVector, std::ostringstream& fileName, vtkStructuredGrid* field);
221  void InitPressureData(int pressure, int prespre, float*& pressureData, float*& prespreData);
223  float* pressureData, float* prespreData, const float* tempgData, const float* densityData);
224  void SetUpVorticityData(float* uData, float* vData, const float* densityData, float* vortData);
226  int var, int& numberOfComponents, float*& varData, int& planeSize, int& rowSize);
227  bool SetUpGlobalData(const std::string& fileName, std::stringstream& inStr);
228  void ProcessZCoords(float* topoData, float* zValues);
229  void ReadBladeHeader(const std::string& fileName, std::stringstream& inStr, int& numColumns);
230  void ReadBladeData(std::stringstream& inStr);
231 
233  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
234  vtkInformationVector* outputVector) override;
235 
236  static void SelectionCallback(
237  vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
238 
239  static void EventCallback(vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
240 
242 
250  vtkInformation* request, vtkInformationVector** inInfo, vtkInformationVector* outInfo) override;
251 
252 private:
253  WindBladeReaderInternal* Internal;
254 
255  vtkWindBladeReader(const vtkWindBladeReader&) = delete;
256  void operator=(const vtkWindBladeReader&) = delete;
257 };
258 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:34
vtkWindBladeReader::SelectionCallback
static void SelectionCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
vtkWindBladeReader::FindVariableOffsets
virtual bool FindVariableOffsets()
vtkWindBladeReader::NUM_BASE_SIDES
static const int NUM_BASE_SIDES
Definition: vtkWindBladeReader.h:90
vtkWindBladeReader
class for reading WindBlade data files
Definition: vtkWindBladeReader.h:49
vtkStructuredGridAlgorithm.h
vtkWindBladeReader::YSpacing
vtkFloatArray * YSpacing
Definition: vtkWindBladeReader.h:130
vtkWindBladeReader::ReadBladeData
void ReadBladeData(std::stringstream &inStr)
vtkWindBladeReader::TimeSteps
double * TimeSteps
Definition: vtkWindBladeReader.h:160
vtkWindBladeReader::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkWindBladeReader::InitPressureData
void InitPressureData(int pressure, int prespre, float *&pressureData, float *&prespreData)
vtkWindBladeReader::Splint
void Splint(float *xa, float *ya, float *y2a, int n, float x, float *y, int)
vtkWindBladeReader::LINE_SIZE
static const int LINE_SIZE
Definition: vtkWindBladeReader.h:91
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkWindBladeReader::BYTES_PER_DATA
static int BYTES_PER_DATA
Definition: vtkWindBladeReader.h:93
vtkStructuredGridAlgorithm
Superclass for algorithms that produce only structured grid as output.
Definition: vtkStructuredGridAlgorithm.h:42
vtkFloatArray
dynamic, self-adjusting array of float
Definition: vtkFloatArray.h:36
vtkWindBladeReader::ReadGlobalData
virtual bool ReadGlobalData()
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:36
vtkWindBladeReader::UseTopographyFile
int UseTopographyFile
Definition: vtkWindBladeReader.h:120
vtkWindBladeReader::LoadBladeData
virtual void LoadBladeData(int timeStep)
vtkWindBladeReader::ZTopographicValues
float * ZTopographicValues
Definition: vtkWindBladeReader.h:132
vtkWindBladeReader::Data
vtkFloatArray ** Data
Definition: vtkWindBladeReader.h:150
vtkWindBladeReader::CreateCoordinates
void CreateCoordinates()
vtkStructuredGrid
topologically regular array of data
Definition: vtkStructuredGrid.h:58
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:54
vtkWindBladeReader::GetPointArrayStatus
int GetPointArrayStatus(const char *name)
vtkWindBladeReader::SetUpGroundData
void SetUpGroundData(vtkInformationVector *outVector)
vtkWindBladeReader::XPosition
vtkFloatArray * XPosition
Definition: vtkWindBladeReader.h:167
vtkWindBladeReader::DRY_AIR_CONSTANT
static float DRY_AIR_CONSTANT
Definition: vtkWindBladeReader.h:88
vtkWindBladeReader::DIMENSION
static int DIMENSION
Definition: vtkWindBladeReader.h:92
vtkWindBladeReader::DisableAllPointArrays
void DisableAllPointArrays()
vtkWindBladeReader::TurbineBladeName
vtkStdString TurbineBladeName
Definition: vtkWindBladeReader.h:177
vtkWindBladeReader::ZMinValue
float ZMinValue
Definition: vtkWindBladeReader.h:133
vtkWindBladeReader::BPoints
vtkPoints * BPoints
Definition: vtkWindBladeReader.h:124
vtkWindBladeReader::SetUpVorticityData
void SetUpVorticityData(float *uData, float *vData, const float *densityData, float *vortData)
vtkWindBladeReader::DataBaseName
vtkStdString DataBaseName
Definition: vtkWindBladeReader.h:153
vtkWindBladeReader::Filename
char * Filename
Definition: vtkWindBladeReader.h:102
vtkWindBladeReader::InitBladeData
void InitBladeData(vtkInformationVector *outVector)
vtkWindBladeReader::VariableName
vtkStdString * VariableName
Definition: vtkWindBladeReader.h:141
vtkWindBladeReader::FLOAT
static int FLOAT
Definition: vtkWindBladeReader.h:96
vtkWindBladeReader::EnableAllPointArrays
void EnableAllPointArrays()
vtkWindBladeReader::ZSpacing
vtkFloatArray * ZSpacing
Definition: vtkWindBladeReader.h:131
vtkWindBladeReader::TimeStepDelta
int TimeStepDelta
Definition: vtkWindBladeReader.h:159
vtkWindBladeReader::TopographyFile
vtkStdString TopographyFile
Definition: vtkWindBladeReader.h:121
vtkWindBladeReader::NumberOfVariables
int NumberOfVariables
Definition: vtkWindBladeReader.h:138
vtkWindBladeReader::SetUpFieldVars
void SetUpFieldVars(vtkStructuredGrid *field)
vtkWindBladeReader::XSpacing
vtkFloatArray * XSpacing
Definition: vtkWindBladeReader.h:129
vtkWindBladeReader::CreateZTopography
virtual void CreateZTopography(float *zdata)
vtkDataArraySelection
Store on/off settings for data arrays for a vtkSource.
Definition: vtkDataArraySelection.h:35
vtkWindBladeReader::DivideByDensity
void DivideByDensity(const char *name)
vtkWindBladeReader::Points
vtkPoints * Points
Definition: vtkWindBladeReader.h:122
vtkWindBladeReader::DivideVariables
vtkStringArray * DivideVariables
Definition: vtkWindBladeReader.h:140
vtkWindBladeReader::NUM_PART_SIDES
static int NUM_PART_SIDES
Definition: vtkWindBladeReader.h:89
vtkWindBladeReader::GetPointArrayName
const char * GetPointArrayName(int index)
vtkWindBladeReader::ProcessZCoords
void ProcessZCoords(float *topoData, float *zValues)
vtkWindBladeReader::GetNumberOfPointArrays
int GetNumberOfPointArrays()
The following methods allow selective reading of solutions fields.
vtkWindBladeReader::RequestData
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
vtkWindBladeReader::RootDirectory
vtkStdString RootDirectory
Definition: vtkWindBladeReader.h:151
vtkWindBladeReader::Fit
float Fit
Definition: vtkWindBladeReader.h:126
vtkWindBladeReader::NumberOfTimeSteps
int NumberOfTimeSteps
Definition: vtkWindBladeReader.h:156
vtkWindBladeReader::ReadBladeHeader
void ReadBladeHeader(const std::string &fileName, std::stringstream &inStr, int &numColumns)
vtkWindBladeReader::GBlockSize
size_t GBlockSize
Definition: vtkWindBladeReader.h:148
vtkWindBladeReader::SetPointArrayStatus
void SetPointArrayStatus(const char *name, int status)
vtkWindBladeReader::NumberOfBladeCells
int NumberOfBladeCells
Definition: vtkWindBladeReader.h:165
vtkWindBladeReader::GetFieldOutput
vtkStructuredGrid * GetFieldOutput()
Get the reader's output.
vtkWindBladeReader::VariableStruct
int * VariableStruct
Definition: vtkWindBladeReader.h:142
vtkWindBladeReader::NumberOfFileVariables
int NumberOfFileVariables
Definition: vtkWindBladeReader.h:136
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:34
vtkWindBladeReader::VariableBasicType
int * VariableBasicType
Definition: vtkWindBladeReader.h:144
vtkWindBladeReader::HubHeight
vtkFloatArray * HubHeight
Definition: vtkWindBladeReader.h:169
vtkIntArray
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:40
vtkX3D::field
@ field
Definition: vtkX3D.h:183
vtkWindBladeReader::VariableOffset
long int * VariableOffset
Definition: vtkWindBladeReader.h:146
vtkWindBladeReader::GPoints
vtkPoints * GPoints
Definition: vtkWindBladeReader.h:123
vtkWindBladeReader::NumberOfBladeTowers
int NumberOfBladeTowers
Definition: vtkWindBladeReader.h:163
vtkWindBladeReader::PointDataArraySelection
vtkDataArraySelection * PointDataArraySelection
Definition: vtkWindBladeReader.h:182
vtkWindBladeReader::vtkWindBladeReader
vtkWindBladeReader()
vtkWindBladeReader::BladeLength
vtkFloatArray * BladeLength
Definition: vtkWindBladeReader.h:171
vtkWindBladeReader::VECTOR
static int VECTOR
Definition: vtkWindBladeReader.h:95
vtkWindBladeReader::NumberOfLinesToSkip
int NumberOfLinesToSkip
Definition: vtkWindBladeReader.h:178
vtkWindBladeReader::SetupBladeData
virtual void SetupBladeData()
vtkWindBladeReader::FillGroundCoordinates
void FillGroundCoordinates()
vtkWindBladeReader::TimeStepLast
int TimeStepLast
Definition: vtkWindBladeReader.h:158
vtkWindBladeReader::TurbineDirectory
vtkStdString TurbineDirectory
Definition: vtkWindBladeReader.h:175
vtkWindBladeReader::LoadVariableData
virtual void LoadVariableData(int var)
vtkX3D::name
@ name
Definition: vtkX3D.h:225
vtkWindBladeReader::SelectionObserver
vtkCallbackCommand * SelectionObserver
Definition: vtkWindBladeReader.h:185
vtkWindBladeReader::GDeform
float GDeform(float sigma, float sigmaMax, int flag)
vtkWindBladeReader::CalculatePressure
virtual void CalculatePressure(int pres, int prespre, int tempg, int density)
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:65
vtkWindBladeReader::AngularVeloc
vtkFloatArray * AngularVeloc
Definition: vtkWindBladeReader.h:170
vtkX3D::string
@ string
Definition: vtkX3D.h:496
vtkWindBladeReader::INTEGER
static int INTEGER
Definition: vtkWindBladeReader.h:97
vtkWindBladeReader::BladeCount
vtkIntArray * BladeCount
Definition: vtkWindBladeReader.h:172
vtkWindBladeReader::SetUpPressureData
void SetUpPressureData(float *pressureData, float *prespreData, const float *tempgData, const float *densityData)
vtkWindBladeReader::InitVariableData
void InitVariableData(int var, int &numberOfComponents, float *&varData, int &planeSize, int &rowSize)
vtkWindBladeReader::FillCoordinates
void FillCoordinates()
vtkWindBladeReader::New
static vtkWindBladeReader * New()
vtkWindBladeReader::NumberOfBladePoints
int NumberOfBladePoints
Definition: vtkWindBladeReader.h:164
vtkWindBladeReader::GetGroundOutput
vtkStructuredGrid * GetGroundOutput()
vtkWindBladeReader::UseTurbineFile
int UseTurbineFile
Definition: vtkWindBladeReader.h:174
vtkWindBladeReader::Compression
float Compression
Definition: vtkWindBladeReader.h:125
vtkWindBladeReader::SetUpGlobalData
bool SetUpGlobalData(const std::string &fileName, std::stringstream &inStr)
vtkWindBladeReader::RequestInformation
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
vtkCallbackCommand
supports function callbacks
Definition: vtkCallbackCommand.h:45
vtkWindBladeReader::TurbineTowerName
vtkStdString TurbineTowerName
Definition: vtkWindBladeReader.h:176
vtkStringArray
a vtkAbstractArray subclass for strings
Definition: vtkStringArray.h:37
vtkWindBladeReader::BlockSize
size_t BlockSize
Definition: vtkWindBladeReader.h:147
vtkWindBladeReader::NumberOfTuples
vtkIdType NumberOfTuples
Definition: vtkWindBladeReader.h:105
vtkWindBladeReader::VariableCompSize
int * VariableCompSize
Definition: vtkWindBladeReader.h:143
vtkWindBladeReader::SCALAR
static int SCALAR
Definition: vtkWindBladeReader.h:94
vtkWindBladeReader::VariableByteCount
int * VariableByteCount
Definition: vtkWindBladeReader.h:145
vtkWindBladeReader::InitFieldData
void InitFieldData(vtkInformationVector *outVector, std::ostringstream &fileName, vtkStructuredGrid *field)
vtkWindBladeReader::ReadDataVariables
void ReadDataVariables(istream &inStr)
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:93
vtkWindBladeReader::TimeStepFirst
int TimeStepFirst
Definition: vtkWindBladeReader.h:157
vtkStdString
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:35
vtkWindBladeReader::YPosition
vtkFloatArray * YPosition
Definition: vtkWindBladeReader.h:168
vtkWindBladeReader::EventCallback
static void EventCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
vtkWindBladeReader::ProcessRequest
vtkTypeBool ProcessRequest(vtkInformation *request, vtkInformationVector **inInfo, vtkInformationVector *outInfo) override
We intercept the requests to check for which port information is being requested for and if there is ...
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkWindBladeReader::DataDirectory
vtkStdString DataDirectory
Definition: vtkWindBladeReader.h:152
vtkWindBladeReader::~vtkWindBladeReader
~vtkWindBladeReader() override
vtkTypeBool
int vtkTypeBool
Definition: vtkABI.h:69
vtkWindBladeReader::FillOutputPortInformation
int FillOutputPortInformation(int, vtkInformation *) override
Fill the output port information objects for this algorithm.
vtkWindBladeReader::Spline
void Spline(float *x, float *y, int n, float yp1, float ypn, float *y2)
vtkWindBladeReader::GetBladeOutput
vtkUnstructuredGrid * GetBladeOutput()
vtkWindBladeReader::NumberOfDerivedVariables
int NumberOfDerivedVariables
Definition: vtkWindBladeReader.h:137
vtkWindBladeReader::CalculateVorticity
virtual void CalculateVorticity(int vort, int uvw, int density)