VTK
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:
55  virtual void PrintSelf(ostream &os, vtkIndent indent);
56 
58 
63  vtkGetMacro(SphericalCoordinates, int);
64  vtkSetMacro(SphericalCoordinates, int);
65  vtkBooleanMacro(SphericalCoordinates, int);
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); }
101  this->SetOutputType(VTK_UNSTRUCTURED_GRID);
102  }
104 
108  static int CanReadFile(const char *filename);
109 
110 protected:
113 
115 
117  double VerticalBias;
118 
120 
121  virtual int RequestDataObject(vtkInformation *request,
122  vtkInformationVector **inputVector,
123  vtkInformationVector *outputVector);
124 
125  virtual int RequestInformation(vtkInformation *request,
126  vtkInformationVector **inputVector,
127  vtkInformationVector *outputVector);
128 
129  virtual int RequestData(vtkInformation *request,
130  vtkInformationVector **inputVector,
131  vtkInformationVector *outputVector);
132 
134 
137  virtual int ReadMetaData(int ncFD);
138  virtual int IsTimeDimension(int ncFD, int dimId);
139  virtual vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId);
141 
143  public:
145  vtkDimensionInfo(int ncFD, int id);
146  const char *GetName() const { return this->Name.c_str(); }
147  enum UnitsEnum {
152  VERTICAL_UNITS
153  };
154  UnitsEnum GetUnits() const { return this->Units; }
155  vtkSmartPointer<vtkDoubleArray> GetCoordinates() {return this->Coordinates;}
156  vtkSmartPointer<vtkDoubleArray> GetBounds() { return this->Bounds; }
157  bool GetHasRegularSpacing() const { return this->HasRegularSpacing; }
158  double GetOrigin() const { return this->Origin; }
159  double GetSpacing() const { return this->Spacing; }
161  return this->SpecialVariables;
162  }
163  protected:
165  int DimId;
170  double Origin, Spacing;
172  int LoadMetaData(int ncFD);
173  };
174  class vtkDimensionInfoVector;
175  friend class vtkDimensionInfoVector;
176  vtkDimensionInfoVector *DimensionInfo;
178 
180  public:
181  vtkDependentDimensionInfo() : Valid(false) { };
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; }
187  return this->GridDimensions;
188  }
190  return this->LongitudeCoordinates;
191  }
193  return this->LatitudeCoordinates;
194  }
196  return this->SpecialVariables;
197  }
198  protected:
199  bool Valid;
200  bool HasBounds;
206  int LoadMetaData(int ncFD, int varId, vtkNetCDFCFReader *parent);
207  int LoadCoordinateVariable(int ncFD, int varId, vtkDoubleArray *coords);
208  int LoadBoundsVariable(int ncFD, int varId, vtkDoubleArray *coords);
209  int LoadUnstructuredBoundsVariable(int ncFD, int varId,
210  vtkDoubleArray *coords);
211  };
213  class vtkDependentDimensionInfoVector;
214  friend class vtkDependentDimensionInfoVector;
215  vtkDependentDimensionInfoVector *DependentDimensionInfo;
216 
217  // Finds the dependent dimension information for the given set of dimensions.
218  // Returns NULL if no information has been recorded.
220 
226  virtual void IdentifySphericalCoordinates(vtkIntArray *dimensions,
227  int &longitudeDim,
228  int &latitudeDim,
229  int &verticalDim);
230 
240  COORDS_SPHERICAL_PSIDED_CELLS
241  };
242 
249 
253  virtual bool DimensionsAreForPointData(vtkIntArray *dimensions);
254 
260  void ExtentForDimensionsAndPiece(int pieceNumber,
261  int numberOfPieces,
262  int ghostLevels,
263  int extent[6]);
264 
268  virtual void GetUpdateExtentForOutput(vtkDataSet *output, int extent[6]);
269 
271 
283  const int extent[6]);
285  const int extent[6]);
287 
289 
297  const int extent[6]);
299  const int extent[6]);
301 
305  void AddStructuredCells(vtkUnstructuredGrid *unstructuredOutput,
306  const int extent[6]);
307 
309 
313  vtkUnstructuredGrid *unstructuredOutput,
314  const int extent[6]);
316  vtkUnstructuredGrid *unstructuredOutput,
317  const int extent[6]);
319 
320 
321 private:
322  vtkNetCDFCFReader(const vtkNetCDFCFReader &) VTK_DELETE_FUNCTION;
323  void operator=(const vtkNetCDFCFReader &) VTK_DELETE_FUNCTION;
324 };
325 
326 #endif //vtkNetCDFCFReader_h
327 
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
dynamic, self-adjusting array of double
topologically and geometrically regular array of data
Definition: vtkImageData.h:46
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:46
vtkSmartPointer< vtkDoubleArray > GetLatitudeCoordinates() const
vtkSmartPointer< vtkStringArray > GetSpecialVariables() const
vtkSmartPointer< vtkIntArray > GridDimensions
int LoadUnstructuredBoundsVariable(int ncFD, int varId, vtkDoubleArray *coords)
vtkSmartPointer< vtkDoubleArray > LongitudeCoordinates
vtkDependentDimensionInfo(int ncFD, int varId, vtkNetCDFCFReader *parent)
vtkSmartPointer< vtkStringArray > SpecialVariables
int LoadMetaData(int ncFD, int varId, vtkNetCDFCFReader *parent)
vtkSmartPointer< vtkDoubleArray > LatitudeCoordinates
vtkSmartPointer< vtkDoubleArray > GetLongitudeCoordinates() const
int LoadCoordinateVariable(int ncFD, int varId, vtkDoubleArray *coords)
int LoadBoundsVariable(int ncFD, int varId, vtkDoubleArray *coords)
vtkSmartPointer< vtkIntArray > GetGridDimensions() const
vtkSmartPointer< vtkDoubleArray > Bounds
vtkSmartPointer< vtkStringArray > SpecialVariables
vtkSmartPointer< vtkDoubleArray > Coordinates
vtkSmartPointer< vtkDoubleArray > GetCoordinates()
vtkSmartPointer< vtkDoubleArray > GetBounds()
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()
virtual int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
virtual void SetOutputType(int type)
vtkDimensionInfoVector * DimensionInfo
void AddUnstructuredSphericalCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
void FakeStructuredCoordinates(vtkStructuredGrid *structuredOutput)
static vtkNetCDFCFReader * New()
void Add2DSphericalCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
void Add1DSphericalCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
virtual void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
virtual int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
static int CanReadFile(const char *filename)
Returns true if the given file can be read.
void Add2DSphericalCoordinates(vtkStructuredGrid *structuredOutput)
void Add1DRectilinearCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
void SetOutputTypeToRectilinear()
virtual int IsTimeDimension(int ncFD, int dimId)
Determines whether the given variable is a time dimension.
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.
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...
vtkDependentDimensionInfo * FindDependentDimensionInfo(vtkIntArray *dims)
void Add1DRectilinearCoordinates(vtkStructuredGrid *structuredOutput)
void SetOutputTypeToStructured()
virtual bool DimensionsAreForPointData(vtkIntArray *dimensions)
Returns false for spherical dimensions, which should use cell data.
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
void Add1DSphericalCoordinates(vtkStructuredGrid *structuredOutput)
void Add2DRectilinearCoordinates(vtkStructuredGrid *structuredOutput)
void Add2DSphericalCoordinates(vtkPoints *points, const int extent[6])
void AddRectilinearCoordinates(vtkRectilinearGrid *rectilinearOutput)
void Add2DRectilinearCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
virtual void GetUpdateExtentForOutput(vtkDataSet *output, int extent[6])
Overridden to retrieve stored extent for unstructured data.
void FakeRectilinearCoordinates(vtkRectilinearGrid *rectilinearOutput)
void Add2DRectilinearCoordinates(vtkPoints *points, const int extent[6])
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])
void AddStructuredCells(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal method for building unstructred cells that match structured cells.
vtkDimensionInfo * GetDimensionInfo(int dimension)
virtual int ReadMetaData(int ncFD)
Interprets the special conventions of COARDS.
void AddUnstructuredRectilinearCoordinates(vtkUnstructuredGrid *unstructuredOutput, const int extent[6])
Internal methods for creating unstructured cells.
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:40
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:49
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types
@ points
Definition: vtkX3D.h:446
@ extent
Definition: vtkX3D.h:345
@ type
Definition: vtkX3D.h:516
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
#define VTK_IMAGE_DATA
Definition: vtkType.h:93
#define VTK_RECTILINEAR_GRID
Definition: vtkType.h:90
#define VTK_UNSTRUCTURED_GRID
Definition: vtkType.h:91
#define VTK_STRUCTURED_GRID
Definition: vtkType.h:89