VTK  9.3.0
vtkNetCDFCFReader.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright 2008 Sandia Corporation
3// SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-LANL-California-USGov
4
18#ifndef vtkNetCDFCFReader_h
19#define vtkNetCDFCFReader_h
20
21#include "vtkIONetCDFModule.h" // For export macro
22#include "vtkNetCDFReader.h"
23
24#include "vtkStdString.h" // Used for ivars.
25
26VTK_ABI_NAMESPACE_BEGIN
27class vtkImageData;
28class vtkPoints;
32
33class VTKIONETCDF_EXPORT vtkNetCDFCFReader : public vtkNetCDFReader
34{
35public:
38 void PrintSelf(ostream& os, vtkIndent indent) override;
39
41
46 vtkGetMacro(SphericalCoordinates, vtkTypeBool);
47 vtkSetMacro(SphericalCoordinates, vtkTypeBool);
48 vtkBooleanMacro(SphericalCoordinates, vtkTypeBool);
50
52
63 vtkGetMacro(VerticalScale, double);
64 vtkSetMacro(VerticalScale, double);
65 vtkGetMacro(VerticalBias, double);
66 vtkSetMacro(VerticalBias, double);
68
70
77 vtkGetMacro(OutputType, int);
78 virtual void SetOutputType(int type);
79 void SetOutputTypeToAutomatic() { this->SetOutputType(-1); }
80 void SetOutputTypeToImage() { this->SetOutputType(VTK_IMAGE_DATA); }
81 void SetOutputTypeToRectilinear() { this->SetOutputType(VTK_RECTILINEAR_GRID); }
82 void SetOutputTypeToStructured() { this->SetOutputType(VTK_STRUCTURED_GRID); }
83 void SetOutputTypeToUnstructured() { this->SetOutputType(VTK_UNSTRUCTURED_GRID); }
85
89 static int CanReadFile(VTK_FILEPATH const char* filename);
90
91protected:
94
96
99
101
103 vtkInformationVector* outputVector) override;
104
106 vtkInformationVector* outputVector) override;
107
109 vtkInformationVector* outputVector) override;
110
112
115 int ReadMetaData(int ncFD) override;
116 int IsTimeDimension(int ncFD, int dimId) override;
117 vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId) override;
119
121 {
122 public:
123 vtkDimensionInfo() = default;
124 vtkDimensionInfo(int ncFD, int id);
125 const char* GetName() const { return this->Name.c_str(); }
127 {
132 VERTICAL_UNITS
133 };
134 UnitsEnum GetUnits() const { return this->Units; }
135 vtkSmartPointer<vtkDoubleArray> GetCoordinates() { return this->Coordinates; }
137 bool GetHasRegularSpacing() const { return this->HasRegularSpacing; }
138 double GetOrigin() const { return this->Origin; }
139 double GetSpacing() const { return this->Spacing; }
140 vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
141
142 protected:
144 int DimId;
149 double Origin, Spacing;
151 int LoadMetaData(int ncFD);
152 };
153 class vtkDimensionInfoVector;
154 friend class vtkDimensionInfoVector;
155 vtkDimensionInfoVector* DimensionInfo;
157
159 {
160 public:
162 : Valid(false)
163 {
164 }
165 vtkDependentDimensionInfo(int ncFD, int varId, vtkNetCDFCFReader* parent);
166 bool GetValid() const { return this->Valid; }
167 bool GetHasBounds() const { return this->HasBounds; }
168 bool GetCellsUnstructured() const { return this->CellsUnstructured; }
169 vtkSmartPointer<vtkIntArray> GetGridDimensions() const { return this->GridDimensions; }
171 {
172 return this->LongitudeCoordinates;
173 }
175 {
176 return this->LatitudeCoordinates;
177 }
178 vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
179
180 protected:
181 bool Valid;
188 int LoadMetaData(int ncFD, int varId, vtkNetCDFCFReader* parent);
189 int LoadCoordinateVariable(int ncFD, int varId, vtkDoubleArray* coords);
190 int LoadBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
191 int LoadUnstructuredBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
192 };
194 class vtkDependentDimensionInfoVector;
195 friend class vtkDependentDimensionInfoVector;
196 vtkDependentDimensionInfoVector* DependentDimensionInfo;
197
198 // Finds the dependent dimension information for the given set of dimensions.
199 // Returns nullptr if no information has been recorded.
201
208 vtkIntArray* dimensions, int& longitudeDim, int& latitudeDim, int& verticalDim);
209
211 {
220 COORDS_SPHERICAL_PSIDED_CELLS
221 };
222
229
233 bool DimensionsAreForPointData(vtkIntArray* dimensions) override;
234
241 int pieceNumber, int numberOfPieces, int ghostLevels, int extent[6]);
242
246 void GetUpdateExtentForOutput(vtkDataSet* output, int extent[6]) override;
247
249
255 void Add1DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
256 void Add2DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
260 void Add1DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
261 void Add2DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
263
265
268 void Add1DSphericalCoordinates(vtkPoints* points, const int extent[6]);
269 void Add2DSphericalCoordinates(vtkPoints* points, const int extent[6]);
272 void Add1DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
273 void Add2DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
275
279 void AddStructuredCells(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
280
282
286 vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
288 vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
290
291private:
292 vtkNetCDFCFReader(const vtkNetCDFCFReader&) = delete;
293 void operator=(const vtkNetCDFCFReader&) = delete;
294};
295
296VTK_ABI_NAMESPACE_END
297#endif // vtkNetCDFCFReader_h
abstract class to specify dataset behavior
Definition vtkDataSet.h:53
dynamic, self-adjusting array of double
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition vtkIntArray.h:35
vtkSmartPointer< vtkDoubleArray > GetLatitudeCoordinates() const
vtkSmartPointer< vtkIntArray > GridDimensions
int LoadUnstructuredBoundsVariable(int ncFD, int varId, vtkDoubleArray *coords)
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
vtkSmartPointer< vtkDoubleArray > LongitudeCoordinates
vtkDependentDimensionInfo(int ncFD, int varId, vtkNetCDFCFReader *parent)
vtkSmartPointer< vtkStringArray > SpecialVariables
int LoadMetaData(int ncFD, int varId, vtkNetCDFCFReader *parent)
vtkSmartPointer< vtkDoubleArray > LatitudeCoordinates
int LoadCoordinateVariable(int ncFD, int varId, vtkDoubleArray *coords)
vtkSmartPointer< vtkDoubleArray > GetLongitudeCoordinates() const
vtkSmartPointer< vtkIntArray > GetGridDimensions() const
int LoadBoundsVariable(int ncFD, int varId, vtkDoubleArray *coords)
vtkSmartPointer< vtkDoubleArray > GetCoordinates()
vtkSmartPointer< vtkDoubleArray > Bounds
vtkSmartPointer< vtkDoubleArray > GetBounds()
vtkSmartPointer< vtkStringArray > SpecialVariables
vtkSmartPointer< vtkDoubleArray > Coordinates
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
Reads netCDF files that follow the CF convention.
void Add1DSphericalCoordinates(vtkPoints *points, const int extent[6])
Internal methods for setting spherical coordinates.
void SetOutputTypeToUnstructured()
Set/get the data type of the output.
vtkDependentDimensionInfo * FindDependentDimensionInfo(vtkIntArray *dims)
virtual void SetOutputType(int type)
Set/get the data type of the output.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
vtkDimensionInfoVector * DimensionInfo
void SetOutputTypeToImage()
Set/get the data type of the output.
void AddUnstructuredSphericalCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal methods for creating unstructured cells.
void FakeStructuredCoordinates(vtkStructuredGrid *structuredOutput)
Internal methods for setting rectilinear coordinates.
void Add2DSphericalCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal methods for setting spherical coordinates.
void Add1DSphericalCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal methods for setting spherical coordinates.
~vtkNetCDFCFReader() override
void Add2DSphericalCoordinates(vtkStructuredGrid *structuredOutput)
Internal methods for setting spherical coordinates.
void Add1DRectilinearCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal methods for setting rectilinear coordinates.
void SetOutputTypeToRectilinear()
Set/get the data type of the output.
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
CoordinateTypesEnum CoordinateType(vtkIntArray *dimensions)
Based on the given dimensions and the current state of the reader, returns how the coordinates should...
void AddRectilinearCoordinates(vtkImageData *imageOutput)
Internal methods for setting rectilinear coordinates.
void GetUpdateExtentForOutput(vtkDataSet *output, int extent[6]) override
Overridden to retrieve stored extent for unstructured data.
vtkDimensionInfo * GetDimensionInfo(int dimension)
vtkTypeBool SphericalCoordinates
static int CanReadFile(VTK_FILEPATH const char *filename)
Returns true if the given file can be read.
int IsTimeDimension(int ncFD, int dimId) override
Interprets the special conventions of COARDS.
void Add1DRectilinearCoordinates(vtkStructuredGrid *structuredOutput)
Internal methods for setting rectilinear coordinates.
void SetOutputTypeToStructured()
Set/get the data type of the output.
static vtkNetCDFCFReader * New()
bool DimensionsAreForPointData(vtkIntArray *dimensions) override
Returns false for spherical dimensions, which should use cell data.
void Add1DSphericalCoordinates(vtkStructuredGrid *structuredOutput)
Internal methods for setting spherical coordinates.
void SetOutputTypeToAutomatic()
Set/get the data type of the output.
void Add2DRectilinearCoordinates(vtkStructuredGrid *structuredOutput)
Internal methods for setting rectilinear coordinates.
void Add2DSphericalCoordinates(vtkPoints *points, const int extent[6])
Internal methods for setting spherical coordinates.
void AddRectilinearCoordinates(vtkRectilinearGrid *rectilinearOutput)
Internal methods for setting rectilinear coordinates.
void Add2DRectilinearCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal methods for setting rectilinear coordinates.
void FakeRectilinearCoordinates(vtkRectilinearGrid *rectilinearOutput)
Internal methods for setting rectilinear coordinates.
void Add2DRectilinearCoordinates(vtkPoints *points, const int extent[6])
Internal methods for setting rectilinear coordinates.
vtkSmartPointer< vtkDoubleArray > GetTimeValues(int ncFD, int dimId) override
Interprets the special conventions of COARDS.
void ExtentForDimensionsAndPiece(int pieceNumber, int numberOfPieces, int ghostLevels, int extent[6])
Convenience function that takes piece information and then returns a set of extents to load based on ...
void Add1DRectilinearCoordinates(vtkPoints *points, const int extent[6])
Internal methods for setting rectilinear coordinates.
void AddStructuredCells(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal method for building unstructred cells that match structured cells.
int ReadMetaData(int ncFD) override
Interprets the special conventions of COARDS.
void AddUnstructuredRectilinearCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal methods for creating unstructured cells.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
vtkDependentDimensionInfoVector * DependentDimensionInfo
virtual void IdentifySphericalCoordinates(vtkIntArray *dimensions, int &longitudeDim, int &latitudeDim, int &verticalDim)
Given the list of dimensions, identify the longitude, latitude, and vertical dimensions.
A superclass for reading netCDF files.
represent and manipulate 3D points
Definition vtkPoints.h:29
a dataset that is topologically regular with variable spacing in the three coordinate directions
Hold a reference to a vtkObjectBase instance.
Wrapper around std::string to keep symbols short.
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_IMAGE_DATA
Definition vtkType.h:71
#define VTK_RECTILINEAR_GRID
Definition vtkType.h:68
#define VTK_UNSTRUCTURED_GRID
Definition vtkType.h:69
#define VTK_STRUCTURED_GRID
Definition vtkType.h:67
#define VTK_FILEPATH