VTK  9.4.20250208
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
92
96 void SetTimeDimensionName(const char* name);
97 void SetLatitudeDimensionName(const char* name);
98 void SetLongitudeDimensionName(const char* name);
99 void SetVerticalDimensionName(const char* name);
101
103
107 const char* GetTimeDimensionName();
112
113protected:
116
118
121
123
125 vtkInformationVector* outputVector) override;
126
128 vtkInformationVector* outputVector) override;
129
131 vtkInformationVector* outputVector) override;
132
134
137 int ReadMetaData(int ncFD) override;
138 int IsTimeDimension(int ncFD, int dimId) override;
139 vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId) override;
141
143 {
144 public:
145 vtkDimensionInfo() = default;
146 vtkDimensionInfo(vtkNetCDFAccessor* accessor, int ncFD, int id,
147 const std::vector<std::string>& dimensionNames);
148
149 const char* GetName() const { return this->Name.c_str(); }
151 {
157 NUMBER_OF_UNITS
158 };
159 UnitsEnum GetUnits() const { return this->Units; }
160 vtkSmartPointer<vtkDoubleArray> GetCoordinates() { return this->Coordinates; }
162 bool GetHasRegularSpacing() const { return this->HasRegularSpacing; }
163 double GetOrigin() const { return this->Origin; }
164 double GetSpacing() const { return this->Spacing; }
165 vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
166 void SetUnitsIfSpecialDimensionOverriden(UnitsEnum unit, const char* name);
167
168 protected:
169 vtkNetCDFAccessor* Accessor = nullptr;
171 int DimId;
176 double Origin, Spacing;
178 std::vector<std::string> SpecialDimensionOverrideNames;
179 int LoadMetaData(int ncFD);
180 };
182 {
183 this->SpecialDimensionOverrideNames[dim] = name;
184 }
186
187 class vtkDimensionInfoVector;
188 friend class vtkDimensionInfoVector;
189 vtkDimensionInfoVector* DimensionInfo;
191
193 {
194 public:
196 : Accessor(accessor)
197 , Valid(false)
198 {
199 }
201 vtkNetCDFAccessor* accessor, int ncFD, int varId, vtkNetCDFCFReader* parent);
202 bool GetValid() const { return this->Valid; }
203 bool GetHasBounds() const { return this->HasBounds; }
204 bool GetCellsUnstructured() const { return this->CellsUnstructured; }
205 vtkSmartPointer<vtkIntArray> GetGridDimensions() const { return this->GridDimensions; }
207 {
208 return this->LongitudeCoordinates;
209 }
211 {
212 return this->LatitudeCoordinates;
213 }
214 vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
215
216 protected:
218 bool Valid;
225 int LoadMetaData(int ncFD, int varId, vtkNetCDFCFReader* parent);
226 int LoadCoordinateVariable(int ncFD, int varId, vtkDoubleArray* coords);
227 int LoadBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
228 int LoadUnstructuredBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
229 };
231 class vtkDependentDimensionInfoVector;
232 friend class vtkDependentDimensionInfoVector;
233 vtkDependentDimensionInfoVector* DependentDimensionInfo;
234
235 // Finds the dependent dimension information for the given set of dimensions.
236 // Returns nullptr if no information has been recorded.
238
245 vtkIntArray* dimensions, int& longitudeDim, int& latitudeDim, int& verticalDim);
246
248 {
257 COORDS_SPHERICAL_PSIDED_CELLS
258 };
259
266
270 bool DimensionsAreForPointData(vtkIntArray* dimensions) override;
271
278 int pieceNumber, int numberOfPieces, int ghostLevels, int extent[6]);
279
283 void GetUpdateExtentForOutput(vtkDataSet* output, int extent[6]) override;
284
286
292 void Add1DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
293 void Add2DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
297 void Add1DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
298 void Add2DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
300
302
305 void Add1DSphericalCoordinates(vtkPoints* points, const int extent[6]);
306 void Add2DSphericalCoordinates(vtkPoints* points, const int extent[6]);
309 void Add1DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
310 void Add2DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
312
316 void AddStructuredCells(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
317
319
323 vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
325 vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
327
328 std::vector<std::string> SpecialDimensionOverrideNames;
329
330private:
331 vtkNetCDFCFReader(const vtkNetCDFCFReader&) = delete;
332 void operator=(const vtkNetCDFCFReader&) = delete;
333};
334
335VTK_ABI_NAMESPACE_END
336#endif // vtkNetCDFCFReader_h
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
dynamic, self-adjusting array of double
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
vtkSmartPointer< vtkDoubleArray > GetLatitudeCoordinates() const
vtkSmartPointer< vtkIntArray > GridDimensions
int LoadUnstructuredBoundsVariable(int ncFD, int varId, vtkDoubleArray *coords)
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
vtkSmartPointer< vtkDoubleArray > LongitudeCoordinates
vtkDependentDimensionInfo(vtkNetCDFAccessor *accessor)
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)
vtkDependentDimensionInfo(vtkNetCDFAccessor *accessor, int ncFD, int varId, vtkNetCDFCFReader *parent)
vtkDimensionInfo(vtkNetCDFAccessor *accessor, int ncFD, int id, const std::vector< std::string > &dimensionNames)
vtkSmartPointer< vtkDoubleArray > GetCoordinates()
std::vector< std::string > SpecialDimensionOverrideNames
vtkSmartPointer< vtkDoubleArray > Bounds
vtkSmartPointer< vtkDoubleArray > GetBounds()
vtkSmartPointer< vtkStringArray > SpecialVariables
vtkSmartPointer< vtkDoubleArray > Coordinates
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
void SetUnitsIfSpecialDimensionOverriden(UnitsEnum unit, const char *name)
Reads netCDF files that follow the CF convention.
const char * GetVerticalDimensionName()
Names for Time, Latitude, Longitude and Vertical.
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.
const char * GetLatitudeDimensionName()
Names for Time, Latitude, Longitude and Vertical.
void SetLatitudeDimensionName(const char *name)
Names for Time, Latitude, Longitude and Vertical which can be set by the user for datasets that don't...
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
const char * GetSpecialDimensionName(vtkDimensionInfo::UnitsEnum dim)
void Add2DSphericalCoordinates(vtkStructuredGrid *structuredOutput)
Internal methods for setting spherical coordinates.
void SetSpecialDimensionOverrideName(vtkDimensionInfo::UnitsEnum dim, const char *name)
const char * GetLongitudeDimensionName()
Names for Time, Latitude, Longitude and Vertical.
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.
const char * GetTimeDimensionName()
Names for Time, Latitude, Longitude and Vertical.
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 SetVerticalDimensionName(const char *name)
Names for Time, Latitude, Longitude and Vertical which can be set by the user for datasets that don't...
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 SetLongitudeDimensionName(const char *name)
Names for Time, Latitude, Longitude and Vertical which can be set by the user for datasets that don't...
void Add2DSphericalCoordinates(vtkPoints *points, const int extent[6])
Internal methods for setting spherical coordinates.
void SetTimeDimensionName(const char *name)
Names for Time, Latitude, Longitude and Vertical which can be set by the user for datasets that don't...
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.
std::vector< std::string > SpecialDimensionOverrideNames
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:139
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