VTK  9.5.20250524
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
134#ifndef vtkImageData_h
135#define vtkImageData_h
136
137#include "vtkCommonDataModelModule.h" // For export macro
138#include "vtkDataSet.h"
139#include "vtkSmartPointer.h" // For vtkSmartPointer ivars
140#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
141
142#include "vtkStructuredData.h" // Needed for inline methods
143
144VTK_ABI_NAMESPACE_BEGIN
145class vtkDataArray;
147class vtkLine;
148class vtkMatrix3x3;
149class vtkMatrix4x4;
150class vtkPixel;
151class vtkPoints;
152class vtkVertex;
153class vtkVoxel;
154
155class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALAUTO vtkImageData : public vtkDataSet
156{
157public:
158 static vtkImageData* New();
160
161 vtkTypeMacro(vtkImageData, vtkDataSet);
162 void PrintSelf(ostream& os, vtkIndent indent) override;
163
168 void CopyStructure(vtkDataSet* ds) override;
169
173 void Initialize() override;
174
178 int GetDataObjectType() VTK_FUTURE_CONST override { return VTK_IMAGE_DATA; }
179
181
188 vtkIdType GetNumberOfCells() override;
189 vtkIdType GetNumberOfPoints() override;
190 vtkPoints* GetPoints() override;
191 double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override;
192 void GetPoint(vtkIdType id, double x[3]) override;
193 vtkCell* GetCell(vtkIdType cellId) override;
194 vtkCell* GetCell(int i, int j, int k) override;
195 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
196 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
198 vtkIdType FindPoint(double x[3]) override;
199 vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
200 double pcoords[3], double* weights) override;
201 vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
202 double tol2, int& subId, double pcoords[3], double* weights) override;
203 vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
204 double pcoords[3], double* weights) override;
205 int GetCellType(vtkIdType cellId) override;
207 void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
208 VTK_SIZEHINT(pts, npts) override;
209 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
210 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
211 {
212 int dimensions[3];
213 this->GetDimensions(dimensions);
214 vtkStructuredData::GetPointCells(ptId, cellIds, dimensions);
215 }
216 void ComputeBounds() override;
217 int GetMaxCellSize() override { return 8; } // voxel is the largest
218 int GetMaxSpatialDimension() override;
219 int GetMinSpatialDimension() override;
220 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
222
229
237
245 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
246
248
254 virtual void BlankPoint(vtkIdType ptId);
255 virtual void UnBlankPoint(vtkIdType ptId);
256 virtual void BlankPoint(int i, int j, int k);
257 virtual void UnBlankPoint(int i, int j, int k);
259
261
267 virtual void BlankCell(vtkIdType ptId);
268 virtual void UnBlankCell(vtkIdType ptId);
269 virtual void BlankCell(int i, int j, int k);
270 virtual void UnBlankCell(int i, int j, int k);
272
278 unsigned char IsPointVisible(vtkIdType ptId);
279
285 unsigned char IsCellVisible(vtkIdType cellId);
286
291 bool HasAnyBlankPoints() override;
292
297 bool HasAnyBlankCells() override;
298
302 vtkGetMacro(DataDescription, int);
303
310 void GetCellDims(int cellDims[3]);
311
315 virtual void SetDimensions(int i, int j, int k);
316
320 virtual void SetDimensions(const int dims[3]);
321
328 virtual int* GetDimensions() VTK_SIZEHINT(3);
329
336 virtual void GetDimensions(int dims[3]);
337#if VTK_ID_TYPE_IMPL != VTK_INT
338 virtual void GetDimensions(vtkIdType dims[3]);
339#endif
340
342
349 virtual int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3]);
351 const double x[3], int ijk[3], double pcoords[3], double tol2);
353
363 virtual void GetVoxelGradient(int i, int j, int k, vtkDataArray* s, vtkDataArray* g);
364
371 virtual void GetPointGradient(int i, int j, int k, vtkDataArray* s, double g[3]);
372
376 virtual int GetDataDimension();
377
381 virtual vtkIdType ComputePointId(int ijk[3])
382 {
383 return vtkStructuredData::ComputePointIdForExtent(this->Extent, ijk);
384 }
385
389 virtual vtkIdType ComputeCellId(int ijk[3])
390 {
391 return vtkStructuredData::ComputeCellIdForExtent(this->Extent, ijk);
392 }
393
395
399 int axis, int min, int max, const int* updateExtent, int* axisUpdateExtent);
400 virtual void GetAxisUpdateExtent(int axis, int& min, int& max, const int* updateExtent);
402
404
415 virtual void SetExtent(int extent[6]);
416 virtual void SetExtent(int x1, int x2, int y1, int y2, int z1, int z2);
417 vtkGetVector6Macro(Extent, int);
419
421
425 virtual double GetScalarTypeMin(vtkInformation* meta_data);
426 virtual double GetScalarTypeMin();
427 virtual double GetScalarTypeMax(vtkInformation* meta_data);
428 virtual double GetScalarTypeMax();
430
432
435 virtual int GetScalarSize(vtkInformation* meta_data);
436 virtual int GetScalarSize();
438
440
452 virtual void GetIncrements(vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
453 virtual void GetIncrements(vtkIdType inc[3]);
454 virtual vtkIdType* GetIncrements(vtkDataArray* scalars) VTK_SIZEHINT(3);
455 virtual void GetIncrements(
456 vtkDataArray* scalars, vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
457 virtual void GetIncrements(vtkDataArray* scalars, vtkIdType inc[3]);
459
461
474 virtual void GetContinuousIncrements(
475 int extent[6], vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
476 virtual void GetContinuousIncrements(
477 vtkDataArray* scalars, int extent[6], vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
479
481
484 virtual void* GetScalarPointerForExtent(int extent[6]);
485 virtual void* GetScalarPointer(int coordinates[3]);
486 virtual void* GetScalarPointer(int x, int y, int z);
487 virtual void* GetScalarPointer();
489
491
494 virtual vtkIdType GetScalarIndexForExtent(int extent[6]);
495 virtual vtkIdType GetScalarIndex(int coordinates[3]);
496 virtual vtkIdType GetScalarIndex(int x, int y, int z);
498
500
503 virtual float GetScalarComponentAsFloat(int x, int y, int z, int component);
504 virtual void SetScalarComponentFromFloat(int x, int y, int z, int component, float v);
505 virtual double GetScalarComponentAsDouble(int x, int y, int z, int component);
506 virtual void SetScalarComponentFromDouble(int x, int y, int z, int component, double v);
508
514 virtual void AllocateScalars(int dataType, int numComponents);
515
522 virtual void AllocateScalars(vtkInformation* pipeline_info);
523
525
531 virtual void CopyAndCastFrom(vtkImageData* inData, int extent[6]);
532 virtual void CopyAndCastFrom(vtkImageData* inData, int x0, int x1, int y0, int y1, int z0, int z1)
533 {
534 int e[6];
535 e[0] = x0;
536 e[1] = x1;
537 e[2] = y0;
538 e[3] = y1;
539 e[4] = z0;
540 e[5] = z1;
541 this->CopyAndCastFrom(inData, e);
542 }
544
550 void Crop(const int* updateExtent) override;
551
560 unsigned long GetActualMemorySize() override;
561
563
567 vtkGetVector3Macro(Spacing, double);
568 virtual void SetSpacing(double i, double j, double k);
569 virtual void SetSpacing(const double ijk[3]);
571
573
581 vtkGetVector3Macro(Origin, double);
582 virtual void SetOrigin(double i, double j, double k);
583 virtual void SetOrigin(const double ijk[3]);
585
587
591 vtkGetObjectMacro(DirectionMatrix, vtkMatrix3x3);
593 virtual void SetDirectionMatrix(const double elements[9]);
594 virtual void SetDirectionMatrix(double e00, double e01, double e02, double e10, double e11,
595 double e12, double e20, double e21, double e22);
597
599
603 vtkGetObjectMacro(IndexToPhysicalMatrix, vtkMatrix4x4);
605
607
618
620
623 virtual void TransformContinuousIndexToPhysicalPoint(double i, double j, double k, double xyz[3]);
624 virtual void TransformContinuousIndexToPhysicalPoint(const double ijk[3], double xyz[3]);
625 virtual void TransformIndexToPhysicalPoint(int i, int j, int k, double xyz[3]);
626 virtual void TransformIndexToPhysicalPoint(const int ijk[3], double xyz[3]);
627 static void TransformContinuousIndexToPhysicalPoint(double i, double j, double k,
628 double const origin[3], double const spacing[3], double const direction[9], double xyz[3]);
630
632
636 vtkGetObjectMacro(PhysicalToIndexMatrix, vtkMatrix4x4);
638
640
651
653
656 virtual void TransformPhysicalPointToContinuousIndex(double x, double y, double z, double ijk[3]);
657 virtual void TransformPhysicalPointToContinuousIndex(const double xyz[3], double ijk[3]);
659
661
664 virtual void TransformPhysicalNormalToContinuousIndex(const double xyz[3], double ijk[3]);
666
671 virtual void TransformPhysicalPlaneToContinuousIndex(double const pplane[4], double iplane[4]);
672
674
678 double const origin[3], double const spacing[3], double const direction[9], double result[16]);
679
681
685 double const origin[3], double const spacing[3], double const direction[9], double result[16]);
687
688 static void SetScalarType(int, vtkInformation* meta_data);
689 static int GetScalarType(vtkInformation* meta_data);
690 static bool HasScalarType(vtkInformation* meta_data);
692 const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
693
695
699 static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
704
709 void CopyInformationFromPipeline(vtkInformation* information) override;
710
716 void CopyInformationToPipeline(vtkInformation* information) override;
717
723 void PrepareForNewData() override;
724
726
729 void ShallowCopy(vtkDataObject* src) override;
730 void DeepCopy(vtkDataObject* src) override;
732
733 //--------------------------------------------------------------------------
734 // Methods that apply to any array (not just scalars).
735 // I am starting to experiment with generalizing imaging filters
736 // to operate on more than just scalars.
737
739
744 void* GetArrayPointerForExtent(vtkDataArray* array, int extent[6]);
745 void* GetArrayPointer(vtkDataArray* array, int coordinates[3]);
747
749
756 vtkIdType GetTupleIndex(vtkDataArray* array, int coordinates[3]);
758
763 void GetArrayIncrements(vtkDataArray* array, vtkIdType increments[3]);
764
771 void ComputeInternalExtent(int* intExt, int* tgtExt, int* bnds);
772
776 int GetExtentType() VTK_FUTURE_CONST override { return VTK_3D_EXTENT; }
777
779
785
786protected:
788 ~vtkImageData() override;
789
790 // The extent of what is currently in the structured grid.
791 // Dimensions is just an array to return a value.
792 // Its contents are out of data until GetDimensions is called.
793 int Dimensions[3];
794 vtkIdType Increments[3];
795
796 // Variables used to define dataset physical orientation
797 double Origin[3];
798 double Spacing[3];
802
803 int Extent[6];
804
808
809 // The first method assumes Active Scalars
810 void ComputeIncrements();
811 // This one is given the number of components of the
812 // scalar field explicitly
813 void ComputeIncrements(int numberOfComponents);
814 void ComputeIncrements(vtkDataArray* scalars);
815
816 // The first method assumes Active Scalars
818 // This one is given the number of components of the
819 // scalar field explicitly
820 void ComputeIncrements(int numberOfComponents, vtkIdType inc[3]);
822
823 // for the index to physical methods
825
830
831private:
832 void InternalImageDataCopy(vtkImageData* src);
833
834 friend class vtkUniformGrid;
835
836 // for the GetPoint method
837 double Point[3];
838
839 int DataDescription;
840 bool DirectionMatrixIsIdentity;
841
842 vtkImageData(const vtkImageData&) = delete;
843 void operator=(const vtkImageData&) = delete;
844};
845
846//----------------------------------------------------------------------------
848{
849 this->ComputeIncrements(this->Increments);
850}
851
852//----------------------------------------------------------------------------
853inline void vtkImageData::ComputeIncrements(int numberOfComponents)
854{
855 this->ComputeIncrements(numberOfComponents, this->Increments);
856}
857
858//----------------------------------------------------------------------------
860{
861 this->ComputeIncrements(scalars, this->Increments);
862}
863
864//----------------------------------------------------------------------------
866{
867 this->GetPoint(id, this->Point);
868 return this->Point;
869}
870
871//----------------------------------------------------------------------------
873{
875}
876
877//----------------------------------------------------------------------------
879{
881}
882
883//----------------------------------------------------------------------------
885{
886 return vtkStructuredData::GetDataDimension(this->DataDescription);
887}
888
889//----------------------------------------------------------------------------
891{
892 return vtkStructuredData::GetDataDimension(this->DataDescription);
893}
894
895//----------------------------------------------------------------------------
897{
898 return vtkStructuredData::GetDataDimension(this->DataDescription);
899}
900VTK_ABI_NAMESPACE_END
901#endif
abstract class to specify cell behavior
Definition vtkCell.h:130
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
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 int GetMaxSpatialDimension()
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
virtual double * GetPoint(vtkIdType ptId)=0
Get point coordinates with ptId such that: 0 <= ptId < NumberOfPoints.
virtual int GetMinSpatialDimension()
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
vtkIdType FindPoint(double x, double y, double z)
Locate the closest point to the global coordinate x.
Definition vtkDataSet.h:352
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:133
topologically and geometrically regular array of data
void Crop(const int *updateExtent) override
Reallocates and copies to set the Extent to updateExtent.
int GetExtentType() VTK_FUTURE_CONST override
The extent type is a 3D extent.
virtual void TransformPhysicalPointToContinuousIndex(const double xyz[3], double ijk[3])
Convert coordinates from physical space (xyz) to index space (ijk).
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.
int GetDataObjectType() VTK_FUTURE_CONST override
Return what type of dataset this is.
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.
vtkPoints * GetPoints() override
Standard vtkDataSet API methods.
vtkIdType GetCellSize(vtkIdType cellId) override
Standard vtkDataSet API methods.
vtkSmartPointer< vtkPoints > StructuredPoints
virtual void UnBlankCell(int i, int j, int k)
Methods for supporting blanking of cells.
void GetCellPoints(vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts, vtkIdList *ptIds) override
Standard vtkDataSet API methods.
int GetMinSpatialDimension() override
Standard vtkDataSet API methods.
virtual void UnBlankCell(vtkIdType ptId)
Methods for supporting blanking of cells.
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 BuildImplicitStructures()
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.
vtkStructuredCellArray * GetCells()
Return the image data connectivity array.
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()
void BuildPoints()
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 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.
void BuildCellTypes()
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 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.
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...
int GetMaxSpatialDimension() override
Standard vtkDataSet API methods.
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.
vtkConstantArray< int > * GetCellTypesArray()
Get the array of all cell types in the image data.
vtkIdType GetNumberOfPoints() override
Standard vtkDataSet API methods.
vtkSmartPointer< vtkStructuredCellArray > StructuredCells
vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
static void ComputePhysicalToIndexMatrix(double const origin[3], double const spacing[3], double const direction[9], double result[16])
Static method to compute the PhysicalToIndexMatrix.
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.
virtual void BlankPoint(int i, int j, int k)
Methods for supporting blanking of cells.
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.
void ApplyIndexToPhysicalMatrix(vtkMatrix4x4 *source)
Set the transformation matrix from the index space to the physical space coordinate system of the dat...
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.
virtual void UnBlankPoint(vtkIdType ptId)
Methods for supporting blanking of cells.
int GetScalarType()
virtual void TransformPhysicalPointToContinuousIndex(double x, double y, double z, double ijk[3])
Convert coordinates from physical space (xyz) to index space (ijk).
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 ApplyPhysicalToIndexMatrix(vtkMatrix4x4 *source)
Get the transformation matrix from the physical space to the index space coordinate system of the dat...
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 object to initial state.
virtual void BlankCell(vtkIdType ptId)
Methods for supporting blanking of cells.
static void ComputeIndexToPhysicalMatrix(double const origin[3], double const spacing[3], double const direction[9], double result[16])
Static method to compute the IndexToPhysicalMatrix.
virtual void UnBlankPoint(int i, int j, int k)
Methods for supporting blanking of cells.
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.
virtual void BlankCell(int i, int j, int k)
Methods for supporting blanking of cells.
vtkSmartPointer< vtkConstantArray< int > > StructuredCellTypes
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 BlankPoint(vtkIdType ptId)
Methods for supporting blanking of cells.
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)
void BuildCells()
virtual int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3], double tol2)
Convenience function computes the structured coordinates for a point x[3].
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 read only array class that wraps an implicit function from integers to any value type supported by ...
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
cell represents a 1D line
Definition vtkLine.h:132
represent and manipulate 3x3 transformation matrices
represent and manipulate 4x4 transformation matrices
a cell that represents an orthogonal quadrilateral
Definition vtkPixel.h:66
represent and manipulate 3D points
Definition vtkPoints.h:139
Hold a reference to a vtkObjectBase instance.
implicit object to represent cell connectivity
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 vtkIdType GetNumberOfCells(const int ext[6], int dataDescription=VTK_EMPTY)
Given the grid extent, this method returns the total number of cells within the extent.
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.
static vtkIdType GetNumberOfPoints(const int ext[6], int dataDescription=VTK_EMPTY)
Given the grid extent, this method returns the total number of points within the extent.
image data with blanking
a cell that represents a 3D point
Definition vtkVertex.h:92
a cell that represents a 3D orthogonal parallelepiped
Definition vtkVoxel.h:80
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
#define VTK_3D_EXTENT
int vtkIdType
Definition vtkType.h:332
@ VTK_IMAGE_DATA
Definition vtkType.h:82
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO
#define max(a, b)