VTK  9.3.20240418
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 
19 #include "vtkDataObjectAlgorithm.h"
20 #include "vtkIONetCDFModule.h" // For export macro
21 
22 #include "vtkSmartPointer.h" // For ivars
23 #include <string> //For std::string
24 
25 VTK_ABI_NAMESPACE_BEGIN
27 class vtkDataSet;
28 class vtkDoubleArray;
29 class vtkIntArray;
30 class vtkStdString;
31 class vtkStringArray;
32 class vtkNetCDFReaderPrivate;
33 
34 class VTKIONETCDF_EXPORT vtkNetCDFReader : public vtkDataObjectAlgorithm
35 {
36 public:
38  static vtkNetCDFReader* New();
39  void PrintSelf(ostream& os, vtkIndent indent) override;
40 
41  virtual void SetFileName(VTK_FILEPATH const char* filename);
43 
49 
50  // // Description:
51  // // Get the data array selection tables used to configure which variables to
52  // // load.
53  // vtkGetObjectMacro(VariableArraySelection, vtkDataArraySelection);
54 
56 
60  virtual const char* GetVariableArrayName(int index);
61  virtual int GetVariableArrayStatus(const char* name);
62  virtual void SetVariableArrayStatus(const char* name, int status);
64 
71 
73 
79  vtkGetObjectMacro(VariableDimensions, vtkStringArray);
81 
89  virtual void SetDimensions(const char* dimensions);
90 
96 
98 
105  vtkGetObjectMacro(AllDimensions, vtkStringArray);
107 
109 
118  vtkGetMacro(ReplaceFillValueWithNan, vtkTypeBool);
119  vtkSetMacro(ReplaceFillValueWithNan, vtkTypeBool);
120  vtkBooleanMacro(ReplaceFillValueWithNan, vtkTypeBool);
122 
124 
129  vtkGetStringMacro(TimeUnits);
130  vtkGetStringMacro(Calendar);
132 
136  std::string QueryArrayUnits(const char* ArrayName);
137 
138 protected:
140  ~vtkNetCDFReader() override;
141 
142  char* FileName;
145 
150 
152 
154 
159 
161 
166 
168 
169  int WholeExtent[6];
170 
172  vtkInformationVector* outputVector) override;
173 
175  vtkInformationVector* outputVector) override;
176 
177  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
178  vtkInformationVector* outputVector) override;
179 
184  vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
185 
190  vtkStdString DescribeDimensions(int ncFD, const int* dimIds, int numDims);
191 
195  virtual int ReadMetaData(int ncFD);
196 
200  virtual int FillVariableDimensions(int ncFD);
201 
209  virtual int IsTimeDimension(int ncFD, int dimId);
210 
218  virtual vtkSmartPointer<vtkDoubleArray> GetTimeValues(int ncFD, int dimId);
219 
226  virtual bool DimensionsAreForPointData(vtkIntArray* vtkNotUsed(dimensions)) { return true; }
227 
234  virtual void GetUpdateExtentForOutput(vtkDataSet* output, int extent[6]);
235 
240  virtual int LoadVariable(int ncFD, const char* varName, double time, vtkDataSet* output);
241 
242 private:
243  vtkNetCDFReader(const vtkNetCDFReader&) = delete;
244  void operator=(const vtkNetCDFReader&) = delete;
245 
246  int UpdateExtent[6];
247  char* TimeUnits;
248  char* Calendar;
249  vtkNetCDFReaderPrivate* Private;
250 };
251 
252 VTK_ABI_NAMESPACE_END
253 #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:166
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
Definition: vtkIntArray.h:144
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 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 int GetVariableArrayStatus(const char *name)
Variable array selection.
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...
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.
virtual const char * GetVariableArrayName(int index)
Variable array selection.
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.
virtual vtkStringArray * GetAllVariableArrayNames()
Convenience method to get a list of variable arrays.
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()
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
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:78
a vtkAbstractArray subclass for strings
record modification and/or execution time
Definition: vtkTimeStamp.h:44
@ time
Definition: vtkX3D.h:497
@ extent
Definition: vtkX3D.h:345
@ name
Definition: vtkX3D.h:219
@ index
Definition: vtkX3D.h:246
@ string
Definition: vtkX3D.h:490
int vtkTypeBool
Definition: vtkABI.h:64
#define VTK_FILEPATH