VTK
vtkNetCDFReader.h
Go to the documentation of this file.
1// -*- c++ -*-
2/*=========================================================================
3
4 Program: Visualization Toolkit
5 Module: vtkNetCDFReader.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
34#ifndef vtkNetCDFReader_h
35#define vtkNetCDFReader_h
36
37#include "vtkIONetCDFModule.h" // For export macro
39
40#include "vtkSmartPointer.h" // For ivars
41#include <string> //For std::string
42
44class vtkDataSet;
45class vtkDoubleArray;
46class vtkIntArray;
47class vtkStdString;
48class vtkStringArray;
49class vtkNetCDFReaderPrivate;
50
51class VTKIONETCDF_EXPORT vtkNetCDFReader : public vtkDataObjectAlgorithm
52{
53public:
56 virtual void PrintSelf(ostream &os, vtkIndent indent);
57
58 virtual void SetFileName(const char *filename);
60
66
67// // Description:
68// // Get the data array selection tables used to configure which variables to
69// // load.
70// vtkGetObjectMacro(VariableArraySelection, vtkDataArraySelection);
71
73
77 virtual const char *GetVariableArrayName(int idx);
78 virtual int GetVariableArrayStatus(const char *name);
79 virtual void SetVariableArrayStatus(const char *name, int status);
81
88
90
96 vtkGetObjectMacro(VariableDimensions, vtkStringArray);
98
106 virtual void SetDimensions(const char *dimensions);
107
109
116 vtkGetObjectMacro(AllDimensions, vtkStringArray);
118
120
129 vtkGetMacro(ReplaceFillValueWithNan, int);
130 vtkSetMacro(ReplaceFillValueWithNan, int);
131 vtkBooleanMacro(ReplaceFillValueWithNan, int);
133
135
143
147 std::string QueryArrayUnits(const char *ArrayName);
148
149protected:
152
153 char *FileName;
156
161
163
165
170
175
177
178 int WholeExtent[6];
179
180 virtual int RequestDataObject(vtkInformation *request,
181 vtkInformationVector **inputVector,
182 vtkInformationVector *outputVector);
183
184 virtual int RequestInformation(vtkInformation *request,
185 vtkInformationVector **inputVector,
186 vtkInformationVector *outputVector);
187
188 virtual int RequestData(vtkInformation *request,
189 vtkInformationVector **inputVector,
190 vtkInformationVector *outputVector);
191
195 static void SelectionModifiedCallback(vtkObject *caller, unsigned long eid,
196 void *clientdata, void *calldata);
197
202 vtkStdString DescribeDimensions(int ncFD, const int *dimIds, int numDims);
203
207 virtual int ReadMetaData(int ncFD);
208
212 virtual int FillVariableDimensions(int ncFD);
213
221 virtual int IsTimeDimension(int ncFD, int dimId);
222
230 virtual vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId);
231
238 virtual bool DimensionsAreForPointData(vtkIntArray *vtkNotUsed(dimensions)) {
239 return true;
240 }
241
248 virtual void GetUpdateExtentForOutput(vtkDataSet *output, int extent[6]);
249
254 virtual int LoadVariable(int ncFD, const char *varName, double time,
255 vtkDataSet *output);
256
257private:
258 vtkNetCDFReader(const vtkNetCDFReader &) VTK_DELETE_FUNCTION;
259 void operator=(const vtkNetCDFReader &) VTK_DELETE_FUNCTION;
260
261 int UpdateExtent[6];
262 char *TimeUnits;
263 char *Calendar;
264 vtkNetCDFReaderPrivate *Private;
265};
266
267#endif //vtkNetCDFReader_h
Store on/off settings for data arrays for a vtkSource.
Superclass for algorithms that produce only data object as output.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
dynamic, self-adjusting array of double
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
A superclass for reading netCDF files.
virtual void SetVariableArrayStatus(const char *name, int status)
static void SelectionModifiedCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
Callback registered with the VariableArraySelection.
int UpdateMetaData()
Update the meta data from the current file.
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
virtual int LoadVariable(int ncFD, const char *varName, double time, vtkDataSet *output)
Load the variable at the given time into the given data set.
std::string QueryArrayUnits(const char *ArrayName)
Get units attached to a particular array in the netcdf file.
vtkSmartPointer< vtkDataArraySelection > VariableArraySelection
virtual void GetUpdateExtentForOutput(vtkDataSet *output, int extent[6])
Retrieves the update extent for the output object.
virtual vtkStringArray * GetAllVariableArrayNames()
Convenience method to get a list of variable arrays.
virtual int GetVariableArrayStatus(const char *name)
virtual int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
vtkSmartPointer< vtkStringArray > AllVariableArrayNames
vtkTimeStamp MetaDataMTime
vtkStdString DescribeDimensions(int ncFD, const int *dimIds, int numDims)
Convenience function for getting a string that describes a set of dimensions.
virtual int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
virtual void SetFileName(const char *filename)
vtkTimeStamp FileNameMTime
virtual bool DimensionsAreForPointData(vtkIntArray *vtkNotUsed(dimensions))
Called internally to determine whether a variable with the given set of dimensions should be loaded a...
virtual int FillVariableDimensions(int ncFD)
Fills the VariableDimensions array.
vtkSmartPointer< vtkIntArray > LoadingDimensions
The dimension ids of the arrays being loaded into the data.
virtual int GetNumberOfVariableArrays()
Variable array selection.
vtkStringArray * VariableDimensions
Placeholder for structure returned from GetVariableDimensions().
virtual int IsTimeDimension(int ncFD, int dimId)
Determines whether the given variable is a time dimension.
virtual int ReadMetaData(int ncFD)
Reads meta data and populates ivars.
virtual const char * GetVariableArrayName(int idx)
static vtkNetCDFReader * New()
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...
virtual void SetDimensions(const char *dimensions)
Loads the grid with the given dimensions.
vtkStringArray * AllDimensions
Placeholder for structure returned from GetAllDimensions().
virtual void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
abstract base class for most VTK objects
Definition: vtkObject.h:60
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:49
a vtkAbstractArray subclass for strings
record modification and/or execution time
Definition: vtkTimeStamp.h:36
@ time
Definition: vtkX3D.h:497
@ extent
Definition: vtkX3D.h:345
@ name
Definition: vtkX3D.h:219
@ string
Definition: vtkX3D.h:490
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkGetStringMacro(ExtensionsString)
Returns a string listing all available extensions.