VTK  9.6.1
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 "vtkCartesianGrid.h"
35#include "vtkCommonDataModelModule.h" // For export macro
36#include "vtkDeprecation.h" // for VTK_DEPRECATED_IN_9_6_0
37#include "vtkSmartPointer.h" // For vtkSmartPointer
38#include "vtkStructuredData.h" // For inline methods
39#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
40
41VTK_ABI_NAMESPACE_BEGIN
42class vtkDataArray;
44class vtkPoints;
45
46class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALAUTO vtkRectilinearGrid : public vtkCartesianGrid
47{
48public:
51
53 void PrintSelf(ostream& os, vtkIndent indent) override;
54
58 int GetDataObjectType() VTK_FUTURE_CONST override { return VTK_RECTILINEAR_GRID; }
59
61
64 void ShallowCopy(vtkDataObject* src) override;
65 void DeepCopy(vtkDataObject* src) override;
67
72 void CopyStructure(vtkDataSet* ds) override;
73
77 void Initialize() override;
78
82
85 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
86 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
87 vtkIdType FindPoint(double x[3]) override;
88 vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
89 double pcoords[3], double* weights) override;
90 void ComputeBounds() override;
92
99 int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3]) override;
100
101 VTK_DEPRECATED_IN_9_6_0("Use the const version instead")
102 int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3])
103 {
104 const double* pt = x;
105 return this->ComputeStructuredCoordinates(pt, ijk, pcoords);
106 };
107
113 vtkIdType ComputePointId(int ijk[3]) override;
114
120 vtkIdType ComputeCellId(int ijk[3]) override;
121
122 using Superclass::GetPoint;
128 void GetPoint(int i, int j, int k, double p[3]);
129
131
135 vtkGetObjectMacro(XCoordinates, vtkDataArray);
137
139
143 vtkGetObjectMacro(YCoordinates, vtkDataArray);
145
147
151 vtkGetObjectMacro(ZCoordinates, vtkDataArray);
153
162 unsigned long GetActualMemorySize() override;
163
169 void Crop(const int* updateExtent) override;
170
172
178
179protected:
182
186
187 void BuildPoints() override;
188
189private:
190 void Cleanup();
191
192 vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
193 void operator=(const vtkRectilinearGrid&) = delete;
194};
195
196//----------------------------------------------------------------------------
201
202//----------------------------------------------------------------------------
207
208VTK_ABI_NAMESPACE_END
209#endif
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 * GetDimensions()
Get dimensions of this structured points dataset.
virtual int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3])=0
Computes the structured coordinates for a point x[3].
virtual vtkIdType ComputeCellId(int ijk[3])=0
Given a location in structured coordinates (i-j-k), return the cell id.
virtual vtkIdType ComputePointId(int ijk[3])=0
Given a location in structured coordinates (i-j-k), return the point id.
vtkCell * GetCell(vtkIdType cellId) override
Standard vtkDataSet API methods.
abstract class to specify cell behavior
Definition vtkCell.h:50
provides thread-safe access to cells
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.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard vtkObject API methods.
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.
vtkIdType ComputeCellId(int ijk[3]) override
Given a location in structured coordinates (i-j-k), return the cell id.
int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3]) override
Computes the structured coordinates for a point x[3].
virtual void SetYCoordinates(vtkDataArray *)
Specify the grid coordinates in the y-direction.
virtual void SetXCoordinates(vtkDataArray *)
Specify the grid coordinates in the x-direction.
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
vtkDataArray * XCoordinates
void BuildPoints() override
Pure abstract method responsible to build and set internal points.
void CopyStructure(vtkDataSet *ds) override
Copy the geometric and topological structure of an input rectilinear grid object.
~vtkRectilinearGrid() override
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 FindPoint(double x[3]) override
Standard vtkDataSet API methods.
vtkDataArray * ZCoordinates
void GetCellBounds(vtkIdType cellId, double bounds[6]) override
Standard vtkDataSet API methods.
static vtkRectilinearGrid * New()
void ComputeBounds() override
Standard vtkDataSet API methods.
void Crop(const int *updateExtent) override
Reallocates and copies to set the Extent to the UpdateExtent.
int GetDataObjectType() VTK_FUTURE_CONST override
Return what type of dataset this is.
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
vtkIdType ComputePointId(int ijk[3]) override
Given a location in structured coordinates (i-j-k), return the point id.
static vtkRectilinearGrid * ExtendedNew()
void GetCell(vtkIdType cellId, vtkGenericCell *cell) override
Standard vtkDataSet API methods.
implicit object to represent cell connectivity
static vtkIdType ComputePointId(const int dim[3], const int ijk[3], int dataDescription=vtkStructuredData::VTK_STRUCTURED_EMPTY)
Given a location in structured coordinates (i-j-k), and the dimensions of the structured dataset,...
static vtkIdType ComputeCellId(const int dim[3], const int ijk[3], int dataDescription=vtkStructuredData::VTK_STRUCTURED_EMPTY)
Given a location in structured coordinates (i-j-k), and the dimensions of the structured dataset,...
#define vtkDataArray
#define VTK_DEPRECATED_IN_9_6_0(reason)
int vtkIdType
Definition vtkType.h:368
@ VTK_RECTILINEAR_GRID
Definition vtkType.h:112
#define VTK_MARSHALAUTO