VTK  9.2.6
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 =========================================================================*/
142 #ifndef vtkRectilinearGrid_h
143 #define vtkRectilinearGrid_h
144 
145 #include "vtkCommonDataModelModule.h" // For export macro
146 #include "vtkDataSet.h"
147 #include "vtkStructuredData.h" // For inline methods
148 
149 class vtkVertex;
150 class vtkLine;
151 class vtkPixel;
152 class vtkVoxel;
153 class vtkDataArray;
154 class vtkPoints;
155 
156 class VTKCOMMONDATAMODEL_EXPORT vtkRectilinearGrid : public vtkDataSet
157 {
158 public:
161 
163  void PrintSelf(ostream& os, vtkIndent indent) override;
164 
168  int GetDataObjectType() override { return VTK_RECTILINEAR_GRID; }
169 
174  void CopyStructure(vtkDataSet* ds) override;
175 
179  void Initialize() override;
180 
182 
185  vtkIdType GetNumberOfCells() override;
186  vtkIdType GetNumberOfPoints() override;
187  double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override;
188  void GetPoint(vtkIdType id, double x[3]) override;
189  vtkCell* GetCell(vtkIdType cellId) override;
190  vtkCell* GetCell(int i, int j, int k) override;
191  void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
192  void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
193  vtkIdType FindPoint(double x, double y, double z) { return this->vtkDataSet::FindPoint(x, y, z); }
194  vtkIdType FindPoint(double x[3]) override;
195  vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
196  double pcoords[3], double* weights) override;
197  vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
198  double tol2, int& subId, double pcoords[3], double* weights) override;
199  vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
200  double pcoords[3], double* weights) override;
201  int GetCellType(vtkIdType cellId) override;
202  vtkIdType GetCellSize(vtkIdType cellId) override;
204  void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override
205  {
206  vtkStructuredData::GetCellPoints(cellId, ptIds, this->DataDescription, this->Dimensions);
207  }
208  void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
209  {
210  vtkStructuredData::GetPointCells(ptId, cellIds, this->Dimensions);
211  }
212  void ComputeBounds() override;
213  int GetMaxCellSize() override { return 8; } // voxel is the largest
214  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
215  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
217 
223  unsigned char IsPointVisible(vtkIdType ptId);
224 
230  unsigned char IsCellVisible(vtkIdType cellId);
231 
236  bool HasAnyBlankPoints() override;
241  bool HasAnyBlankCells() override;
242 
249  void GetCellDims(int cellDims[3]);
250 
255  void GetPoints(vtkPoints* pnts);
256 
258 
262  void SetDimensions(int i, int j, int k);
263  void SetDimensions(const int dim[3]);
265 
267 
270  vtkGetVectorMacro(Dimensions, int, 3);
272 
276  int GetDataDimension();
277 
284  int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3]);
285 
289  vtkIdType ComputePointId(int ijk[3]);
290 
294  vtkIdType ComputeCellId(int ijk[3]);
295 
301  void GetPoint(const int i, const int j, const int k, double p[3]);
302 
304 
308  vtkGetObjectMacro(XCoordinates, vtkDataArray);
310 
312 
316  vtkGetObjectMacro(YCoordinates, vtkDataArray);
318 
320 
324  vtkGetObjectMacro(ZCoordinates, vtkDataArray);
326 
328 
333  void SetExtent(int extent[6]);
334  void SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax);
335  vtkGetVector6Macro(Extent, int);
337 
346  unsigned long GetActualMemorySize() override;
347 
349 
352  void ShallowCopy(vtkDataObject* src) override;
353  void DeepCopy(vtkDataObject* src) override;
355 
359  int GetExtentType() override { return VTK_3D_EXTENT; }
360 
366  void Crop(const int* updateExtent) override;
367 
369 
375 
377 
380  static void SetScalarType(int, vtkInformation* meta_data);
381  static int GetScalarType(vtkInformation* meta_data);
382  static bool HasScalarType(vtkInformation* meta_data);
384  const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
386 
388 
392  static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
397 
398 protected:
401 
402  // for the GetCell method
407 
408  int Dimensions[3];
410 
411  int Extent[6];
412 
416 
417  // Hang on to some space for returning points when GetPoint(id) is called.
418  double PointReturn[3];
419 
420 private:
421  void Cleanup();
422 
423 private:
424  vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
425  void operator=(const vtkRectilinearGrid&) = delete;
426 };
427 
428 //----------------------------------------------------------------------------
430 {
431  vtkIdType nCells = 1;
432  int i;
433 
434  for (i = 0; i < 3; i++)
435  {
436  if (this->Dimensions[i] <= 0)
437  {
438  return 0;
439  }
440  if (this->Dimensions[i] > 1)
441  {
442  nCells *= (this->Dimensions[i] - 1);
443  }
444  }
445 
446  return nCells;
447 }
448 
449 //----------------------------------------------------------------------------
451 {
452  return static_cast<vtkIdType>(this->Dimensions[0]) * this->Dimensions[1] * this->Dimensions[2];
453 }
454 
455 //----------------------------------------------------------------------------
457 {
459 }
460 
461 //----------------------------------------------------------------------------
463 {
465 }
466 
467 //----------------------------------------------------------------------------
469 {
470  return vtkStructuredData::ComputeCellId(this->Dimensions, ijk);
471 }
472 
473 #endif
abstract class to specify cell behavior
Definition: vtkCell.h:150
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
general representation of visualization data
abstract class to specify dataset behavior
Definition: vtkDataSet.h:172
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:341
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:143
a simple class to control print indentation
Definition: vtkIndent.h:119
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
cell represents a 1D line
Definition: vtkLine.h:143
a cell that represents an orthogonal quadrilateral
Definition: vtkPixel.h:77
represent and manipulate 3D points
Definition: vtkPoints.h:149
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 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, 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:103
a cell that represents a 3D orthogonal parallelepiped
Definition: vtkVoxel.h:91
@ 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(...)