VTK  9.4.20251007
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
30
31#ifndef vtkRectilinearGrid_h
32#define vtkRectilinearGrid_h
33
34#include "vtkCommonDataModelModule.h" // For export macro
35#include "vtkDataSet.h"
36#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_3_0
37#include "vtkSmartPointer.h" // For vtkSmartPointer
38#include "vtkStructuredData.h" // For inline methods
39
40VTK_ABI_NAMESPACE_BEGIN
41class vtkDataArray;
43class vtkPoints;
44
45class VTKCOMMONDATAMODEL_EXPORT vtkRectilinearGrid : public vtkDataSet
46{
47public:
50
52 void PrintSelf(ostream& os, vtkIndent indent) override;
53
57 int GetDataObjectType() override { return VTK_RECTILINEAR_GRID; }
58
63 void CopyStructure(vtkDataSet* ds) override;
64
68 void Initialize() override;
69
71
74 vtkIdType GetNumberOfCells() override;
76 vtkPoints* GetPoints() 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[3]) override;
84 vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
85 double pcoords[3], double* weights) override;
86 vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
87 double tol2, int& subId, double pcoords[3], double* weights) override;
88 vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
89 double pcoords[3], double* weights) override;
90 int GetCellType(vtkIdType cellId) override;
92 void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
93 VTK_SIZEHINT(pts, npts) override;
94 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
95 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
96 {
97 vtkStructuredData::GetPointCells(ptId, cellIds, this->Dimensions);
98 }
99 void ComputeBounds() override;
100 int GetMaxCellSize() override { return 8; } // voxel is the largest
101 int GetMaxSpatialDimension() override;
102 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
103 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
105
112
120
122
128 virtual void BlankPoint(vtkIdType ptId);
129 virtual void UnBlankPoint(vtkIdType ptId);
130 virtual void BlankPoint(int i, int j, int k);
131 virtual void UnBlankPoint(int i, int j, int k);
133
135
141 virtual void BlankCell(vtkIdType ptId);
142 virtual void UnBlankCell(vtkIdType ptId);
143 virtual void BlankCell(int i, int j, int k);
144 virtual void UnBlankCell(int i, int j, int k);
146
152 unsigned char IsPointVisible(vtkIdType ptId);
153
159 unsigned char IsCellVisible(vtkIdType cellId);
160
165 bool HasAnyBlankPoints() override;
170 bool HasAnyBlankCells() override;
171
175 vtkGetMacro(DataDescription, int);
176
183 void GetCellDims(int cellDims[3]);
184
189 VTK_DEPRECATED_IN_9_3_0("Use vtkPoints* GetPoints() instead.")
190 void GetPoints(vtkPoints* points);
191
193
197 void SetDimensions(int i, int j, int k);
198 void SetDimensions(const int dim[3]);
200
202
205 vtkGetVectorMacro(Dimensions, int, 3);
207
211 int GetDataDimension();
212
219 int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3]);
220
224 vtkIdType ComputePointId(int ijk[3]);
225
229 vtkIdType ComputeCellId(int ijk[3]);
230
236 void GetPoint(int i, int j, int k, double p[3]);
237
239
243 vtkGetObjectMacro(XCoordinates, vtkDataArray);
245
247
251 vtkGetObjectMacro(YCoordinates, vtkDataArray);
253
255
259 vtkGetObjectMacro(ZCoordinates, vtkDataArray);
261
263
268 void SetExtent(int extent[6]);
269 void SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax);
270 vtkGetVector6Macro(Extent, int);
272
281 unsigned long GetActualMemorySize() override;
282
284
287 void ShallowCopy(vtkDataObject* src) override;
288 void DeepCopy(vtkDataObject* src) override;
290
294 int GetExtentType() override { return VTK_3D_EXTENT; }
295
301 void Crop(const int* updateExtent) override;
302
304
310
312
315 static void SetScalarType(int, vtkInformation* meta_data);
316 static int GetScalarType(vtkInformation* meta_data);
317 static bool HasScalarType(vtkInformation* meta_data);
319 const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
321
323
327 static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
332
333protected:
336
339
340 int Extent[6];
341
345
346 // Hang on to some space for returning points when GetPoint(id) is called.
347 double Point[3];
348
352
357
358private:
359 void Cleanup();
360
361 vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
362 void operator=(const vtkRectilinearGrid&) = delete;
363};
364
365//----------------------------------------------------------------------------
367{
368 this->GetPoint(id, this->Point);
369 return this->Point;
370}
371
372//----------------------------------------------------------------------------
377
378//----------------------------------------------------------------------------
383
384//----------------------------------------------------------------------------
389
390//----------------------------------------------------------------------------
395
396//----------------------------------------------------------------------------
401
402//----------------------------------------------------------------------------
407
408VTK_ABI_NAMESPACE_END
409#endif
abstract class to specify cell behavior
Definition vtkCell.h:51
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 int GetMaxSpatialDimension()
Get the maximum spatial dimensionality of the data which is the maximum dimension of all cells.
virtual double * GetPoint(vtkIdType ptId)=0
Get point coordinates with ptId such that: 0 <= ptId < NumberOfPoints.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:24
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition vtkPoints.h:30
static vtkRectilinearGrid * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
static bool HasScalarType(vtkInformation *meta_data)
Set/Get the scalar data type for the points.
virtual void UnBlankCell(vtkIdType ptId)
Methods for supporting blanking of cells.
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.
virtual void BlankCell(int i, int j, int k)
Methods for supporting blanking of cells.
virtual void BlankPoint(int i, int j, int k)
Methods for supporting blanking of cells.
virtual void BlankCell(vtkIdType ptId)
Methods for supporting blanking of cells.
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.
virtual void UnBlankCell(int i, int j, int k)
Methods for supporting blanking of cells.
vtkIdType ComputeCellId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the cell id.
int GetExtentType() override
Structured extent.
vtkConstantArray< int > * GetCellTypesArray()
Get the array of all cell types in the rectilinear grid.
void DeepCopy(vtkDataObject *src) override
Shallow and Deep copy.
static vtkRectilinearGrid * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
vtkDataArray * YCoordinates
virtual void UnBlankPoint(int i, int j, int k)
Methods for supporting blanking of cells.
vtkCell * GetCell(vtkIdType cellId) override
Standard vtkDataSet API methods.
vtkIdType GetNumberOfCells() override
Standard vtkDataSet API methods.
void SetExtent(int extent[6])
Different ways to set the extent of the data array.
int GetMaxSpatialDimension() override
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
void GetCellPoints(vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts, vtkIdList *ptIds) override
Standard vtkDataSet API methods.
const char * GetScalarTypeAsString()
Set/Get the scalar data type for the points.
vtkSmartPointer< vtkStructuredCellArray > StructuredCells
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.
vtkSmartPointer< vtkPoints > StructuredPoints
void BuildImplicitStructures()
static void SetScalarType(int, vtkInformation *meta_data)
Set/Get the scalar data type for the points.
~vtkRectilinearGrid() override
virtual void BlankPoint(vtkIdType ptId)
Methods for supporting blanking of cells.
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.
virtual void UnBlankPoint(vtkIdType ptId)
Methods for supporting blanking of cells.
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.
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.
int GetCellType(vtkIdType cellId) override
Standard vtkDataSet API methods.
vtkSmartPointer< vtkConstantArray< int > > StructuredCellTypes
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.
double * GetPoint(vtkIdType ptId) override
Standard vtkDataSet API methods.
vtkDataArray * ZCoordinates
void GetCellBounds(vtkIdType cellId, double bounds[6]) override
Standard vtkDataSet API methods.
vtkCell * GetCell(int i, int j, int k) override
Standard vtkDataSet API methods.
static vtkRectilinearGrid * New()
vtkStructuredCellArray * GetCells()
Return the rectilinear grid connectivity array.
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.
vtkPoints * GetPoints() override
Standard vtkDataSet API methods.
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.
vtkCell * FindAndGetCell(double x[3], vtkCell *cell, 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].
int GetMaxCellSize() override
Standard vtkDataSet API methods.
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
static vtkRectilinearGrid * ExtendedNew()
void GetCell(vtkIdType cellId, vtkGenericCell *cell) override
Standard vtkDataSet API methods.
Hold a reference to a vtkObjectBase instance.
implicit object to represent cell connectivity
static vtkIdType GetNumberOfCells(const int ext[6], int dataDescription=VTK_EMPTY)
Given the grid extent, this method returns the total number of cells within the extent.
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 GetNumberOfPoints(const int ext[6], int dataDescription=VTK_EMPTY)
Given the grid extent, this method returns the total number of points within the extent.
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,...
#define vtkDataArray
vtkImplicitArray< vtkConstantImplicitBackend< T > > vtkConstantArray
A utility alias for wrapping constant functions in implicit arrays.
#define VTK_3D_EXTENT
#define VTK_DEPRECATED_IN_9_3_0(reason)
int vtkIdType
Definition vtkType.h:315
#define VTK_RECTILINEAR_GRID
Definition vtkType.h:68
#define VTK_SIZEHINT(...)