VTK  9.3.1
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 
26 VTK_ABI_NAMESPACE_BEGIN
27 class vtkImageData;
28 class vtkPoints;
29 class vtkRectilinearGrid;
30 class vtkStructuredGrid;
32 
33 class VTKIONETCDF_EXPORT vtkNetCDFCFReader : public vtkNetCDFReader
34 {
35 public:
37  static vtkNetCDFCFReader* New();
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 
91 protected:
93  ~vtkNetCDFCFReader() override;
94 
96 
97  double VerticalScale;
98  double VerticalBias;
99 
101 
102  int RequestDataObject(vtkInformation* request, vtkInformationVector** inputVector,
103  vtkInformationVector* outputVector) override;
104 
105  int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
106  vtkInformationVector* outputVector) override;
107 
108  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
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; }
136  vtkSmartPointer<vtkDoubleArray> GetBounds() { return this->Bounds; }
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;
156  vtkDimensionInfo* GetDimensionInfo(int dimension);
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;
182  bool HasBounds;
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.
200  vtkDependentDimensionInfo* FindDependentDimensionInfo(vtkIntArray* dims);
201 
207  virtual void IdentifySphericalCoordinates(
208  vtkIntArray* dimensions, int& longitudeDim, int& latitudeDim, int& verticalDim);
209 
211  {
220  COORDS_SPHERICAL_PSIDED_CELLS
221  };
222 
228  CoordinateTypesEnum CoordinateType(vtkIntArray* dimensions);
229 
233  bool DimensionsAreForPointData(vtkIntArray* dimensions) override;
234 
240  void ExtentForDimensionsAndPiece(
241  int pieceNumber, int numberOfPieces, int ghostLevels, int extent[6]);
242 
246  void GetUpdateExtentForOutput(vtkDataSet* output, int extent[6]) override;
247 
249 
252  void AddRectilinearCoordinates(vtkImageData* imageOutput);
253  void AddRectilinearCoordinates(vtkRectilinearGrid* rectilinearOutput);
254  void FakeRectilinearCoordinates(vtkRectilinearGrid* rectilinearOutput);
255  void Add1DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
256  void Add2DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
257  void Add1DRectilinearCoordinates(vtkStructuredGrid* structuredOutput);
258  void Add2DRectilinearCoordinates(vtkStructuredGrid* structuredOutput);
259  void FakeStructuredCoordinates(vtkStructuredGrid* structuredOutput);
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]);
270  void Add1DSphericalCoordinates(vtkStructuredGrid* structuredOutput);
271  void Add2DSphericalCoordinates(vtkStructuredGrid* structuredOutput);
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 
285  void AddUnstructuredRectilinearCoordinates(
286  vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
287  void AddUnstructuredSphericalCoordinates(
288  vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
290 
291 private:
292  vtkNetCDFCFReader(const vtkNetCDFCFReader&) = delete;
293  void operator=(const vtkNetCDFCFReader&) = delete;
294 };
295 
296 VTK_ABI_NAMESPACE_END
297 #endif // vtkNetCDFCFReader_h
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:28
vtkSmartPointer< vtkDoubleArray > Coordinates
vtkSmartPointer< vtkDoubleArray > GetCoordinates()
a dataset that is topologically regular with variable spacing in the three coordinate directions ...
#define VTK_IMAGE_DATA
Definition: vtkType.h:71
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
vtkSmartPointer< vtkDoubleArray > LatitudeCoordinates
void SetOutputTypeToStructured()
Set/get the data type of the output.
#define VTK_RECTILINEAR_GRID
Definition: vtkType.h:68
Store vtkAlgorithm input/output information.
Reads netCDF files that follow the CF convention.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:52
virtual bool DimensionsAreForPointData(vtkIntArray *vtkNotUsed(dimensions))
Called internally to determine whether a variable with the given set of dimensions should be loaded a...
vtkDimensionInfoVector * DimensionInfo
vtkSmartPointer< vtkStringArray > SpecialVariables
void SetOutputTypeToUnstructured()
Set/get the data type of the output.
void SetOutputTypeToImage()
Set/get the data type of the output.
virtual int ReadMetaData(int ncFD)
Reads meta data and populates ivars.
A superclass for reading netCDF files.
dynamic, self-adjusting array of double
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
int vtkTypeBool
Definition: vtkABI.h:64
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:34
void SetOutputTypeToAutomatic()
Set/get the data type of the output.
static vtkNetCDFReader * New()
virtual void GetUpdateExtentForOutput(vtkDataSet *output, int extent[6])
Retrieves the update extent for the output object.
a simple class to control print indentation
Definition: vtkIndent.h:28
vtkSmartPointer< vtkDoubleArray > GetLatitudeCoordinates() const
topologically and geometrically regular array of data
Definition: vtkImageData.h:42
dataset represents arbitrary combinations of all possible cell types
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
virtual int IsTimeDimension(int ncFD, int dimId)
Determines whether the given variable is a time dimension.
vtkSmartPointer< vtkStringArray > SpecialVariables
vtkTypeBool SphericalCoordinates
vtkSmartPointer< vtkDoubleArray > GetBounds()
virtual vtkSmartPointer< vtkDoubleArray > GetTimeValues(int ncFD, int dimId)
Given a dimension already determined to be a time dimension (via a call to IsTimeDimension) returns a...
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSmartPointer< vtkDoubleArray > LongitudeCoordinates
#define VTK_FILEPATH
vtkDependentDimensionInfoVector * DependentDimensionInfo
topologically regular array of data
void SetOutputTypeToRectilinear()
Set/get the data type of the output.
Store zero or more vtkInformation instances.
vtkSmartPointer< vtkDoubleArray > GetLongitudeCoordinates() const
vtkSmartPointer< vtkDoubleArray > Bounds
vtkSmartPointer< vtkIntArray > GridDimensions
#define VTK_STRUCTURED_GRID
Definition: vtkType.h:67
represent and manipulate 3D points
Definition: vtkPoints.h:28
#define VTK_UNSTRUCTURED_GRID
Definition: vtkType.h:69
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
vtkSmartPointer< vtkIntArray > GetGridDimensions() const