VTK  9.3.20240918
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
130#ifndef vtkRectilinearGrid_h
131#define vtkRectilinearGrid_h
132
133#include "vtkCommonDataModelModule.h" // For export macro
134#include "vtkDataSet.h"
135#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_3_0
136#include "vtkSmartPointer.h" // For vtkSmartPointer
137#include "vtkStructuredData.h" // For inline methods
138
139VTK_ABI_NAMESPACE_BEGIN
140class vtkDataArray;
142class vtkPoints;
143
144class VTKCOMMONDATAMODEL_EXPORT vtkRectilinearGrid : public vtkDataSet
145{
146public:
149
151 void PrintSelf(ostream& os, vtkIndent indent) override;
152
156 int GetDataObjectType() override { return VTK_RECTILINEAR_GRID; }
157
162 void CopyStructure(vtkDataSet* ds) override;
163
167 void Initialize() override;
168
170
173 vtkIdType GetNumberOfCells() override;
174 vtkIdType GetNumberOfPoints() override;
175 vtkPoints* GetPoints() override;
176 double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override;
177 void GetPoint(vtkIdType id, double x[3]) override;
178 vtkCell* GetCell(vtkIdType cellId) override;
179 vtkCell* GetCell(int i, int j, int k) override;
180 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
181 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
182 vtkIdType FindPoint(double x[3]) override;
183 vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
184 double pcoords[3], double* weights) override;
185 vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
186 double tol2, int& subId, double pcoords[3], double* weights) override;
187 vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
188 double pcoords[3], double* weights) override;
189 int GetCellType(vtkIdType cellId) override;
191 void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
192 VTK_SIZEHINT(pts, npts) override;
193 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
194 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
195 {
196 vtkStructuredData::GetPointCells(ptId, cellIds, this->Dimensions);
197 }
198 void ComputeBounds() override;
199 int GetMaxCellSize() override { return 8; } // voxel is the largest
200 int GetMaxSpatialDimension() override;
201 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
202 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
204
211
219
221
227 virtual void BlankPoint(vtkIdType ptId);
228 virtual void UnBlankPoint(vtkIdType ptId);
229 virtual void BlankPoint(int i, int j, int k);
230 virtual void UnBlankPoint(int i, int j, int k);
232
234
240 virtual void BlankCell(vtkIdType ptId);
241 virtual void UnBlankCell(vtkIdType ptId);
242 virtual void BlankCell(int i, int j, int k);
243 virtual void UnBlankCell(int i, int j, int k);
245
251 unsigned char IsPointVisible(vtkIdType ptId);
252
258 unsigned char IsCellVisible(vtkIdType cellId);
259
264 bool HasAnyBlankPoints() override;
269 bool HasAnyBlankCells() override;
270
274 vtkGetMacro(DataDescription, int);
275
282 void GetCellDims(int cellDims[3]);
283
288 VTK_DEPRECATED_IN_9_3_0("Use vtkPoints* GetPoints() instead.")
289 void GetPoints(vtkPoints* points);
290
292
296 void SetDimensions(int i, int j, int k);
297 void SetDimensions(const int dim[3]);
299
301
304 vtkGetVectorMacro(Dimensions, int, 3);
306
310 int GetDataDimension();
311
318 int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3]);
319
323 vtkIdType ComputePointId(int ijk[3]);
324
328 vtkIdType ComputeCellId(int ijk[3]);
329
335 void GetPoint(int i, int j, int k, double p[3]);
336
338
341 virtual void SetXCoordinates(vtkDataArray*);
342 vtkGetObjectMacro(XCoordinates, vtkDataArray);
344
346
349 virtual void SetYCoordinates(vtkDataArray*);
350 vtkGetObjectMacro(YCoordinates, vtkDataArray);
352
354
357 virtual void SetZCoordinates(vtkDataArray*);
358 vtkGetObjectMacro(ZCoordinates, vtkDataArray);
360
362
367 void SetExtent(int extent[6]);
368 void SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax);
369 vtkGetVector6Macro(Extent, int);
371
380 unsigned long GetActualMemorySize() override;
381
383
386 void ShallowCopy(vtkDataObject* src) override;
387 void DeepCopy(vtkDataObject* src) override;
389
393 int GetExtentType() override { return VTK_3D_EXTENT; }
394
400 void Crop(const int* updateExtent) override;
401
403
409
411
414 static void SetScalarType(int, vtkInformation* meta_data);
415 static int GetScalarType(vtkInformation* meta_data);
416 static bool HasScalarType(vtkInformation* meta_data);
418 const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
420
422
426 static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
431
432protected:
435
436 int Dimensions[3];
438
439 int Extent[6];
440
444
445 // Hang on to some space for returning points when GetPoint(id) is called.
446 double Point[3];
447
451
456
457private:
458 void Cleanup();
459
460 vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
461 void operator=(const vtkRectilinearGrid&) = delete;
462};
463
464//----------------------------------------------------------------------------
466{
467 this->GetPoint(id, this->Point);
468 return this->Point;
469}
470
471//----------------------------------------------------------------------------
473{
475}
476
477//----------------------------------------------------------------------------
479{
481}
482
483//----------------------------------------------------------------------------
485{
487}
488
489//----------------------------------------------------------------------------
491{
493}
494
495//----------------------------------------------------------------------------
497{
499}
500
501//----------------------------------------------------------------------------
503{
505}
506
507VTK_ABI_NAMESPACE_END
508#endif
void GetPoint(int i, int j, int k, double pnt[3])
abstract class to specify cell behavior
Definition vtkCell.h:130
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
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:133
A read only array class that wraps an implicit function from integers to any value type supported by ...
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition vtkPoints.h:139
a dataset that is topologically regular with variable spacing in the three coordinate directions
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.
vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
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 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.
vtkConstantArray< int > * GetCellTypesArray()
Get the array of all cell types in the rectilinear grid.
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.
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.
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 GetMaxCellSize() override
Standard vtkDataSet API methods.
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 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(...)