VTK  9.3.1
vtkRectilinearGrid.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
31 #ifndef vtkRectilinearGrid_h
32 #define vtkRectilinearGrid_h
33 
34 #include "vtkCommonDataModelModule.h" // For export macro
35 #include "vtkDataSet.h"
36 #include "vtkStructuredData.h" // For inline methods
37 
38 VTK_ABI_NAMESPACE_BEGIN
39 class vtkVertex;
40 class vtkLine;
41 class vtkPixel;
42 class vtkVoxel;
43 class vtkDataArray;
44 class vtkPoints;
45 
46 class VTKCOMMONDATAMODEL_EXPORT vtkRectilinearGrid : public vtkDataSet
47 {
48 public:
49  static vtkRectilinearGrid* New();
50  static vtkRectilinearGrid* ExtendedNew();
51 
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
58  int GetDataObjectType() override { return VTK_RECTILINEAR_GRID; }
59 
64  void CopyStructure(vtkDataSet* ds) override;
65 
69  void Initialize() override;
70 
72 
75  vtkIdType GetNumberOfCells() override;
76  vtkIdType GetNumberOfPoints() override;
77  double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override;
78  void GetPoint(vtkIdType id, double x[3]) override;
79  vtkCell* GetCell(vtkIdType cellId) override;
80  vtkCell* GetCell(int i, int j, int k) override;
81  void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
82  void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
83  vtkIdType FindPoint(double x, double y, double z) { return this->vtkDataSet::FindPoint(x, y, z); }
84  vtkIdType FindPoint(double x[3]) override;
85  vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
86  double pcoords[3], double* weights) override;
87  vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
88  double tol2, int& subId, double pcoords[3], double* weights) override;
89  vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
90  double pcoords[3], double* weights) override;
91  int GetCellType(vtkIdType cellId) override;
92  vtkIdType GetCellSize(vtkIdType cellId) override;
94  void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override
95  {
96  vtkStructuredData::GetCellPoints(cellId, ptIds, this->DataDescription, this->Dimensions);
97  }
98  void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
99  {
100  vtkStructuredData::GetPointCells(ptId, cellIds, this->Dimensions);
101  }
102  void ComputeBounds() override;
103  int GetMaxCellSize() override { return 8; } // voxel is the largest
104  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
105  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
107 
113  unsigned char IsPointVisible(vtkIdType ptId);
114 
120  unsigned char IsCellVisible(vtkIdType cellId);
121 
126  bool HasAnyBlankPoints() override;
131  bool HasAnyBlankCells() override;
132 
139  void GetCellDims(int cellDims[3]);
140 
145  void GetPoints(vtkPoints* pnts);
146 
148 
152  void SetDimensions(int i, int j, int k);
153  void SetDimensions(const int dim[3]);
155 
157 
160  vtkGetVectorMacro(Dimensions, int, 3);
162 
166  int GetDataDimension();
167 
174  int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3]);
175 
179  vtkIdType ComputePointId(int ijk[3]);
180 
184  vtkIdType ComputeCellId(int ijk[3]);
185 
191  void GetPoint(int i, int j, int k, double p[3]);
192 
194 
197  virtual void SetXCoordinates(vtkDataArray*);
198  vtkGetObjectMacro(XCoordinates, vtkDataArray);
200 
202 
205  virtual void SetYCoordinates(vtkDataArray*);
206  vtkGetObjectMacro(YCoordinates, vtkDataArray);
208 
210 
213  virtual void SetZCoordinates(vtkDataArray*);
214  vtkGetObjectMacro(ZCoordinates, vtkDataArray);
216 
218 
223  void SetExtent(int extent[6]);
224  void SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax);
225  vtkGetVector6Macro(Extent, int);
227 
236  unsigned long GetActualMemorySize() override;
237 
239 
242  void ShallowCopy(vtkDataObject* src) override;
243  void DeepCopy(vtkDataObject* src) override;
245 
249  int GetExtentType() override { return VTK_3D_EXTENT; }
250 
256  void Crop(const int* updateExtent) override;
257 
259 
263  static vtkRectilinearGrid* GetData(vtkInformationVector* v, int i = 0);
265 
267 
270  static void SetScalarType(int, vtkInformation* meta_data);
271  static int GetScalarType(vtkInformation* meta_data);
272  static bool HasScalarType(vtkInformation* meta_data);
273  int GetScalarType();
274  const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
276 
278 
282  static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
283  static int GetNumberOfScalarComponents(vtkInformation* meta_data);
284  static bool HasNumberOfScalarComponents(vtkInformation* meta_data);
285  int GetNumberOfScalarComponents();
287 
288 protected:
290  ~vtkRectilinearGrid() override;
291 
292  // for the GetCell method
297 
298  int Dimensions[3];
300 
301  int Extent[6];
302 
306 
307  // Hang on to some space for returning points when GetPoint(id) is called.
308  double PointReturn[3];
309 
310 private:
311  void Cleanup();
312 
313  vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
314  void operator=(const vtkRectilinearGrid&) = delete;
315 };
316 
317 //----------------------------------------------------------------------------
319 {
320  vtkIdType nCells = 1;
321  int i;
322 
323  for (i = 0; i < 3; i++)
324  {
325  if (this->Dimensions[i] <= 0)
326  {
327  return 0;
328  }
329  if (this->Dimensions[i] > 1)
330  {
331  nCells *= (this->Dimensions[i] - 1);
332  }
333  }
334 
335  return nCells;
336 }
337 
338 //----------------------------------------------------------------------------
340 {
341  return static_cast<vtkIdType>(this->Dimensions[0]) * this->Dimensions[1] * this->Dimensions[2];
342 }
343 
344 //----------------------------------------------------------------------------
346 {
348 }
349 
350 //----------------------------------------------------------------------------
352 {
354 }
355 
356 //----------------------------------------------------------------------------
358 {
359  return vtkStructuredData::ComputeCellId(this->Dimensions, ijk);
360 }
361 
362 VTK_ABI_NAMESPACE_END
363 #endif
int GetDataDimension()
Return the dimensionality of the data.
vtkIdType FindPoint(double x, double y, double z)
Standard vtkDataSet API methods.
int GetDataObjectType() override
Return what type of dataset this is.
a dataset that is topologically regular with variable spacing in the three coordinate directions ...
virtual vtkCell * FindAndGetCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
Locate the cell that contains a point and return the cell.
virtual vtkIdType GetNumberOfCells()=0
Determine the number of cells composing the dataset.
vtkIdType ComputeCellId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the cell id.
int GetMaxCellSize() override
Standard vtkDataSet API methods.
static vtkDataObject * New()
virtual vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)=0
Locate cell based on global coordinate x and tolerance squared.
#define VTK_RECTILINEAR_GRID
Definition: vtkType.h:68
Store vtkAlgorithm input/output information.
virtual vtkIdType GetNumberOfPoints()=0
Determine the number of points composing the dataset.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:52
a cell that represents a 3D point
Definition: vtkVertex.h:22
virtual vtkIdType GetCellSize(vtkIdType cellId)
Get the size of cell with cellId such that: 0 <= cellId < NumberOfCells.
#define VTK_3D_EXTENT
Definition: vtkDataObject.h:51
static int GetDataDimension(int dataDescription)
Return the topological dimension of the data (e.g., 0, 1, 2, or 3D).
a cell that represents an orthogonal quadrilateral
Definition: vtkPixel.h:26
virtual void ComputeBounds()
Compute the data bounding box from data points.
int GetExtentType() override
Structured extent.
vtkDataArray * XCoordinates
int vtkIdType
Definition: vtkType.h:315
void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds) override
Standard vtkDataSet API methods.
virtual bool HasAnyBlankCells()
Returns 1 if there are any blanking cells 0 otherwise.
Definition: vtkDataSet.h:473
static vtkDataSet * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
provides thread-safe access to cells
vtkIdType ComputePointId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the point id.
virtual void Crop(const int *updateExtent)
This method crops the data object (if necessary) so that the extent matches the update extent...
static vtkIdType ComputePointId(const int dim[3], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the dimensions of the structured dataset...
const char * GetScalarTypeAsString()
Set/Get the scalar data type for the points.
cell represents a 1D line
Definition: vtkLine.h:22
abstract class to specify cell behavior
Definition: vtkCell.h:49
a cell that represents a 3D orthogonal parallelepiped
Definition: vtkVoxel.h:30
unsigned long GetActualMemorySize() override
Return the actual size of the data in kibibytes (1024 bytes).
a simple class to control print indentation
Definition: vtkIndent.h:28
vtkIdType GetNumberOfPoints() override
Standard vtkDataSet API methods.
static void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds, int dataDescription, int dim[3])
Get the points defining a cell.
vtkDataArray * ZCoordinates
virtual bool HasAnyBlankPoints()
Returns 1 if there are any blanking points 0 otherwise.
Definition: vtkDataSet.h:479
list of point or cell ids
Definition: vtkIdList.h:22
virtual void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)=0
Topological inquiry to get points defining cell.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:44
void DeepCopy(vtkDataObject *src) override
Shallow and Deep copy.
virtual void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds)
Topological inquiry to get all cells using list of points exclusive of cell specified (e...
#define VTK_SIZEHINT(...)
void Initialize() override
Restore data object to initial state.
void GetPointCells(vtkIdType ptId, vtkIdList *cellIds) override
Standard vtkDataSet API methods.
static vtkIdType ComputeCellId(const int dim[3], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the dimensions of the structured dataset...
vtkDataArray * YCoordinates
vtkIdType GetNumberOfCells() override
Standard vtkDataSet API methods.
static void GetPointCells(vtkIdType ptId, vtkIdList *cellIds, VTK_FUTURE_CONST int dim[3])
Get the cells using a point.
virtual void CopyStructure(vtkDataSet *ds)=0
Copy the geometric and topological structure of an object.
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
Store zero or more vtkInformation instances.
virtual double * GetPoint(vtkIdType ptId)=0
Get point coordinates with ptId such that: 0 <= ptId < NumberOfPoints.
general representation of visualization data
Definition: vtkDataObject.h:54
vtkIdType FindPoint(double x, double y, double z)
Locate the closest point to the global coordinate x.
Definition: vtkDataSet.h:233
virtual vtkCell * GetCell(vtkIdType cellId)=0
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
represent and manipulate 3D points
Definition: vtkPoints.h:28
virtual void GetCellBounds(vtkIdType cellId, double bounds[6])
Get the bounds of the cell with cellId such that: 0 <= cellId < NumberOfCells.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual int GetCellType(vtkIdType cellId)=0
Get type of cell with cellId such that: 0 <= cellId < NumberOfCells.