VTK  9.3.20240418
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 
144 VTK_ABI_NAMESPACE_BEGIN
145 class vtkDataArray;
147 class vtkLine;
148 class vtkMatrix3x3;
149 class vtkMatrix4x4;
150 class vtkPixel;
151 class vtkPoints;
152 class vtkVertex;
153 class vtkVoxel;
154 
155 class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALAUTO vtkImageData : public vtkDataSet
156 {
157 public:
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() 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;
197  vtkIdType FindPoint(double x[3]) override;
198  vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
199  double pcoords[3], double* weights) override;
200  vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
201  double tol2, int& subId, double pcoords[3], double* weights) override;
202  vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
203  double pcoords[3], double* weights) override;
204  int GetCellType(vtkIdType cellId) override;
205  vtkIdType GetCellSize(vtkIdType cellId) override;
206  void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
207  VTK_SIZEHINT(pts, npts) override;
208  void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
209  void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override
210  {
211  int dimensions[3];
212  this->GetDimensions(dimensions);
213  vtkStructuredData::GetPointCells(ptId, cellIds, dimensions);
214  }
215  void ComputeBounds() override;
216  int GetMaxCellSize() override { return 8; } // voxel is the largest
217  int GetMaxSpatialDimension() override;
218  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
220 
227 
235 
243  void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc);
244 
246 
252  virtual void BlankPoint(vtkIdType ptId);
253  virtual void UnBlankPoint(vtkIdType ptId);
254  virtual void BlankPoint(int i, int j, int k);
255  virtual void UnBlankPoint(int i, int j, int k);
257 
259 
265  virtual void BlankCell(vtkIdType ptId);
266  virtual void UnBlankCell(vtkIdType ptId);
267  virtual void BlankCell(int i, int j, int k);
268  virtual void UnBlankCell(int i, int j, int k);
270 
276  unsigned char IsPointVisible(vtkIdType ptId);
277 
283  unsigned char IsCellVisible(vtkIdType cellId);
284 
289  bool HasAnyBlankPoints() override;
290 
295  bool HasAnyBlankCells() override;
296 
300  vtkGetMacro(DataDescription, int);
301 
308  void GetCellDims(int cellDims[3]);
309 
313  virtual void SetDimensions(int i, int j, int k);
314 
318  virtual void SetDimensions(const int dims[3]);
319 
326  virtual int* GetDimensions() VTK_SIZEHINT(3);
327 
334  virtual void GetDimensions(int dims[3]);
335 #if VTK_ID_TYPE_IMPL != VTK_INT
336  virtual void GetDimensions(vtkIdType dims[3]);
337 #endif
338 
345  virtual int ComputeStructuredCoordinates(const double x[3], int ijk[3], double pcoords[3]);
346 
356  virtual void GetVoxelGradient(int i, int j, int k, vtkDataArray* s, vtkDataArray* g);
357 
364  virtual void GetPointGradient(int i, int j, int k, vtkDataArray* s, double g[3]);
365 
369  virtual int GetDataDimension();
370 
374  virtual vtkIdType ComputePointId(int ijk[3])
375  {
376  return vtkStructuredData::ComputePointIdForExtent(this->Extent, ijk);
377  }
378 
382  virtual vtkIdType ComputeCellId(int ijk[3])
383  {
384  return vtkStructuredData::ComputeCellIdForExtent(this->Extent, ijk);
385  }
386 
388 
391  virtual void SetAxisUpdateExtent(
392  int axis, int min, int max, const int* updateExtent, int* axisUpdateExtent);
393  virtual void GetAxisUpdateExtent(int axis, int& min, int& max, const int* updateExtent);
395 
397 
408  virtual void SetExtent(int extent[6]);
409  virtual void SetExtent(int x1, int x2, int y1, int y2, int z1, int z2);
410  vtkGetVector6Macro(Extent, int);
412 
414 
418  virtual double GetScalarTypeMin(vtkInformation* meta_data);
419  virtual double GetScalarTypeMin();
420  virtual double GetScalarTypeMax(vtkInformation* meta_data);
421  virtual double GetScalarTypeMax();
423 
425 
428  virtual int GetScalarSize(vtkInformation* meta_data);
429  virtual int GetScalarSize();
431 
433 
445  virtual void GetIncrements(vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
446  virtual void GetIncrements(vtkIdType inc[3]);
447  virtual vtkIdType* GetIncrements(vtkDataArray* scalars) VTK_SIZEHINT(3);
448  virtual void GetIncrements(
449  vtkDataArray* scalars, vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
450  virtual void GetIncrements(vtkDataArray* scalars, vtkIdType inc[3]);
452 
454 
467  virtual void GetContinuousIncrements(
468  int extent[6], vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
469  virtual void GetContinuousIncrements(
470  vtkDataArray* scalars, int extent[6], vtkIdType& incX, vtkIdType& incY, vtkIdType& incZ);
472 
474 
477  virtual void* GetScalarPointerForExtent(int extent[6]);
478  virtual void* GetScalarPointer(int coordinates[3]);
479  virtual void* GetScalarPointer(int x, int y, int z);
480  virtual void* GetScalarPointer();
482 
484 
487  virtual vtkIdType GetScalarIndexForExtent(int extent[6]);
488  virtual vtkIdType GetScalarIndex(int coordinates[3]);
489  virtual vtkIdType GetScalarIndex(int x, int y, int z);
491 
493 
496  virtual float GetScalarComponentAsFloat(int x, int y, int z, int component);
497  virtual void SetScalarComponentFromFloat(int x, int y, int z, int component, float v);
498  virtual double GetScalarComponentAsDouble(int x, int y, int z, int component);
499  virtual void SetScalarComponentFromDouble(int x, int y, int z, int component, double v);
501 
507  virtual void AllocateScalars(int dataType, int numComponents);
508 
515  virtual void AllocateScalars(vtkInformation* pipeline_info);
516 
518 
524  virtual void CopyAndCastFrom(vtkImageData* inData, int extent[6]);
525  virtual void CopyAndCastFrom(vtkImageData* inData, int x0, int x1, int y0, int y1, int z0, int z1)
526  {
527  int e[6];
528  e[0] = x0;
529  e[1] = x1;
530  e[2] = y0;
531  e[3] = y1;
532  e[4] = z0;
533  e[5] = z1;
534  this->CopyAndCastFrom(inData, e);
535  }
537 
543  void Crop(const int* updateExtent) override;
544 
553  unsigned long GetActualMemorySize() override;
554 
556 
560  vtkGetVector3Macro(Spacing, double);
561  virtual void SetSpacing(double i, double j, double k);
562  virtual void SetSpacing(const double ijk[3]);
564 
566 
574  vtkGetVector3Macro(Origin, double);
575  virtual void SetOrigin(double i, double j, double k);
576  virtual void SetOrigin(const double ijk[3]);
578 
580 
584  vtkGetObjectMacro(DirectionMatrix, vtkMatrix3x3);
586  virtual void SetDirectionMatrix(const double elements[9]);
587  virtual void SetDirectionMatrix(double e00, double e01, double e02, double e10, double e11,
588  double e12, double e20, double e21, double e22);
590 
592 
596  vtkGetObjectMacro(IndexToPhysicalMatrix, vtkMatrix4x4);
598 
600 
603  virtual void TransformContinuousIndexToPhysicalPoint(double i, double j, double k, double xyz[3]);
604  virtual void TransformContinuousIndexToPhysicalPoint(const double ijk[3], double xyz[3]);
605  virtual void TransformIndexToPhysicalPoint(int i, int j, int k, double xyz[3]);
606  virtual void TransformIndexToPhysicalPoint(const int ijk[3], double xyz[3]);
607  static void TransformContinuousIndexToPhysicalPoint(double i, double j, double k,
608  double const origin[3], double const spacing[3], double const direction[9], double xyz[3]);
610 
612 
616  vtkGetObjectMacro(PhysicalToIndexMatrix, vtkMatrix4x4);
618 
620 
623  virtual void TransformPhysicalPointToContinuousIndex(double x, double y, double z, double ijk[3]);
624  virtual void TransformPhysicalPointToContinuousIndex(const double xyz[3], double ijk[3]);
626 
628  double const origin[3], double const spacing[3], double const direction[9], double result[16]);
629 
631 
634  virtual void TransformPhysicalNormalToContinuousIndex(const double xyz[3], double ijk[3]);
636 
641  virtual void TransformPhysicalPlaneToContinuousIndex(double const pplane[4], double iplane[4]);
642 
643  static void SetScalarType(int, vtkInformation* meta_data);
644  static int GetScalarType(vtkInformation* meta_data);
645  static bool HasScalarType(vtkInformation* meta_data);
647  const char* GetScalarTypeAsString() { return vtkImageScalarTypeNameMacro(this->GetScalarType()); }
648 
650 
654  static void SetNumberOfScalarComponents(int n, vtkInformation* meta_data);
659 
664  void CopyInformationFromPipeline(vtkInformation* information) override;
665 
671  void CopyInformationToPipeline(vtkInformation* information) override;
672 
678  void PrepareForNewData() override;
679 
681 
684  void ShallowCopy(vtkDataObject* src) override;
685  void DeepCopy(vtkDataObject* src) override;
687 
688  //--------------------------------------------------------------------------
689  // Methods that apply to any array (not just scalars).
690  // I am starting to experiment with generalizing imaging filters
691  // to operate on more than just scalars.
692 
694 
700  void* GetArrayPointer(vtkDataArray* array, int coordinates[3]);
702 
704 
711  vtkIdType GetTupleIndex(vtkDataArray* array, int coordinates[3]);
713 
718  void GetArrayIncrements(vtkDataArray* array, vtkIdType increments[3]);
719 
726  void ComputeInternalExtent(int* intExt, int* tgtExt, int* bnds);
727 
731  int GetExtentType() override { return VTK_3D_EXTENT; }
732 
734 
738  static vtkImageData* GetData(vtkInformationVector* v, int i = 0);
740 
741 protected:
743  ~vtkImageData() override;
744 
745  // The extent of what is currently in the structured grid.
746  // Dimensions is just an array to return a value.
747  // Its contents are out of data until GetDimensions is called.
748  int Dimensions[3];
749  vtkIdType Increments[3];
750 
751  // Variables used to define dataset physical orientation
752  double Origin[3];
753  double Spacing[3];
757 
758  int Extent[6];
759 
763 
764  // The first method assumes Active Scalars
765  void ComputeIncrements();
766  // This one is given the number of components of the
767  // scalar field explicitly
768  void ComputeIncrements(int numberOfComponents);
769  void ComputeIncrements(vtkDataArray* scalars);
770 
771  // The first method assumes Acitive Scalars
773  // This one is given the number of components of the
774  // scalar field explicitly
775  void ComputeIncrements(int numberOfComponents, vtkIdType inc[3]);
776  void ComputeIncrements(vtkDataArray* scalars, vtkIdType inc[3]);
777 
778  // for the index to physical methods
780 
782  void BuildPoints();
783  void BuildCells();
785 
786 private:
787  void InternalImageDataCopy(vtkImageData* src);
788 
789  friend class vtkUniformGrid;
790 
791  // for the GetPoint method
792  double Point[3];
793 
794  int DataDescription;
795  bool DirectionMatrixIsIdentity;
796 
797  vtkImageData(const vtkImageData&) = delete;
798  void operator=(const vtkImageData&) = delete;
799 };
800 
801 //----------------------------------------------------------------------------
803 {
804  this->ComputeIncrements(this->Increments);
805 }
806 
807 //----------------------------------------------------------------------------
808 inline void vtkImageData::ComputeIncrements(int numberOfComponents)
809 {
810  this->ComputeIncrements(numberOfComponents, this->Increments);
811 }
812 
813 //----------------------------------------------------------------------------
815 {
816  this->ComputeIncrements(scalars, this->Increments);
817 }
818 
819 //----------------------------------------------------------------------------
821 {
822  this->GetPoint(id, this->Point);
823  return this->Point;
824 }
825 
826 //----------------------------------------------------------------------------
828 {
830 }
831 
832 //----------------------------------------------------------------------------
834 {
836 }
837 
838 //----------------------------------------------------------------------------
840 {
841  return vtkStructuredData::GetDataDimension(this->DataDescription);
842 }
843 
844 //----------------------------------------------------------------------------
846 {
847  return vtkStructuredData::GetDataDimension(this->DataDescription);
848 }
849 VTK_ABI_NAMESPACE_END
850 #endif
abstract class to specify cell behavior
Definition: vtkCell.h:130
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:155
general representation of visualization data
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
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 spatial dimensionality of the data which is the maximum dimension of all cells.
virtual double * GetPoint(vtkIdType ptId)=0
Get point coordinates with ptId such that: 0 <= ptId < NumberOfPoints.
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:133
topologically and geometrically regular array of data
Definition: vtkImageData.h:156
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).
virtual int GetDataDimension()
Return the dimensionality of the data.
Definition: vtkImageData.h:839
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.
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.
vtkIdType GetCellSize(vtkIdType cellId) override
Standard vtkDataSet API methods.
vtkCell * FindAndGetCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
vtkSmartPointer< vtkPoints > StructuredPoints
Definition: vtkImageData.h:760
virtual void UnBlankCell(int i, int j, int k)
Methods for supporting blanking of cells.
virtual vtkIdType * GetIncrements()
Different ways to get the increments for moving around the data.
void GetCellPoints(vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts, vtkIdList *ptIds) 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.
Definition: vtkImageData.h:374
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.
vtkMatrix4x4 * IndexToPhysicalMatrix
Definition: vtkImageData.h:755
vtkConstantArray< int > * GetCellTypesArray()
Get the array of all cell types in the image data.
static vtkImageData * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
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()
Definition: vtkImageData.h:802
void BuildPoints()
vtkPoints * GetPoints() override
Standard vtkDataSet API methods.
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.
vtkCell * GetCell(vtkIdType cellId) override
Standard vtkDataSet API methods.
void ComputeBounds() override
Standard vtkDataSet API methods.
double * GetPoint(vtkIdType ptId) override
Standard vtkDataSet API methods.
Definition: vtkImageData.h:820
static vtkImageData * New()
void BuildCellTypes()
virtual double GetScalarTypeMin(vtkInformation *meta_data)
These returns the minimum and maximum values the ScalarType can hold without overflowing.
static vtkImageData * ExtendedNew()
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.
virtual int * GetDimensions()
Get dimensions of this structured points dataset.
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.
Definition: vtkImageData.h:845
virtual double GetScalarTypeMax(vtkInformation *meta_data)
These returns the minimum and maximum values the ScalarType can hold without overflowing.
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])
const char * GetScalarTypeAsString()
Definition: vtkImageData.h:647
int GetNumberOfScalarComponents()
Set/Get the number of scalar components for points.
vtkIdType GetNumberOfPoints() override
Standard vtkDataSet API methods.
Definition: vtkImageData.h:827
vtkSmartPointer< vtkStructuredCellArray > StructuredCells
Definition: vtkImageData.h:761
vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Standard vtkDataSet API methods.
vtkCell * GetCell(int i, int j, int k) 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.
int GetExtentType() override
The extent type is a 3D extent.
Definition: vtkImageData.h:731
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]
Definition: vtkImageData.h:749
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.
Definition: vtkImageData.h:833
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).
int GetDataObjectType() override
Return what type of dataset this is.
Definition: vtkImageData.h:178
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
Definition: vtkImageData.h:754
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.
static vtkImageData * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
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])
virtual void UnBlankPoint(int i, int j, int k)
Methods for supporting blanking of cells.
vtkStructuredCellArray * GetCells()
Return the image data connectivity array.
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.
Definition: vtkImageData.h:216
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.
Definition: vtkImageData.h:209
virtual void BlankCell(int i, int j, int k)
Methods for supporting blanking of cells.
vtkSmartPointer< vtkConstantArray< int > > StructuredCellTypes
Definition: vtkImageData.h:762
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
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)
void * GetArrayPointerForExtent(vtkDataArray *array, int extent[6])
These are convenience methods for getting a pointer from any filed array.
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()
vtkMatrix4x4 * PhysicalToIndexMatrix
Definition: vtkImageData.h:756
virtual vtkIdType ComputeCellId(int ijk[3])
Given a location in structured coordinates (i-j-k), return the cell id.
Definition: vtkImageData.h:382
void * GetArrayPointer(vtkDataArray *array, int coordinates[3])
These are convenience methods for getting a pointer from any filed array.
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
Definition: vtkMatrix3x3.h:56
represent and manipulate 4x4 transformation matrices
Definition: vtkMatrix4x4.h:141
a cell that represents an orthogonal quadrilateral
Definition: vtkPixel.h:66
represent and manipulate 3D points
Definition: vtkPoints.h:139
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
@ component
Definition: vtkX3D.h:175
@ info
Definition: vtkX3D.h:376
@ direction
Definition: vtkX3D.h:260
@ extent
Definition: vtkX3D.h:345
@ spacing
Definition: vtkX3D.h:481
#define VTK_3D_EXTENT
int vtkIdType
Definition: vtkType.h:315
#define VTK_IMAGE_DATA
Definition: vtkType.h:71
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO
#define max(a, b)