VTK  9.3.0
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
40 #ifndef vtkRectilinearGrid_h
41 #define vtkRectilinearGrid_h
42 
43 #include "vtkCommonDataModelModule.h" // For export macro
44 #include "vtkDataSet.h"
45 #include "vtkStructuredData.h" // For inline methods
46 
47 VTK_ABI_NAMESPACE_BEGIN
48 class vtkVertex;
49 class vtkLine;
50 class vtkPixel;
51 class vtkVoxel;
52 class vtkDataArray;
53 class vtkPoints;
54 
55 class VTKCOMMONDATAMODEL_EXPORT vtkRectilinearGrid : public vtkDataSet
56 {
57 public:
60 
62  void PrintSelf(ostream& os, vtkIndent indent) override;
63 
67  int GetDataObjectType() override { return VTK_RECTILINEAR_GRID; }
68 
73  void CopyStructure(vtkDataSet* ds) override;
74 
78  void Initialize() override;
79 
81 
84  vtkIdType GetNumberOfCells() override;
85  vtkIdType GetNumberOfPoints() override;
86  double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override;
87  void GetPoint(vtkIdType id, double x[3]) override;
88  vtkCell* GetCell(vtkIdType cellId) override;
89  vtkCell* GetCell(int i, int j, int k) override;
90  void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
91  void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
92  vtkIdType FindPoint(double x, double y, double z) { return this->vtkDataSet::FindPoint(x, y, z); }
93  vtkIdType FindPoint(double x[3]) override;
94  vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
95  double pcoords[3], double* weights) override;
96  vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
97  double tol2, int& subId, double pcoords[3], double* weights) override;
98  vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
99  double pcoords[3], double* weights) override;
100  int GetCellType(vtkIdType cellId) override;
101  vtkIdType GetCellSize(vtkIdType cellId) override;
103  void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override
104  {
105  vtkStructuredData::GetCellPoints(cellId, ptIds, this->DataDescription, this->Dimensions);
106  }
107  void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
108  {
109  vtkStructuredData::GetPointCells(ptId, cellIds, this->Dimensions);
110  }
111  void ComputeBounds() override;
112  int GetMaxCellSize() override { return 8; } // voxel is the largest
113  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
114  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
116 
122  unsigned char IsPointVisible(vtkIdType ptId);
123 
129  unsigned char IsCellVisible(vtkIdType cellId);
130 
135  bool HasAnyBlankPoints() override;
140  bool HasAnyBlankCells() override;
141 
148  void GetCellDims(int cellDims[3]);
149 
154  void GetPoints(vtkPoints* pnts);
155 
157 
161  void SetDimensions(int i, int j, int k);
162  void SetDimensions(const int dim[3]);
164 
166 
169  vtkGetVectorMacro(Dimensions, int, 3);
171 
175  int GetDataDimension();
176 
183  int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3]);
184 
188  vtkIdType ComputePointId(int ijk[3]);
189 
193  vtkIdType ComputeCellId(int ijk[3]);
194 
200  void GetPoint(int i, int j, int k, double p[3]);
201 
203 
207  vtkGetObjectMacro(XCoordinates, vtkDataArray);
209 
211 
215  vtkGetObjectMacro(YCoordinates, vtkDataArray);
217 
219 
223  vtkGetObjectMacro(ZCoordinates, vtkDataArray);
225 
227 
232  void SetExtent(int extent[6]);
233  void SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax);
234  vtkGetVector6Macro(Extent, int);
236 
245  unsigned long GetActualMemorySize() override;
246 
248 
251  void ShallowCopy(vtkDataObject* src) override;
252  void DeepCopy(vtkDataObject* src) override;
254 
258  int GetExtentType() override { return VTK_3D_EXTENT; }
259 
265  void Crop(const int* updateExtent) override;
266 
268 
274 
276 
279  static void SetScalarType(int, vtkInformation* meta_data);
280  static int GetScalarType(vtkInformation* meta_data);
281  static bool HasScalarType(vtkInformation* meta_data);
283  const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
285 
287 
291  static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
296 
297 protected:
300 
301  // for the GetCell method
306 
307  int Dimensions[3];
309 
310  int Extent[6];
311 
315 
316  // Hang on to some space for returning points when GetPoint(id) is called.
317  double PointReturn[3];
318 
319 private:
320  void Cleanup();
321 
322  vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
323  void operator=(const vtkRectilinearGrid&) = delete;
324 };
325 
326 //----------------------------------------------------------------------------
328 {
329  vtkIdType nCells = 1;
330  int i;
331 
332  for (i = 0; i < 3; i++)
333  {
334  if (this->Dimensions[i] <= 0)
335  {
336  return 0;
337  }
338  if (this->Dimensions[i] > 1)
339  {
340  nCells *= (this->Dimensions[i] - 1);
341  }
342  }
343 
344  return nCells;
345 }
346 
347 //----------------------------------------------------------------------------
349 {
350  return static_cast<vtkIdType>(this->Dimensions[0]) * this->Dimensions[1] * this->Dimensions[2];
351 }
352 
353 //----------------------------------------------------------------------------
355 {
357 }
358 
359 //----------------------------------------------------------------------------
361 {
363 }
364 
365 //----------------------------------------------------------------------------
367 {
368  return vtkStructuredData::ComputeCellId(this->Dimensions, ijk);
369 }
370 
371 VTK_ABI_NAMESPACE_END
372 #endif
abstract class to specify cell behavior
Definition: vtkCell.h:59
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
general representation of visualization data
Definition: vtkDataObject.h:64
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
virtual vtkIdType GetNumberOfPoints()=0
Determine the number of points composing the dataset.
virtual vtkIdType GetNumberOfCells()=0
Determine the number of cells composing the dataset.
virtual void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)=0
Topological inquiry to get points defining cell.
vtkIdType FindPoint(double x, double y, double z)
Locate the closest point to the global coordinate x.
Definition: vtkDataSet.h:242
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:32
a simple class to control print indentation
Definition: vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
cell represents a 1D line
Definition: vtkLine.h:32
a cell that represents an orthogonal quadrilateral
Definition: vtkPixel.h:36
represent and manipulate 3D points
Definition: vtkPoints.h:38
a dataset that is topologically regular with variable spacing in the three coordinate directions
static bool HasScalarType(vtkInformation *meta_data)
Set/Get the scalar data type for the points.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Initialize() override
Restore object to initial state.
unsigned long GetActualMemorySize() override
Return the actual size of the data in kibibytes (1024 bytes).
vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
virtual void SetZCoordinates(vtkDataArray *)
Specify the grid coordinates in the z-direction.
void GetPointCells(vtkIdType ptId, vtkIdList *cellIds) override
Standard vtkDataSet API methods.
vtkIdType GetNumberOfPoints() override
Standard vtkDataSet API methods.
bool HasAnyBlankCells() override
Returns 1 if there is any visibility constraint on the cells, 0 otherwise.
unsigned char IsCellVisible(vtkIdType cellId)
Return non-zero value if specified point is visible.
const char * GetScalarTypeAsString()
Set/Get the scalar data type for the points.
bool HasAnyBlankPoints() override
Returns 1 if there is any visibility constraint on the points, 0 otherwise.
virtual void SetYCoordinates(vtkDataArray *)
Specify the grid coordinates in the y-direction.
virtual void SetXCoordinates(vtkDataArray *)
Specify the grid coordinates in the x-direction.
vtkCell * GetCell(vtkIdType cellId) override
Standard vtkDataSet API methods.
vtkIdType ComputeCellId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the cell id.
int GetExtentType() override
Structured extent.
void DeepCopy(vtkDataObject *src) override
Shallow and Deep copy.
vtkDataArray * YCoordinates
vtkIdType GetNumberOfCells() override
Standard vtkDataSet API methods.
void SetExtent(int extent[6])
Different ways to set the extent of the data array.
vtkIdType FindPoint(double x, double y, double z)
Standard vtkDataSet API methods.
int GetDataObjectType() override
Return what type of dataset this is.
static void SetNumberOfScalarComponents(int n, vtkInformation *meta_data)
Set/Get the number of scalar components for points.
static int GetNumberOfScalarComponents(vtkInformation *meta_data)
Set/Get the number of scalar components for points.
vtkDataArray * XCoordinates
static vtkRectilinearGrid * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
void CopyStructure(vtkDataSet *ds) override
Copy the geometric and topological structure of an input rectilinear grid object.
static int GetScalarType(vtkInformation *meta_data)
Set/Get the scalar data type for the points.
void SetDimensions(int i, int j, int k)
Set dimensions of rectilinear grid dataset.
void GetPoints(vtkPoints *pnts)
Given a user-supplied vtkPoints container object, this method fills in all the points of the Rectilin...
static void SetScalarType(int, vtkInformation *meta_data)
Set/Get the scalar data type for the points.
~vtkRectilinearGrid() override
void GetCellDims(int cellDims[3])
Given the node dimensions of this grid instance, this method computes the node dimensions.
unsigned char IsPointVisible(vtkIdType ptId)
Return non-zero value if specified point is visible.
vtkIdType ComputePointId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the point id.
void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds) override
Standard vtkDataSet API methods.
static bool HasNumberOfScalarComponents(vtkInformation *meta_data)
Set/Get the number of scalar components for points.
void GetPoint(int i, int j, int k, double p[3])
Given the IJK-coordinates of the point, it returns the corresponding xyz-coordinates.
vtkIdType GetCellSize(vtkIdType cellId) override
Standard vtkDataSet API methods.
vtkIdType FindPoint(double x[3]) override
Standard vtkDataSet API methods.
int GetNumberOfScalarComponents()
Set/Get the number of scalar components for points.
double * GetPoint(vtkIdType ptId) override
Standard vtkDataSet API methods.
int GetCellType(vtkIdType cellId) override
Standard vtkDataSet API methods.
void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds, int *seedLoc)
Standard vtkDataSet API methods.
void GetPoint(vtkIdType id, double x[3]) override
Standard vtkDataSet API methods.
vtkDataArray * ZCoordinates
static vtkRectilinearGrid * ExtendedNew()
void GetCellBounds(vtkIdType cellId, double bounds[6]) override
Standard vtkDataSet API methods.
static vtkRectilinearGrid * New()
vtkCell * GetCell(int i, int j, int k) override
Standard vtkDataSet API methods.
int GetScalarType()
Set/Get the scalar data type for the points.
int GetDataDimension()
Return the dimensionality of the data.
void ComputeBounds() override
Standard vtkDataSet API methods.
void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds) override
Standard vtkDataSet API methods.
void Crop(const int *updateExtent) override
Reallocates and copies to set the Extent to the UpdateExtent.
vtkIdType FindCell(double x[3], vtkCell *cell, vtkGenericCell *gencell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3])
Convenience function computes the structured coordinates for a point x[3].
static vtkRectilinearGrid * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
int GetMaxCellSize() override
Standard vtkDataSet API methods.
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
void SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax)
Different ways to set the extent of the data array.
vtkCell * FindAndGetCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
void GetCell(vtkIdType cellId, vtkGenericCell *cell) override
Standard vtkDataSet API methods.
void SetDimensions(const int dim[3])
Set dimensions of rectilinear grid dataset.
static void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds, int dataDescription, int dim[3])
Get the points defining a cell.
static int GetDataDimension(int dataDescription)
Return the topological dimension of the data (e.g., 0, 1, 2, or 3D).
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,...
static void GetPointCells(vtkIdType ptId, vtkIdList *cellIds, VTK_FUTURE_CONST int dim[3])
Get the cells using a point.
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,...
a cell that represents a 3D point
Definition: vtkVertex.h:32
a cell that represents a 3D orthogonal parallelepiped
Definition: vtkVoxel.h:40
@ info
Definition: vtkX3D.h:376
@ extent
Definition: vtkX3D.h:345
#define VTK_3D_EXTENT
Definition: vtkDataObject.h:60
int vtkIdType
Definition: vtkType.h:315
#define VTK_RECTILINEAR_GRID
Definition: vtkType.h:68
#define VTK_SIZEHINT(...)