VTK  9.2.6
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 
36 #ifndef vtkNetCDFCFReader_h
37 #define vtkNetCDFCFReader_h
38 
39 #include "vtkIONetCDFModule.h" // For export macro
40 #include "vtkNetCDFReader.h"
41 
42 #include "vtkStdString.h" // Used for ivars.
43 
44 class vtkImageData;
45 class vtkPoints;
46 class vtkRectilinearGrid;
47 class vtkStructuredGrid;
49 
50 class VTKIONETCDF_EXPORT vtkNetCDFCFReader : public vtkNetCDFReader
51 {
52 public:
54  static vtkNetCDFCFReader* New();
55  void PrintSelf(ostream& os, vtkIndent indent) override;
56 
58 
63  vtkGetMacro(SphericalCoordinates, vtkTypeBool);
64  vtkSetMacro(SphericalCoordinates, vtkTypeBool);
65  vtkBooleanMacro(SphericalCoordinates, vtkTypeBool);
67 
69 
80  vtkGetMacro(VerticalScale, double);
81  vtkSetMacro(VerticalScale, double);
82  vtkGetMacro(VerticalBias, double);
83  vtkSetMacro(VerticalBias, double);
85 
87 
94  vtkGetMacro(OutputType, int);
95  virtual void SetOutputType(int type);
96  void SetOutputTypeToAutomatic() { this->SetOutputType(-1); }
97  void SetOutputTypeToImage() { this->SetOutputType(VTK_IMAGE_DATA); }
98  void SetOutputTypeToRectilinear() { this->SetOutputType(VTK_RECTILINEAR_GRID); }
99  void SetOutputTypeToStructured() { this->SetOutputType(VTK_STRUCTURED_GRID); }
100  void SetOutputTypeToUnstructured() { this->SetOutputType(VTK_UNSTRUCTURED_GRID); }
102 
106  static int CanReadFile(VTK_FILEPATH const char* filename);
107 
108 protected:
110  ~vtkNetCDFCFReader() override;
111 
113 
115  double VerticalBias;
116 
118 
119  int RequestDataObject(vtkInformation* request, vtkInformationVector** inputVector,
120  vtkInformationVector* outputVector) override;
121 
122  int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
123  vtkInformationVector* outputVector) override;
124 
125  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
126  vtkInformationVector* outputVector) override;
127 
129 
132  int ReadMetaData(int ncFD) override;
133  int IsTimeDimension(int ncFD, int dimId) override;
134  vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId) override;
136 
138  {
139  public:
140  vtkDimensionInfo() = default;
141  vtkDimensionInfo(int ncFD, int id);
142  const char* GetName() const { return this->Name.c_str(); }
144  {
149  VERTICAL_UNITS
150  };
151  UnitsEnum GetUnits() const { return this->Units; }
152  vtkSmartPointer<vtkDoubleArray> GetCoordinates() { return this->Coordinates; }
153  vtkSmartPointer<vtkDoubleArray> GetBounds() { return this->Bounds; }
154  bool GetHasRegularSpacing() const { return this->HasRegularSpacing; }
155  double GetOrigin() const { return this->Origin; }
156  double GetSpacing() const { return this->Spacing; }
157  vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
158 
159  protected:
161  int DimId;
166  double Origin, Spacing;
168  int LoadMetaData(int ncFD);
169  };
170  class vtkDimensionInfoVector;
171  friend class vtkDimensionInfoVector;
172  vtkDimensionInfoVector* DimensionInfo;
173  vtkDimensionInfo* GetDimensionInfo(int dimension);
174 
176  {
177  public:
179  : Valid(false)
180  {
181  }
182  vtkDependentDimensionInfo(int ncFD, int varId, vtkNetCDFCFReader* parent);
183  bool GetValid() const { return this->Valid; }
184  bool GetHasBounds() const { return this->HasBounds; }
185  bool GetCellsUnstructured() const { return this->CellsUnstructured; }
186  vtkSmartPointer<vtkIntArray> GetGridDimensions() const { return this->GridDimensions; }
188  {
189  return this->LongitudeCoordinates;
190  }
192  {
193  return this->LatitudeCoordinates;
194  }
195  vtkSmartPointer<vtkStringArray> GetSpecialVariables() const { return this->SpecialVariables; }
196 
197  protected:
198  bool Valid;
199  bool HasBounds;
205  int LoadMetaData(int ncFD, int varId, vtkNetCDFCFReader* parent);
206  int LoadCoordinateVariable(int ncFD, int varId, vtkDoubleArray* coords);
207  int LoadBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
208  int LoadUnstructuredBoundsVariable(int ncFD, int varId, vtkDoubleArray* coords);
209  };
211  class vtkDependentDimensionInfoVector;
212  friend class vtkDependentDimensionInfoVector;
213  vtkDependentDimensionInfoVector* DependentDimensionInfo;
214 
215  // Finds the dependent dimension information for the given set of dimensions.
216  // Returns nullptr if no information has been recorded.
217  vtkDependentDimensionInfo* FindDependentDimensionInfo(vtkIntArray* dims);
218 
224  virtual void IdentifySphericalCoordinates(
225  vtkIntArray* dimensions, int& longitudeDim, int& latitudeDim, int& verticalDim);
226 
228  {
237  COORDS_SPHERICAL_PSIDED_CELLS
238  };
239 
245  CoordinateTypesEnum CoordinateType(vtkIntArray* dimensions);
246 
250  bool DimensionsAreForPointData(vtkIntArray* dimensions) override;
251 
257  void ExtentForDimensionsAndPiece(
258  int pieceNumber, int numberOfPieces, int ghostLevels, int extent[6]);
259 
263  void GetUpdateExtentForOutput(vtkDataSet* output, int extent[6]) override;
264 
266 
269  void AddRectilinearCoordinates(vtkImageData* imageOutput);
270  void AddRectilinearCoordinates(vtkRectilinearGrid* rectilinearOutput);
271  void FakeRectilinearCoordinates(vtkRectilinearGrid* rectilinearOutput);
272  void Add1DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
273  void Add2DRectilinearCoordinates(vtkPoints* points, const int extent[6]);
274  void Add1DRectilinearCoordinates(vtkStructuredGrid* structuredOutput);
275  void Add2DRectilinearCoordinates(vtkStructuredGrid* structuredOutput);
276  void FakeStructuredCoordinates(vtkStructuredGrid* structuredOutput);
277  void Add1DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
278  void Add2DRectilinearCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
280 
282 
285  void Add1DSphericalCoordinates(vtkPoints* points, const int extent[6]);
286  void Add2DSphericalCoordinates(vtkPoints* points, const int extent[6]);
287  void Add1DSphericalCoordinates(vtkStructuredGrid* structuredOutput);
288  void Add2DSphericalCoordinates(vtkStructuredGrid* structuredOutput);
289  void Add1DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
290  void Add2DSphericalCoordinates(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
292 
296  void AddStructuredCells(vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
297 
299 
302  void AddUnstructuredRectilinearCoordinates(
303  vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
304  void AddUnstructuredSphericalCoordinates(
305  vtkUnstructuredGrid* unstructuredOutput, const int extent[6]);
307 
308 private:
309  vtkNetCDFCFReader(const vtkNetCDFCFReader&) = delete;
310  void operator=(const vtkNetCDFCFReader&) = delete;
311 };
312 
313 #endif // vtkNetCDFCFReader_h
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:38
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:83
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:80
Store vtkAlgorithm input/output information.
Reads netCDF files that follow the CF convention.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
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:69
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:45
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:39
vtkSmartPointer< vtkDoubleArray > GetLatitudeCoordinates() const
topologically and geometrically regular array of data
Definition: vtkImageData.h:53
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:79
represent and manipulate 3D points
Definition: vtkPoints.h:39
#define VTK_UNSTRUCTURED_GRID
Definition: vtkType.h:81
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
vtkSmartPointer< vtkIntArray > GetGridDimensions() const