VTK  9.1.0
vtkRectilinearGrid.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkRectilinearGrid.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
136 #ifndef vtkRectilinearGrid_h
137 #define vtkRectilinearGrid_h
138 
139 #include "vtkCommonDataModelModule.h" // For export macro
140 #include "vtkDataSet.h"
141 #include "vtkStructuredData.h" // For inline methods
142 
143 class vtkVertex;
144 class vtkLine;
145 class vtkPixel;
146 class vtkVoxel;
147 class vtkDataArray;
148 class vtkPoints;
149 
150 class VTKCOMMONDATAMODEL_EXPORT vtkRectilinearGrid : public vtkDataSet
151 {
152 public:
155 
157  void PrintSelf(ostream& os, vtkIndent indent) override;
158 
162  int GetDataObjectType() override { return VTK_RECTILINEAR_GRID; }
163 
168  void CopyStructure(vtkDataSet* ds) override;
169 
173  void Initialize() override;
174 
176 
179  vtkIdType GetNumberOfCells() override;
180  vtkIdType GetNumberOfPoints() override;
181  double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override;
182  void GetPoint(vtkIdType id, double x[3]) override;
183  vtkCell* GetCell(vtkIdType cellId) override;
184  vtkCell* GetCell(int i, int j, int k) override;
185  void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
186  void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
187  vtkIdType FindPoint(double x, double y, double z) { return this->vtkDataSet::FindPoint(x, y, z); }
188  vtkIdType FindPoint(double x[3]) override;
189  vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
190  double pcoords[3], double* weights) override;
191  vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
192  double tol2, int& subId, double pcoords[3], double* weights) override;
193  vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
194  double pcoords[3], double* weights) override;
195  int GetCellType(vtkIdType cellId) override;
196  void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override
197  {
198  vtkStructuredData::GetCellPoints(cellId, ptIds, this->DataDescription, this->Dimensions);
199  }
200  void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
201  {
202  vtkStructuredData::GetPointCells(ptId, cellIds, this->Dimensions);
203  }
204  void ComputeBounds() override;
205  int GetMaxCellSize() override { return 8; } // voxel is the largest
206  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
207  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
209 
215  unsigned char IsPointVisible(vtkIdType ptId);
216 
222  unsigned char IsCellVisible(vtkIdType cellId);
223 
228  bool HasAnyBlankPoints() override;
233  bool HasAnyBlankCells() override;
234 
241  void GetCellDims(int cellDims[3]);
242 
247  void GetPoints(vtkPoints* pnts);
248 
250 
254  void SetDimensions(int i, int j, int k);
255  void SetDimensions(const int dim[3]);
257 
259 
262  vtkGetVectorMacro(Dimensions, int, 3);
264 
268  int GetDataDimension();
269 
276  int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3]);
277 
281  vtkIdType ComputePointId(int ijk[3]);
282 
286  vtkIdType ComputeCellId(int ijk[3]);
287 
293  void GetPoint(const int i, const int j, const int k, double p[3]);
294 
296 
300  vtkGetObjectMacro(XCoordinates, vtkDataArray);
302 
304 
308  vtkGetObjectMacro(YCoordinates, vtkDataArray);
310 
312 
316  vtkGetObjectMacro(ZCoordinates, vtkDataArray);
318 
320 
325  void SetExtent(int extent[6]);
326  void SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax);
327  vtkGetVector6Macro(Extent, int);
329 
338  unsigned long GetActualMemorySize() override;
339 
341 
344  void ShallowCopy(vtkDataObject* src) override;
345  void DeepCopy(vtkDataObject* src) override;
347 
351  int GetExtentType() override { return VTK_3D_EXTENT; }
352 
358  void Crop(const int* updateExtent) override;
359 
361 
367 
369 
372  static void SetScalarType(int, vtkInformation* meta_data);
373  static int GetScalarType(vtkInformation* meta_data);
374  static bool HasScalarType(vtkInformation* meta_data);
376  const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
378 
380 
384  static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
389 
390 protected:
393 
394  // for the GetCell method
399 
400  int Dimensions[3];
402 
403  int Extent[6];
404 
408 
409  // Hang on to some space for returning points when GetPoint(id) is called.
410  double PointReturn[3];
411 
412 private:
413  void Cleanup();
414 
415 private:
416  vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
417  void operator=(const vtkRectilinearGrid&) = delete;
418 };
419 
420 //----------------------------------------------------------------------------
422 {
423  vtkIdType nCells = 1;
424  int i;
425 
426  for (i = 0; i < 3; i++)
427  {
428  if (this->Dimensions[i] <= 0)
429  {
430  return 0;
431  }
432  if (this->Dimensions[i] > 1)
433  {
434  nCells *= (this->Dimensions[i] - 1);
435  }
436  }
437 
438  return nCells;
439 }
440 
441 //----------------------------------------------------------------------------
443 {
444  return static_cast<vtkIdType>(this->Dimensions[0]) * this->Dimensions[1] * this->Dimensions[2];
445 }
446 
447 //----------------------------------------------------------------------------
449 {
451 }
452 
453 //----------------------------------------------------------------------------
455 {
457 }
458 
459 //----------------------------------------------------------------------------
461 {
462  return vtkStructuredData::ComputeCellId(this->Dimensions, ijk);
463 }
464 
465 #endif
abstract class to specify cell behavior
Definition: vtkCell.h:147
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
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.
vtkIdType FindPoint(double x, double y, double z)
Locate the closest point to the global coordinate x.
Definition: vtkDataSet.h:310
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:140
a simple class to control print indentation
Definition: vtkIndent.h:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
cell represents a 1D line
Definition: vtkLine.h:140
a cell that represents an orthogonal quadrilateral
Definition: vtkPixel.h:74
represent and manipulate 3D points
Definition: vtkPoints.h:143
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.
void GetPoint(const int i, const int j, const int k, double p[3])
Given the IJK-coordinates of the point, it returns the corresponding xyz-coordinates.
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.
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, 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:100
a cell that represents a 3D orthogonal parallelepiped
Definition: vtkVoxel.h:88
@ info
Definition: vtkX3D.h:382
@ extent
Definition: vtkX3D.h:351
#define VTK_3D_EXTENT
int vtkIdType
Definition: vtkType.h:332
#define VTK_RECTILINEAR_GRID
Definition: vtkType.h:80
#define VTK_SIZEHINT(...)