VTK  9.5.20250906
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 "vtkCartesianGrid.h"
134#include "vtkCommonDataModelModule.h" // For export macro
135#include "vtkDeprecation.h" // for VTK_DEPRECATED_IN_9_6_0
136#include "vtkSmartPointer.h" // For vtkSmartPointer
137#include "vtkStructuredData.h" // For inline methods
138#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
139
140VTK_ABI_NAMESPACE_BEGIN
141class vtkDataArray;
143class vtkPoints;
144
145class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALAUTO vtkRectilinearGrid : public vtkCartesianGrid
146{
147public:
150
152 void PrintSelf(ostream& os, vtkIndent indent) override;
153
157 int GetDataObjectType() VTK_FUTURE_CONST override { return VTK_RECTILINEAR_GRID; }
158
160
163 void ShallowCopy(vtkDataObject* src) override;
164 void DeepCopy(vtkDataObject* src) override;
166
171 void CopyStructure(vtkDataSet* ds) override;
172
176 void Initialize() override;
177
178 using Superclass::FindCell;
179 using Superclass::GetCell;
181
184 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
185 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
186 vtkIdType FindPoint(double x[3]) override;
187 vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
188 double pcoords[3], double* weights) override;
189 void ComputeBounds() override;
191
198 int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3]) override;
199
200 VTK_DEPRECATED_IN_9_6_0("Use the const version instead")
201 int ComputeStructuredCoordinates(double x[3], int ijk[3], double pcoords[3])
202 {
203 const double* pt = x;
204 return this->ComputeStructuredCoordinates(pt, ijk, pcoords);
205 };
206
212 vtkIdType ComputePointId(int ijk[3]) override;
213
219 vtkIdType ComputeCellId(int ijk[3]) override;
220
221 using Superclass::GetPoint;
227 void GetPoint(int i, int j, int k, double p[3]);
228
230
234 vtkGetObjectMacro(XCoordinates, vtkDataArray);
236
238
242 vtkGetObjectMacro(YCoordinates, vtkDataArray);
244
246
250 vtkGetObjectMacro(ZCoordinates, vtkDataArray);
252
261 unsigned long GetActualMemorySize() override;
262
268 void Crop(const int* updateExtent) override;
269
271
277
278protected:
281
285
286 void BuildPoints() override;
287
288private:
289 void Cleanup();
290
291 vtkRectilinearGrid(const vtkRectilinearGrid&) = delete;
292 void operator=(const vtkRectilinearGrid&) = delete;
293};
294
295//----------------------------------------------------------------------------
297{
299}
300
301//----------------------------------------------------------------------------
303{
305}
306
307VTK_ABI_NAMESPACE_END
308#endif
Abstract API for vtkImageData and vtkRectilinearGrid.
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.
abstract class to specify cell behavior
Definition vtkCell.h:129
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
provides thread-safe access to cells
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.
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 VTK_DEPRECATED_IN_9_6_0(reason)
int vtkIdType
Definition vtkType.h:332
@ VTK_RECTILINEAR_GRID
Definition vtkType.h:79
#define VTK_MARSHALAUTO