VTK  9.3.0
vtkImageData.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
34#ifndef vtkImageData_h
35#define vtkImageData_h
36
37#include "vtkCommonDataModelModule.h" // For export macro
38#include "vtkDataSet.h"
39
40#include "vtkStructuredData.h" // Needed for inline methods
41
42VTK_ABI_NAMESPACE_BEGIN
43class vtkDataArray;
44class vtkLine;
45class vtkMatrix3x3;
46class vtkMatrix4x4;
47class vtkPixel;
48class vtkVertex;
49class vtkVoxel;
50
51class VTKCOMMONDATAMODEL_EXPORT vtkImageData : public vtkDataSet
52{
53public:
54 static vtkImageData* New();
56
57 vtkTypeMacro(vtkImageData, vtkDataSet);
58 void PrintSelf(ostream& os, vtkIndent indent) override;
59
64 void CopyStructure(vtkDataSet* ds) override;
65
69 int GetDataObjectType() override { return VTK_IMAGE_DATA; }
70
72
81 double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override;
82 void GetPoint(vtkIdType id, double x[3]) override;
83 vtkCell* GetCell(vtkIdType cellId) override;
84 vtkCell* GetCell(int i, int j, int k) override;
85 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
86 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
87 virtual vtkIdType FindPoint(double x, double y, double z)
88 {
89 return this->vtkDataSet::FindPoint(x, y, z);
90 }
91 vtkIdType FindPoint(double x[3]) override;
92 vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
93 double pcoords[3], double* weights) override;
94 vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
95 double tol2, int& subId, double pcoords[3], double* weights) override;
96 vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
97 double pcoords[3], double* weights) override;
98 int GetCellType(vtkIdType cellId) override;
101 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override
102 {
103 int dimensions[3];
104 this->GetDimensions(dimensions);
105 vtkStructuredData::GetCellPoints(cellId, ptIds, this->DataDescription, dimensions);
106 }
107 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
108 {
109 int dimensions[3];
110 this->GetDimensions(dimensions);
111 vtkStructuredData::GetPointCells(ptId, cellIds, dimensions);
112 }
113 void ComputeBounds() override;
114 int GetMaxCellSize() override { return 8; } // voxel is the largest
115 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
117
125 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
126
130 void Initialize() override;
131
137 unsigned char IsPointVisible(vtkIdType ptId);
138
144 unsigned char IsCellVisible(vtkIdType cellId);
145
150 bool HasAnyBlankPoints() override;
155 bool HasAnyBlankCells() override;
156
163 void GetCellDims(int cellDims[3]);
164
168 virtual void SetDimensions(int i, int j, int k);
169
173 virtual void SetDimensions(const int dims[3]);
174
181 virtual int* GetDimensions() VTK_SIZEHINT(3);
182
189 virtual void GetDimensions(int dims[3]);
190#if VTK_ID_TYPE_IMPL != VTK_INT
191 virtual void GetDimensions(vtkIdType dims[3]);
192#endif
193
200 virtual int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3]);
201
211 virtual void GetVoxelGradient(int i, int j, int k, vtkDataArray* s, vtkDataArray* g);
212
219 virtual void GetPointGradient(int i, int j, int k, vtkDataArray* s, double g[3]);
220
224 virtual int GetDataDimension();
225
229 virtual vtkIdType ComputePointId(int ijk[3])
230 {
231 return vtkStructuredData::ComputePointIdForExtent(this->Extent, ijk);
232 }
233
237 virtual vtkIdType ComputeCellId(int ijk[3])
238 {
239 return vtkStructuredData::ComputeCellIdForExtent(this->Extent, ijk);
240 }
241
243
247 int axis, int min, int max, const int* updateExtent, int* axisUpdateExtent);
248 virtual void GetAxisUpdateExtent(int axis, int& min, int& max, const int* updateExtent);
250
252
263 virtual void SetExtent(int extent[6]);
264 virtual void SetExtent(int x1, int x2, int y1, int y2, int z1, int z2);
265 vtkGetVector6Macro(Extent, int);
267
269
273 virtual double GetScalarTypeMin(vtkInformation* meta_data);
274 virtual double GetScalarTypeMin();
275 virtual double GetScalarTypeMax(vtkInformation* meta_data);
276 virtual double GetScalarTypeMax();
278
280
283 virtual int GetScalarSize(vtkInformation* meta_data);
284 virtual int GetScalarSize();
286
288
300 virtual void GetIncrements(vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
301 virtual void GetIncrements(vtkIdType inc[3]);
302 virtual vtkIdType* GetIncrements(vtkDataArray* scalars) VTK_SIZEHINT(3);
303 virtual void GetIncrements(
304 vtkDataArray* scalars, vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
305 virtual void GetIncrements(vtkDataArray* scalars, vtkIdType inc[3]);
307
309
322 virtual void GetContinuousIncrements(
323 int extent[6], vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
324 virtual void GetContinuousIncrements(
325 vtkDataArray* scalars, int extent[6], vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
327
329
332 virtual void* GetScalarPointerForExtent(int extent[6]);
333 virtual void* GetScalarPointer(int coordinates[3]);
334 virtual void* GetScalarPointer(int x, int y, int z);
335 virtual void* GetScalarPointer();
337
339
342 virtual vtkIdType GetScalarIndexForExtent(int extent[6]);
343 virtual vtkIdType GetScalarIndex(int coordinates[3]);
344 virtual vtkIdType GetScalarIndex(int x, int y, int z);
346
348
351 virtual float GetScalarComponentAsFloat(int x, int y, int z, int component);
352 virtual void SetScalarComponentFromFloat(int x, int y, int z, int component, float v);
353 virtual double GetScalarComponentAsDouble(int x, int y, int z, int component);
354 virtual void SetScalarComponentFromDouble(int x, int y, int z, int component, double v);
356
362 virtual void AllocateScalars(int dataType, int numComponents);
363
370 virtual void AllocateScalars(vtkInformation* pipeline_info);
371
373
379 virtual void CopyAndCastFrom(vtkImageData* inData, int extent[6]);
380 virtual void CopyAndCastFrom(vtkImageData* inData, int x0, int x1, int y0, int y1, int z0, int z1)
381 {
382 int e[6];
383 e[0] = x0;
384 e[1] = x1;
385 e[2] = y0;
386 e[3] = y1;
387 e[4] = z0;
388 e[5] = z1;
389 this->CopyAndCastFrom(inData, e);
390 }
392
398 void Crop(const int* updateExtent) override;
399
408 unsigned long GetActualMemorySize() override;
409
411
415 vtkGetVector3Macro(Spacing, double);
416 virtual void SetSpacing(double i, double j, double k);
417 virtual void SetSpacing(const double ijk[3]);
419
421
429 vtkGetVector3Macro(Origin, double);
430 virtual void SetOrigin(double i, double j, double k);
431 virtual void SetOrigin(const double ijk[3]);
433
435
439 vtkGetObjectMacro(DirectionMatrix, vtkMatrix3x3);
441 virtual void SetDirectionMatrix(const double elements[9]);
442 virtual void SetDirectionMatrix(double e00, double e01, double e02, double e10, double e11,
443 double e12, double e20, double e21, double e22);
445
447
451 vtkGetObjectMacro(IndexToPhysicalMatrix, vtkMatrix4x4);
453
455
458 virtual void TransformContinuousIndexToPhysicalPoint(double i, double j, double k, double xyz[3]);
459 virtual void TransformContinuousIndexToPhysicalPoint(const double ijk[3], double xyz[3]);
460 virtual void TransformIndexToPhysicalPoint(int i, int j, int k, double xyz[3]);
461 virtual void TransformIndexToPhysicalPoint(const int ijk[3], double xyz[3]);
462 static void TransformContinuousIndexToPhysicalPoint(double i, double j, double k,
463 double const origin[3], double const spacing[3], double const direction[9], double xyz[3]);
465
467
471 vtkGetObjectMacro(PhysicalToIndexMatrix, vtkMatrix4x4);
473
475
478 virtual void TransformPhysicalPointToContinuousIndex(double x, double y, double z, double ijk[3]);
479 virtual void TransformPhysicalPointToContinuousIndex(const double xyz[3], double ijk[3]);
481
483 double const origin[3], double const spacing[3], double const direction[9], double result[16]);
484
486
489 virtual void TransformPhysicalNormalToContinuousIndex(const double xyz[3], double ijk[3]);
491
496 virtual void TransformPhysicalPlaneToContinuousIndex(double const pplane[4], double iplane[4]);
497
498 static void SetScalarType(int, vtkInformation* meta_data);
499 static int GetScalarType(vtkInformation* meta_data);
500 static bool HasScalarType(vtkInformation* meta_data);
502 const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
503
505
509 static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
514
519 void CopyInformationFromPipeline(vtkInformation* information) override;
520
526 void CopyInformationToPipeline(vtkInformation* information) override;
527
533 void PrepareForNewData() override;
534
536
539 void ShallowCopy(vtkDataObject* src) override;
540 void DeepCopy(vtkDataObject* src) override;
542
543 //--------------------------------------------------------------------------
544 // Methods that apply to any array (not just scalars).
545 // I am starting to experiment with generalizing imaging filters
546 // to operate on more than just scalars.
547
549
554 void* GetArrayPointerForExtent(vtkDataArray* array, int extent[6]);
555 void* GetArrayPointer(vtkDataArray* array, int coordinates[3]);
557
559
566 vtkIdType GetTupleIndex(vtkDataArray* array, int coordinates[3]);
568
573 void GetArrayIncrements(vtkDataArray* array, vtkIdType increments[3]);
574
581 void ComputeInternalExtent(int* intExt, int* tgtExt, int* bnds);
582
586 int GetExtentType() override { return VTK_3D_EXTENT; }
587
589
595
596protected:
598 ~vtkImageData() override;
599
600 // The extent of what is currently in the structured grid.
601 // Dimensions is just an array to return a value.
602 // Its contents are out of data until GetDimensions is called.
603 int Dimensions[3];
604 vtkIdType Increments[3];
605
606 // Variables used to define dataset physical orientation
607 double Origin[3];
608 double Spacing[3];
612
613 int Extent[6];
614
615 // The first method assumes Active Scalars
616 void ComputeIncrements();
617 // This one is given the number of components of the
618 // scalar field explicitly
619 void ComputeIncrements(int numberOfComponents);
620 void ComputeIncrements(vtkDataArray* scalars);
621
622 // The first method assumes Acitive Scalars
624 // This one is given the number of components of the
625 // scalar field explicitly
626 void ComputeIncrements(int numberOfComponents, vtkIdType inc[3]);
628
629 // for the index to physical methods
631
632 // Cell utilities
635 bool GetIJKMinForCellId(vtkIdType cellId, int ijkMin[3]);
636 bool GetIJKMaxForIJKMin(int ijkMin[3], int ijkMax[3]);
637 void AddPointsToCellTemplate(vtkCell* cell, int ijkMin[3], int ijkMax[3]);
638
640
641 void SetDataDescription(int desc);
642 int GetDataDescription() { return this->DataDescription; }
643
644private:
645 void InternalImageDataCopy(vtkImageData* src);
646
647 friend class vtkUniformGrid;
648
649 // for the GetCell method
650 vtkVertex* Vertex;
651 vtkLine* Line;
652 vtkPixel* Pixel;
653 vtkVoxel* Voxel;
654
655 // for the GetPoint method
656 double Point[3];
657
658 int DataDescription;
659
660 vtkImageData(const vtkImageData&) = delete;
661 void operator=(const vtkImageData&) = delete;
662};
663
664//----------------------------------------------------------------------------
666{
667 this->ComputeIncrements(this->Increments);
668}
669
670//----------------------------------------------------------------------------
671inline void vtkImageData::ComputeIncrements(int numberOfComponents)
672{
673 this->ComputeIncrements(numberOfComponents, this->Increments);
674}
675
676//----------------------------------------------------------------------------
678{
679 this->ComputeIncrements(scalars, this->Increments);
680}
681
682//----------------------------------------------------------------------------
684{
685 this->GetPoint(id, this->Point);
686 return this->Point;
687}
688
689//----------------------------------------------------------------------------
691{
692 const int* extent = this->Extent;
693 vtkIdType dims[3];
694 dims[0] = extent[1] - extent[0] + 1;
695 dims[1] = extent[3] - extent[2] + 1;
696 dims[2] = extent[5] - extent[4] + 1;
697
698 return dims[0] * dims[1] * dims[2];
699}
700
701//----------------------------------------------------------------------------
703{
704 return vtkStructuredData::GetDataDimension(this->DataDescription);
705}
706
707VTK_ABI_NAMESPACE_END
708#endif
abstract class to specify cell behavior
Definition vtkCell.h:59
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:62
virtual vtkIdType GetNumberOfPoints()=0
Determine the number of points composing the dataset.
virtual double * GetPoint(vtkIdType ptId)=0
Get point coordinates with ptId such that: 0 <= ptId < NumberOfPoints.
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:242
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:32
topologically and geometrically regular array of data
void Crop(const int *updateExtent) override
Reallocates and copies to set the Extent to updateExtent.
virtual void TransformPhysicalPointToContinuousIndex(const double xyz[3], double ijk[3])
Convert coordinates from physical space (xyz) to index space (ijk).
bool GetCellTemplateForDataDescription(vtkGenericCell *cell)
virtual int GetDataDimension()
Return the dimensionality of the data.
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.
virtual vtkIdType * GetIncrements()
Different ways to get the increments for moving around the data.
void GetArrayIncrements(vtkDataArray *array, vtkIdType increments[3])
Since various arrays have different number of components, the will have different increments.
virtual void TransformContinuousIndexToPhysicalPoint(const double ijk[3], double xyz[3])
Convert coordinates from index space (ijk) to physical space (xyz).
void CopyInformationToPipeline(vtkInformation *information) override
Copy information from this data object to the pipeline information.
bool GetIJKMaxForIJKMin(int ijkMin[3], int ijkMax[3])
vtkIdType GetCellSize(vtkIdType cellId) override
Standard vtkDataSet API methods.
vtkCell * GetCellTemplateForDataDescription()
bool HasAnyBlankCells() override
Returns 1 if there is any visibility constraint on the cells, 0 otherwise.
virtual vtkIdType ComputePointId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the point id.
void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds) override
Standard vtkDataSet API methods.
virtual void SetDirectionMatrix(vtkMatrix3x3 *m)
Set/Get the direction transform of the dataset.
vtkMatrix4x4 * IndexToPhysicalMatrix
virtual int * GetDimensions()
Get dimensions of this structured points dataset.
void ComputeInternalExtent(int *intExt, int *tgtExt, int *bnds)
Given how many pixel are required on a side for boundary conditions (in bnds), the target extent to t...
virtual void SetDimensions(int i, int j, int k)
Same as SetExtent(0, i-1, 0, j-1, 0, k-1)
virtual double GetScalarTypeMin()
These returns the minimum and maximum values the ScalarType can hold without overflowing.
virtual void SetDirectionMatrix(const double elements[9])
Set/Get the direction transform of the dataset.
void ComputeIncrements()
vtkCell * FindAndGetCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
static vtkImageData * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
virtual void TransformIndexToPhysicalPoint(const int ijk[3], double xyz[3])
Convert coordinates from index space (ijk) to physical space (xyz).
static int GetScalarType(vtkInformation *meta_data)
void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds, int *seedLoc)
Get cell neighbors around cell located at seedloc, except cell of id cellId.
void AddPointsToCellTemplate(vtkCell *cell, int ijkMin[3], int ijkMax[3])
void ComputeBounds() override
Standard vtkDataSet API methods.
double * GetPoint(vtkIdType ptId) override
Standard vtkDataSet API methods.
void * GetArrayPointerForExtent(vtkDataArray *array, int extent[6])
These are convenience methods for getting a pointer from any filed array.
vtkCell * GetCell(vtkIdType cellId) override
Standard vtkDataSet API methods.
bool GetIJKMinForCellId(vtkIdType cellId, int ijkMin[3])
static vtkImageData * ExtendedNew()
virtual double GetScalarTypeMin(vtkInformation *meta_data)
These returns the minimum and maximum values the ScalarType can hold without overflowing.
void GetCellBounds(vtkIdType cellId, double bounds[6]) override
Standard vtkDataSet API methods.
virtual void TransformPhysicalNormalToContinuousIndex(const double xyz[3], double ijk[3])
Convert normal from physical space (xyz) to index space (ijk).
void SetDataDescription(int desc)
void GetPoint(vtkIdType id, double x[3]) override
Standard vtkDataSet API methods.
void GetCellDims(int cellDims[3])
Given the node dimensions of this grid instance, this method computes the node dimensions.
vtkTimeStamp ExtentComputeTime
int GetDataDescription()
virtual int GetScalarSize(vtkInformation *meta_data)
Get the size of the scalar type in bytes.
virtual void SetExtent(int extent[6])
Set/Get the extent.
static int GetNumberOfScalarComponents(vtkInformation *meta_data)
Set/Get the number of scalar components for points.
virtual void TransformPhysicalPlaneToContinuousIndex(double const pplane[4], double iplane[4])
Convert a plane from physical to a continuous index.
void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds) override
Standard vtkDataSet API methods.
virtual void GetVoxelGradient(int i, int j, int k, vtkDataArray *s, vtkDataArray *g)
Given structured coordinates (i,j,k) for a voxel cell, compute the eight gradient values for the voxe...
virtual double GetScalarTypeMax(vtkInformation *meta_data)
These returns the minimum and maximum values the ScalarType can hold without overflowing.
static vtkImageData * New()
vtkIdType GetTupleIndex(vtkDataArray *array, int coordinates[3])
Given a data array and a coordinate, return the index of the tuple in the array corresponding to that...
void ComputeIncrements(int numberOfComponents, vtkIdType inc[3])
int GetNumberOfScalarComponents()
Set/Get the number of scalar components for points.
vtkIdType GetNumberOfPoints() override
Standard vtkDataSet API methods.
vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
void ComputeIncrements(vtkDataArray *scalars, vtkIdType inc[3])
virtual void GetAxisUpdateExtent(int axis, int &min, int &max, const int *updateExtent)
Set / Get the extent on just one axis.
void * GetArrayPointer(vtkDataArray *array, int coordinates[3])
These are convenience methods for getting a pointer from any filed array.
int GetExtentType() override
The extent type is a 3D extent.
virtual int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3])
Convenience function computes the structured coordinates for a point x[3].
virtual void SetSpacing(double i, double j, double k)
Set the spacing (width,height,length) of the cubical cells that compose the data set.
static void SetNumberOfScalarComponents(int n, vtkInformation *meta_data)
Set/Get the number of scalar components for points.
virtual void SetDirectionMatrix(double e00, double e01, double e02, double e10, double e11, double e12, double e20, double e21, double e22)
Set/Get the direction transform of the dataset.
vtkIdType Increments[3]
int GetCellType(vtkIdType cellId) override
Standard vtkDataSet API methods.
virtual void SetSpacing(const double ijk[3])
Set the spacing (width,height,length) of the cubical cells that compose the data set.
static void SetScalarType(int, vtkInformation *meta_data)
vtkIdType GetNumberOfCells() override
Standard vtkDataSet API methods.
virtual void SetOrigin(const double ijk[3])
Set/Get the origin of the dataset.
void DeepCopy(vtkDataObject *src) override
Shallow and Deep copy.
virtual int GetScalarSize()
Get the size of the scalar type in bytes.
int GetScalarType()
virtual void TransformPhysicalPointToContinuousIndex(double x, double y, double z, double ijk[3])
Convert coordinates from physical space (xyz) to index space (ijk).
int GetDataObjectType() override
Return what type of dataset this is.
const char * GetScalarTypeAsString()
void CopyInformationFromPipeline(vtkInformation *information) override
Override these to handle origin, spacing, scalar type, and scalar number of components.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static void TransformContinuousIndexToPhysicalPoint(double i, double j, double k, double const origin[3], double const spacing[3], double const direction[9], double xyz[3])
Convert coordinates from index space (ijk) to physical space (xyz).
void ComputeIncrements(vtkIdType inc[3])
void ComputeTransforms()
vtkMatrix3x3 * DirectionMatrix
vtkIdType FindPoint(double x[3]) override
Standard vtkDataSet API methods.
unsigned long GetActualMemorySize() override
Return the actual size of the data in kibibytes (1024 bytes).
void GetCell(vtkIdType cellId, vtkGenericCell *cell) override
Standard vtkDataSet API methods.
void Initialize() override
Restore data object to initial state.
static void ComputeIndexToPhysicalMatrix(double const origin[3], double const spacing[3], double const direction[9], double result[16])
virtual void TransformContinuousIndexToPhysicalPoint(double i, double j, double k, double xyz[3])
Convert coordinates from index space (ijk) to physical space (xyz).
virtual void SetOrigin(double i, double j, double k)
Set/Get the origin of the dataset.
virtual double GetScalarTypeMax()
These returns the minimum and maximum values the ScalarType can hold without overflowing.
virtual void TransformIndexToPhysicalPoint(int i, int j, int k, double xyz[3])
Convert coordinates from index space (ijk) to physical space (xyz).
int GetMaxCellSize() override
Standard vtkDataSet API methods.
unsigned char IsPointVisible(vtkIdType ptId)
Return non-zero value if specified point is visible.
void GetPointCells(vtkIdType ptId, vtkIdList *cellIds) override
Standard vtkDataSet API methods.
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
vtkCell * GetCell(int i, int j, int k) override
Standard vtkDataSet API methods.
virtual void SetAxisUpdateExtent(int axis, int min, int max, const int *updateExtent, int *axisUpdateExtent)
Set / Get the extent on just one axis.
static bool HasNumberOfScalarComponents(vtkInformation *meta_data)
Set/Get the number of scalar components for points.
static bool HasScalarType(vtkInformation *meta_data)
virtual void GetPointGradient(int i, int j, int k, vtkDataArray *s, double g[3])
Given structured coordinates (i,j,k) for a point in a structured point dataset, compute the gradient ...
~vtkImageData() override
unsigned char IsCellVisible(vtkIdType cellId)
Return non-zero value if specified point is visible.
void PrepareForNewData() override
make the output data ready for new data to be inserted.
virtual void SetExtent(int x1, int x2, int y1, int y2, int z1, int z2)
Set/Get the extent.
void CopyStructure(vtkDataSet *ds) override
Copy the geometric and topological structure of an input image data object.
bool HasAnyBlankPoints() override
Returns 1 if there is any visibility constraint on the points, 0 otherwise.
virtual void SetDimensions(const int dims[3])
Same as SetExtent(0, dims[0]-1, 0, dims[1]-1, 0, dims[2]-1)
virtual vtkIdType FindPoint(double x, double y, double z)
Standard vtkDataSet API methods.
vtkMatrix4x4 * PhysicalToIndexMatrix
virtual vtkIdType ComputeCellId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the cell id.
static vtkImageData * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
a simple class to control print indentation
Definition vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
cell represents a 1D line
Definition vtkLine.h:32
represent and manipulate 3x3 transformation matrices
represent and manipulate 4x4 transformation matrices
a cell that represents an orthogonal quadrilateral
Definition vtkPixel.h:36
static vtkIdType ComputePointIdForExtent(const int extent[6], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the extent of the structured 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 ComputeCellIdForExtent(const int extent[6], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the extent of the structured dataset,...
static void GetPointCells(vtkIdType ptId, vtkIdList *cellIds, VTK_FUTURE_CONST int dim[3])
Get the cells using a point.
record modification and/or execution time
image data with blanking
a cell that represents a 3D point
Definition vtkVertex.h:32
a cell that represents a 3D orthogonal parallelepiped
Definition vtkVoxel.h:40
#define VTK_3D_EXTENT
int vtkIdType
Definition vtkType.h:315
#define VTK_IMAGE_DATA
Definition vtkType.h:71
#define VTK_SIZEHINT(...)
#define max(a, b)