VTK  9.1.0
vtkNetCDFCFReader.h
Go to the documentation of this file.
1// -*- c++ -*-
2/*=========================================================================
3
4 Program: Visualization Toolkit
5 Module: vtkNetCDFCFReader.h
6
7 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
8 All rights reserved.
9 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
10
11 This software is distributed WITHOUT ANY WARRANTY; without even
12 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
13 PURPOSE. See the above copyright notice for more information.
14
15=========================================================================*/
16
17/*-------------------------------------------------------------------------
18 Copyright 2008 Sandia Corporation.
19 Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
20 the U.S. Government retains certain rights in this software.
21-------------------------------------------------------------------------*/
22
33#ifndef vtkNetCDFCFReader_h
34#define vtkNetCDFCFReader_h
35
36#include "vtkIONetCDFModule.h" // For export macro
37#include "vtkNetCDFReader.h"
38
39#include "vtkStdString.h" // Used for ivars.
40
41class vtkImageData;
42class vtkPoints;
46
47class VTKIONETCDF_EXPORT vtkNetCDFCFReader : public vtkNetCDFReader
48{
49public:
52 void PrintSelf(ostream& os, vtkIndent indent) override;
53
55
60 vtkGetMacro(SphericalCoordinates, vtkTypeBool);
61 vtkSetMacro(SphericalCoordinates, vtkTypeBool);
62 vtkBooleanMacro(SphericalCoordinates, vtkTypeBool);
64
66
77 vtkGetMacro(VerticalScale, double);
78 vtkSetMacro(VerticalScale, double);
79 vtkGetMacro(VerticalBias, double);
80 vtkSetMacro(VerticalBias, double);
82
84
91 vtkGetMacro(OutputType, int);
92 virtual void SetOutputType(int type);
93 void SetOutputTypeToAutomatic() { this->SetOutputType(-1); }
94 void SetOutputTypeToImage() { this->SetOutputType(VTK_IMAGE_DATA); }
95 void SetOutputTypeToRectilinear() { this->SetOutputType(VTK_RECTILINEAR_GRID); }
96 void SetOutputTypeToStructured() { this->SetOutputType(VTK_STRUCTURED_GRID); }
97 void SetOutputTypeToUnstructured() { this->SetOutputType(VTK_UNSTRUCTURED_GRID); }
99
103 static int CanReadFile(VTK_FILEPATH const char* filename);
104
105protected:
108
110
113
115
117 vtkInformationVector* outputVector) override;
118
120 vtkInformationVector* outputVector) override;
121
123 vtkInformationVector* outputVector) override;
124
126
129 int ReadMetaData(int ncFD) override;
130 int IsTimeDimension(int ncFD, int dimId) override;
131 vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId) override;
133
135 {
136 public:
137 vtkDimensionInfo() = default;
138 vtkDimensionInfo(int ncFD, int id);
139 const char* GetName() const { return this->Name.c_str(); }
141 {
146 VERTICAL_UNITS
147 };
148 UnitsEnum GetUnits() const { return this->Units; }
149 vtkSmartPointer<vtkDoubleArray> GetCoordinates() { return this->Coordinates; }
151 bool GetHasRegularSpacing() const { return this->HasRegularSpacing; }
152 double GetOrigin() const { return this->Origin; }
153 double GetSpacing() const { return this->Spacing; }
154 vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
155
156 protected:
158 int DimId;
163 double Origin, Spacing;
165 int LoadMetaData(int ncFD);
166 };
167 class vtkDimensionInfoVector;
168 friend class vtkDimensionInfoVector;
169 vtkDimensionInfoVector* DimensionInfo;
171
173 {
174 public:
176 : Valid(false)
177 {
178 }
179 vtkDependentDimensionInfo(int ncFD, int varId, vtkNetCDFCFReader* parent);
180 bool GetValid() const { return this->Valid; }
181 bool GetHasBounds() const { return this->HasBounds; }
182 bool GetCellsUnstructured() const { return this->CellsUnstructured; }
183 vtkSmartPointer<vtkIntArray> GetGridDimensions() const { return this->GridDimensions; }
185 {
186 return this->LongitudeCoordinates;
187 }
189 {
190 return this->LatitudeCoordinates;
191 }
192 vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
193
194 protected:
195 bool Valid;
202 int LoadMetaData(int ncFD, int varId, vtkNetCDFCFReader* parent);
203 int LoadCoordinateVariable(int ncFD, int varId, vtkDoubleArray* coords);
204 int LoadBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
205 int LoadUnstructuredBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
206 };
208 class vtkDependentDimensionInfoVector;
209 friend class vtkDependentDimensionInfoVector;
210 vtkDependentDimensionInfoVector* DependentDimensionInfo;
211
212 // Finds the dependent dimension information for the given set of dimensions.
213 // Returns nullptr if no information has been recorded.
215
222 vtkIntArray* dimensions, int& longitudeDim, int& latitudeDim, int& verticalDim);
223
225 {
234 COORDS_SPHERICAL_PSIDED_CELLS
235 };
236
243
247 bool DimensionsAreForPointData(vtkIntArray* dimensions) override;
248
255 int pieceNumber, int numberOfPieces, int ghostLevels, int extent[6]);
256
260 void GetUpdateExtentForOutput(vtkDataSet* output, int extent[6]) override;
261
263
274 void Add1DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
275 void Add2DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
277
279
286 void Add1DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
287 void Add2DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
289
293 void AddStructuredCells(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
294
296
300 vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
302 vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
304
305private:
306 vtkNetCDFCFReader(const vtkNetCDFCFReader&) = delete;
307 void operator=(const vtkNetCDFCFReader&) = delete;
308};
309
310#endif // vtkNetCDFCFReader_h
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
dynamic, self-adjusting array of double
topologically and geometrically regular array of data
Definition: vtkImageData.h:157
a simple class to control print indentation
Definition: vtkIndent.h:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:149
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:143
a dataset that is topologically regular with variable spacing in the three coordinate directions
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:105
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types
@ points
Definition: vtkX3D.h:452
@ extent
Definition: vtkX3D.h:351
@ type
Definition: vtkX3D.h:522
int vtkTypeBool
Definition: vtkABI.h:69
#define VTK_IMAGE_DATA
Definition: vtkType.h:83
#define VTK_RECTILINEAR_GRID
Definition: vtkType.h:80
#define VTK_UNSTRUCTURED_GRID
Definition: vtkType.h:81
#define VTK_STRUCTURED_GRID
Definition: vtkType.h:79
#define VTK_FILEPATH