VTK  9.3.0
vtkNetCDFUGRIDReader.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
3 
4 #ifndef vtkNetCDFUGRIDReader_h
5 #define vtkNetCDFUGRIDReader_h
6 
7 #include "vtkDataArraySelection.h" // For vtkSmartPointer downcast
8 #include "vtkIONetCDFModule.h" // For export macro
9 #include "vtkSmartPointer.h" // For vtkSmartPointer
11 
12 #include <cstdlib> // For std::size_t
13 #include <vector> // For std::vector
14 
15 VTK_ABI_NAMESPACE_BEGIN
16 
27 class VTKIONETCDF_EXPORT vtkNetCDFUGRIDReader : public vtkUnstructuredGridAlgorithm
28 {
29 public:
32  void PrintSelf(ostream& os, vtkIndent indent) override;
33 
35 
41 
43 
47  vtkGetObjectMacro(PointDataArraySelection, vtkDataArraySelection);
48  vtkGetObjectMacro(CellDataArraySelection, vtkDataArraySelection);
50 
52 
58 
60 
63  const char* GetPointArrayName(int index);
64  const char* GetCellArrayName(int index);
66 
68 
71  int GetPointArrayStatus(const char* name);
72  int GetCellArrayStatus(const char* name);
73  void SetPointArrayStatus(const char* name, int status);
74  void SetCellArrayStatus(const char* name, int status);
76 
78 
87  vtkGetMacro(ReplaceFillValueWithNan, bool);
88  vtkSetMacro(ReplaceFillValueWithNan, bool);
89  vtkBooleanMacro(ReplaceFillValueWithNan, bool);
91 
92 protected:
95 
99  vtkInformation*, vtkInformationVector**, vtkInformationVector* outputVector) override;
100 
101  bool Open();
102  bool ParseHeader();
103  bool FillArraySelection(const std::vector<int>& ids, vtkDataArraySelection* selection);
106  bool FillArrays(vtkUnstructuredGrid* output, std::size_t timeStep);
107  void Close();
108 
109  bool CheckError(int error);
112  std::string GetAttributeName(int var, int att);
114  vtkSmartPointer<vtkDataArray> GetArrayData(int var, std::size_t time, std::size_t size);
115 
116 private:
117  char* FileName = nullptr;
118 
119  int NcId = -1;
120  int MeshVarId = -1;
121  int FaceVarId = -1;
122  int FaceFillValue = -1;
123  int FaceStartIndex = 0;
124  int NodeXVarId = -1;
125  int NodeYVarId = -1;
126  int NodeType = -1;
127  std::size_t NodeCount = 0;
128  std::size_t FaceCount = 0;
129  std::size_t NodesPerFace = 0;
130  std::size_t FaceStride = 0;
131  std::size_t NodesPerFaceStride = 0;
132  bool ReplaceFillValueWithNan = false;
133  std::vector<int> NodeArrayVarIds; // data variables linked to nodes (points)
134  std::vector<int> FaceArrayVarIds; // data variables linked to face (cells)
135  std::vector<double> TimeSteps;
136 
137  vtkSmartPointer<vtkDataArraySelection> PointDataArraySelection;
138  vtkSmartPointer<vtkDataArraySelection> CellDataArraySelection;
139 
141  void operator=(const vtkNetCDFUGRIDReader&) = delete;
142 };
143 
144 VTK_ABI_NAMESPACE_END
145 
146 #endif
Store on/off settings for data arrays, etc.
a simple class to control print indentation
Definition: vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
vtkSetFilePathMacro(FileName)
Get/Set the file name of the file.
const char * GetPointArrayName(int index)
Get the name of the point or cell with the given index in the input.
int GetNumberOfCellArrays()
Get the number of point or cell arrays available in the input.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
bool FillArrays(vtkUnstructuredGrid *output, std::size_t timeStep)
bool FillArraySelection(const std::vector< int > &ids, vtkDataArraySelection *selection)
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *outputVector) override
This is called by the superclass.
int GetCellArrayStatus(const char *name)
Get/Set whether the point or cell with the given name is to be read.
const char * GetCellArrayName(int index)
Get the name of the point or cell with the given index in the input.
void SetPointArrayStatus(const char *name, int status)
Get/Set whether the point or cell with the given name is to be read.
vtkGetFilePathMacro(FileName)
Get/Set the file name of the file.
int GetPointArrayStatus(const char *name)
Get/Set whether the point or cell with the given name is to be read.
~vtkNetCDFUGRIDReader() override
int GetNumberOfPointArrays()
Get the number of point or cell arrays available in the input.
void SetCellArrayStatus(const char *name, int status)
Get/Set whether the point or cell with the given name is to be read.
bool CheckError(int error)
std::string GetVariableName(int var)
bool FillCells(vtkUnstructuredGrid *output)
static vtkNetCDFUGRIDReader * New()
std::string GetDimensionName(int dim)
std::string GetAttributeName(int var, int att)
bool FillPoints(vtkUnstructuredGrid *output)
std::string GetAttributeString(int var, std::string name)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
vtkSmartPointer< vtkDataArray > GetArrayData(int var, std::size_t time, std::size_t size)
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
@ time
Definition: vtkX3D.h:497
@ name
Definition: vtkX3D.h:219
@ size
Definition: vtkX3D.h:253
@ index
Definition: vtkX3D.h:246
@ string
Definition: vtkX3D.h:490