VTK  9.4.20241226
vtkNetCDFReader.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
16#ifndef vtkNetCDFReader_h
17#define vtkNetCDFReader_h
18
20#include "vtkIONetCDFModule.h" // For export macro
21#include "vtkNetCDFAccessor.h" // For netcdf/xarray accessor
22
23#include "vtkSmartPointer.h" // For ivars
24#include <string> //For std::string
25
26VTK_ABI_NAMESPACE_BEGIN
28class vtkDataSet;
29class vtkDoubleArray;
30class vtkIntArray;
31class vtkStdString;
32class vtkStringArray;
33class vtkNetCDFReaderPrivate;
34
35class VTKIONETCDF_EXPORT vtkNetCDFReader : public vtkDataObjectAlgorithm
36{
37public:
40 void PrintSelf(ostream& os, vtkIndent indent) override;
41
42 virtual void SetFileName(VTK_FILEPATH const char* filename);
45 vtkSetObjectMacro(Accessor, vtkNetCDFAccessor);
46 vtkGetObjectMacro(Accessor, vtkNetCDFAccessor);
48
54
55 // // Description:
56 // // Get the data array selection tables used to configure which variables to
57 // // load.
58 // vtkGetObjectMacro(VariableArraySelection, vtkDataArraySelection);
59
61
65 virtual const char* GetVariableArrayName(int index);
66 virtual int GetVariableArrayStatus(const char* name);
67 virtual void SetVariableArrayStatus(const char* name, int status);
69
76
78
84 vtkGetObjectMacro(VariableDimensions, vtkStringArray);
86
94 virtual void SetDimensions(const char* dimensions);
95
101
103
110 vtkGetObjectMacro(AllDimensions, vtkStringArray);
112
114
123 vtkGetMacro(ReplaceFillValueWithNan, vtkTypeBool);
124 vtkSetMacro(ReplaceFillValueWithNan, vtkTypeBool);
125 vtkBooleanMacro(ReplaceFillValueWithNan, vtkTypeBool);
127
129
134 vtkGetStringMacro(TimeUnits);
135 vtkGetStringMacro(Calendar);
137
141 std::string QueryArrayUnits(const char* ArrayName);
142
143protected:
146
147 char* FileName;
150
155
157
159
164
165 std::string CurrentDimensions;
166
171
173
174 int WholeExtent[6];
175
177 vtkInformationVector* outputVector) override;
178
180 vtkInformationVector* outputVector) override;
181
183 vtkInformationVector* outputVector) override;
184
189 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
190
195 vtkStdString DescribeDimensions(int ncFD, const int* dimIds, int numDims);
196
200 virtual int ReadMetaData(int ncFD);
201
205 virtual int FillVariableDimensions(int ncFD);
206
214 virtual int IsTimeDimension(int ncFD, int dimId);
215
223 virtual vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId);
224
231 virtual bool DimensionsAreForPointData(vtkIntArray* vtkNotUsed(dimensions)) { return true; }
232
239 virtual void GetUpdateExtentForOutput(vtkDataSet* output, int extent[6]);
240
245 virtual int LoadVariable(int ncFD, const char* varName, double time, vtkDataSet* output);
246
248
249private:
250 vtkNetCDFReader(const vtkNetCDFReader&) = delete;
251 void operator=(const vtkNetCDFReader&) = delete;
252
253 int UpdateExtent[6];
254 char* TimeUnits;
255 char* Calendar;
256 vtkNetCDFReaderPrivate* Private;
257};
258
259VTK_ABI_NAMESPACE_END
260#endif // vtkNetCDFReader_h
Store on/off settings for data arrays, etc.
Superclass for algorithms that produce only data object as output.
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
dynamic, self-adjusting array of double
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
A superclass for reading netCDF files.
virtual void SetVariableArrayStatus(const char *name, int status)
Variable array selection.
static void SelectionModifiedCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
Callback registered with the VariableArraySelection.
vtkTypeBool ReplaceFillValueWithNan
int UpdateMetaData()
Update the meta data from the current file.
virtual const char * GetVariableArrayName(int index)
Variable array selection.
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
int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
bool ComputeArraySelection()
Enables arrays in VariableArraySelection depending on Dimensions.
vtkGetFilePathMacro(FileName)
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)
Variable array selection.
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.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkTimeStamp FileNameMTime
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
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.
~vtkNetCDFReader() override
static vtkNetCDFReader * New()
vtkNetCDFAccessor * Accessor
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.
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
vtkStringArray * AllDimensions
Placeholder for structure returned from GetAllDimensions().
std::string CurrentDimensions
virtual void SetFileName(VTK_FILEPATH const char *filename)
abstract base class for most VTK objects
Definition vtkObject.h:162
Hold a reference to a vtkObjectBase instance.
Wrapper around std::string to keep symbols short.
a vtkAbstractArray subclass for strings
record modification and/or execution time
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_FILEPATH