ITK  6.0.0
Insight Toolkit
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Attributes | Friends | List of all members
itk::Image< TPixel, VImageDimension > Class Template Reference

#include <itkImage.h>

Detailed Description

template<typename TPixel, unsigned int VImageDimension = 2>
class itk::Image< TPixel, VImageDimension >

Templated n-dimensional image class.

Images are templated over a pixel type (modeling the dependent variables), and a dimension (number of independent variables). The container for the pixel data is the ImportImageContainer.

Within the pixel container, images are modeled as arrays, defined by a start index and a size.

The superclass of Image, ImageBase, defines the geometry of the image in terms of where the image sits in physical space, how the image is oriented in physical space, the size of a pixel, and the extent of the image itself. ImageBase provides the methods to convert between the index and physical space coordinate frames.

Pixels can be accessed directly using the SetPixel() and GetPixel() methods or can be accessed via iterators that define the region of the image they traverse.

The pixel type may be one of the native types; a Insight-defined class type such as Vector; or a user-defined type. Note that depending on the type of pixel that you use, the process objects (i.e., those filters processing data objects) may not operate on the image and/or pixel type. This becomes apparent at compile-time because operator overloading (for the pixel type) is not supported.

The data in an image is arranged in a 1D array as [][][][slice][row][col] with the column index varying most rapidly. The Index type reverses the order so that with Index[0] = col, Index[1] = row, Index[2] = slice, ...

See also
ImageBase
ImageContainerInterface
ITK Sphinx Examples:
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=
  • <a href=

Definition at line 88 of file itkImage.h.

+ Inheritance diagram for itk::Image< TPixel, VImageDimension >:
+ Collaboration diagram for itk::Image< TPixel, VImageDimension >:

Classes

struct  Rebind
 

Public Types

using AccessorFunctorType = DefaultPixelAccessorFunctor< Self >
 
using AccessorType = DefaultPixelAccessor< PixelType >
 
using ConstPointer = SmartPointer< const Self >
 
using ConstWeakPointer = WeakPointer< const Self >
 
using InternalPixelType = TPixel
 
using IOPixelType = PixelType
 
using NeighborhoodAccessorFunctorType = NeighborhoodAccessorFunctor< Self >
 
using PixelContainer = ImportImageContainer< SizeValueType, PixelType >
 
using PixelContainerConstPointer = typename PixelContainer::ConstPointer
 
using PixelContainerPointer = typename PixelContainer::Pointer
 
using PixelType = TPixel
 
using Pointer = SmartPointer< Self >
 
template<typename UPixelType , unsigned int VUImageDimension = VImageDimension>
using RebindImageType = itk::Image< UPixelType, VUImageDimension >
 
using Self = Image
 
using Superclass = ImageBase< VImageDimension >
 
using ValueType = TPixel
 
- Public Types inherited from itk::ImageBase< 2 >
using ConstPointer = SmartPointer< const Self >
 
using DirectionType = Matrix< SpacePrecisionType, VImageDimension, VImageDimension >
 
using ImageDimensionType = unsigned int
 
using IndexType = Index< VImageDimension >
 
using IndexValueType = typename IndexType::IndexValueType
 
using OffsetType = Offset< VImageDimension >
 
using OffsetValueType = typename OffsetType::OffsetValueType
 
using Pointer = SmartPointer< Self >
 
using PointType = Point< PointValueType, VImageDimension >
 
using PointValueType = SpacePrecisionType
 
using RegionType = ImageRegion< VImageDimension >
 
using Self = ImageBase
 
using SizeType = Size< VImageDimension >
 
using SizeValueType = typename SizeType::SizeValueType
 
using SpacingType = Vector< SpacingValueType, VImageDimension >
 
using SpacingValueType = SpacePrecisionType
 
using Superclass = DataObject
 
- Public Types inherited from itk::DataObject
using ConstPointer = SmartPointer< const Self >
 
using DataObjectIdentifierType = std::string
 
using DataObjectPointerArraySizeType = std::vector< Pointer >::size_type
 
using Pointer = SmartPointer< Self >
 
using Self = DataObject
 
using Superclass = Object
 
- Public Types inherited from itk::Object
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = Object
 
using Superclass = LightObject
 
- Public Types inherited from itk::LightObject
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = LightObject
 

Public Member Functions

void Allocate (bool initializePixels=false) override
 
void FillBuffer (const TPixel &value)
 
const char * GetNameOfClass () const override
 
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor ()
 
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor () const
 
unsigned int GetNumberOfComponentsPerPixel () const override
 
TPixel & GetPixel (const IndexType &index)
 
const TPixel & GetPixel (const IndexType &index) const
 
AccessorType GetPixelAccessor ()
 
const AccessorType GetPixelAccessor () const
 
PixelContainerGetPixelContainer ()
 
const PixelContainerGetPixelContainer () const
 
virtual void Graft (const Self *image)
 
void Initialize () override
 
TPixel & operator[] (const IndexType &index)
 
const TPixel & operator[] (const IndexType &index) const
 
void SetPixel (const IndexType &index, const TPixel &value)
 
void SetPixelContainer (PixelContainer *container)
 
virtual TPixel * GetBufferPointer ()
 
virtual const TPixel * GetBufferPointer () const
 
- Public Member Functions inherited from itk::ImageBase< 2 >
virtual void Allocate (bool initialize=false)
 
void AllocateInitialized ()
 
OffsetValueType ComputeOffset (const IndexType &ind) const
 
void CopyInformation (const DataObject *data) override
 
virtual const RegionTypeGetBufferedRegion () const
 
virtual const DirectionTypeGetDirection () const
 
virtual const DirectionTypeGetInverseDirection () const
 
virtual const RegionTypeGetLargestPossibleRegion () const
 
const char * GetNameOfClass () const override
 
virtual const PointTypeGetOrigin () const
 
virtual const RegionTypeGetRequestedRegion () const
 
virtual const SpacingTypeGetSpacing () const
 
virtual void Graft (const Self *image)
 
void Initialize () override
 
bool IsCongruentImageGeometry (const ImageBase *otherImage, double coordinateTolerance, double directionTolerance) const
 
bool IsSameImageGeometryAs (const ImageBase *otherImage, double coordinateTolerance=DefaultImageCoordinateTolerance, double directionTolerance=DefaultImageDirectionTolerance) const
 
bool RequestedRegionIsOutsideOfTheBufferedRegion () override
 
virtual void SetBufferedRegion (const RegionType &region)
 
virtual void SetDirection (const DirectionType &direction)
 
virtual void SetLargestPossibleRegion (const RegionType &region)
 
virtual void SetRegions (const SizeType &size)
 
void SetRequestedRegion (const DataObject *data) override
 
virtual void SetRequestedRegion (const RegionType &region)
 
void SetRequestedRegionToLargestPossibleRegion () override
 
void TransformLocalVectorToPhysicalVector (const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
 
ContinuousIndex< TIndexRep, VImageDimension > TransformPhysicalPointToContinuousIndex (const Point< TCoordRep, VImageDimension > &point) const
 
bool TransformPhysicalPointToContinuousIndex (const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const
 
bool TransformPhysicalPointToIndex (const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
 
void TransformPhysicalVectorToLocalVector (const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
 
void UpdateOutputData () override
 
void UpdateOutputInformation () override
 
bool VerifyRequestedRegion () override
 
virtual void SetOrigin (PointType _arg)
 
virtual void SetOrigin (const double origin[VImageDimension])
 
virtual void SetOrigin (const float origin[VImageDimension])
 
virtual void SetRegions (const RegionType &region)
 
const OffsetValueTypeGetOffsetTable () const
 
IndexType ComputeIndex (OffsetValueType offset) const
 
virtual void SetSpacing (const SpacingType &spacing)
 
virtual void SetSpacing (const double spacing[VImageDimension])
 
virtual void SetSpacing (const float spacing[VImageDimension])
 
IndexType TransformPhysicalPointToIndex (const Point< TCoordRep, VImageDimension > &point) const
 
void TransformContinuousIndexToPhysicalPoint (const ContinuousIndex< TIndexRep, VImageDimension > &index, Point< TCoordRep, VImageDimension > &point) const
 
Point< TCoordRep, VImageDimension > TransformContinuousIndexToPhysicalPoint (const ContinuousIndex< TIndexRep, VImageDimension > &index) const
 
void TransformIndexToPhysicalPoint (const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
 
Point< TCoordRep, VImageDimension > TransformIndexToPhysicalPoint (const IndexType &index) const
 
TVector TransformLocalVectorToPhysicalVector (const TVector &inputGradient) const
 
TVector TransformPhysicalVectorToLocalVector (const TVector &inputGradient) const
 
virtual void SetNumberOfComponentsPerPixel (unsigned int)
 
- Public Member Functions inherited from itk::DataObject
virtual void CopyInformation (const DataObject *)
 
virtual void DataHasBeenGenerated ()
 
void DisconnectPipeline ()
 
bool GetDataReleased () const
 
const char * GetNameOfClass () const override
 
virtual const bool & GetReleaseDataFlag () const
 
SmartPointer< ProcessObjectGetSource () const
 
DataObjectPointerArraySizeType GetSourceOutputIndex () const
 
const DataObjectIdentifierTypeGetSourceOutputName () const
 
virtual ModifiedTimeType GetUpdateMTime () const
 
virtual void Graft (const DataObject *)
 
virtual void Initialize ()
 
virtual void PrepareForNewData ()
 
virtual void PropagateRequestedRegion ()
 
void ReleaseData ()
 
virtual void ReleaseDataFlagOn ()
 
virtual bool RequestedRegionIsOutsideOfTheBufferedRegion ()
 
virtual void ResetPipeline ()
 
void SetReleaseDataFlag (bool flag)
 
virtual void SetRequestedRegion (const DataObject *)
 
virtual void SetRequestedRegionToLargestPossibleRegion ()
 
bool ShouldIReleaseData () const
 
virtual void Update ()
 
virtual void UpdateOutputData ()
 
virtual void UpdateOutputInformation ()
 
void UpdateSource () const
 
virtual bool VerifyRequestedRegion ()
 
void SetPipelineMTime (ModifiedTimeType time)
 
virtual const ModifiedTimeTypeGetPipelineMTime () const
 
virtual void SetRealTimeStamp (RealTimeStamp _arg)
 
virtual const RealTimeStampGetRealTimeStamp () const
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *cmd) const
 
unsigned long AddObserver (const EventObject &event, std::function< void(const EventObject &)> function) const
 
LightObject::Pointer CreateAnother () const override
 
virtual void DebugOff () const
 
virtual void DebugOn () const
 
CommandGetCommand (unsigned long tag)
 
bool GetDebug () const
 
MetaDataDictionaryGetMetaDataDictionary ()
 
const MetaDataDictionaryGetMetaDataDictionary () const
 
virtual ModifiedTimeType GetMTime () const
 
const char * GetNameOfClass () const override
 
virtual const TimeStampGetTimeStamp () const
 
bool HasObserver (const EventObject &event) const
 
void InvokeEvent (const EventObject &)
 
void InvokeEvent (const EventObject &) const
 
virtual void Modified () const
 
void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag) const
 
void SetDebug (bool debugFlag) const
 
void SetReferenceCount (int) override
 
void UnRegister () const noexcept override
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
void SetMetaDataDictionary (MetaDataDictionary &&rrhs)
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
Pointer Clone () const
 
virtual Pointer CreateAnother () const
 
virtual void Delete ()
 
virtual const char * GetNameOfClass () const
 
virtual int GetReferenceCount () const
 
void Print (std::ostream &os, Indent indent=0) const
 
virtual void Register () const
 
virtual void SetReferenceCount (int)
 
virtual void UnRegister () const noexcept
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::ImageBase< 2 >
static constexpr unsigned int GetImageDimension ()
 
static Pointer New ()
 
- Static Public Member Functions inherited from itk::DataObject
static bool GetGlobalReleaseDataFlag ()
 
static void GlobalReleaseDataFlagOff ()
 
static void GlobalReleaseDataFlagOn ()
 
static Pointer New ()
 
static void SetGlobalReleaseDataFlag (bool val)
 
- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool val)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 

Protected Member Functions

void ComputeIndexToPhysicalPointMatrices () override
 
void Graft (const DataObject *data) override
 
 Image ()=default
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~Image () override=default
 
- Protected Member Functions inherited from itk::ImageBase< 2 >
virtual void ComputeIndexToPhysicalPointMatrices ()
 
void ComputeOffsetTable ()
 
void Graft (const DataObject *data) override
 
 ImageBase ()=default
 
virtual void InitializeBufferedRegion ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~ImageBase () override=default
 
OffsetValueType FastComputeOffset (const IndexType &ind) const
 
IndexType FastComputeIndex (OffsetValueType offset) const
 
- Protected Member Functions inherited from itk::DataObject
 DataObject ()
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
virtual void PropagateResetPipeline ()
 
 ~DataObject () override
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
virtual void SetTimeStamp (const TimeStamp &timeStamp)
 
 ~Object () override
 
- Protected Member Functions inherited from itk::LightObject
virtual LightObject::Pointer InternalClone () const
 
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintSelf (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 

Private Attributes

PixelContainerPointer m_Buffer { PixelContainer::New() }
 

Friends

template<typename TEqualityComparable >
std::enable_if_t< std::is_same_v< TEqualityComparable, TPixel >, bool > operator!= (const Image< TEqualityComparable, VImageDimension > &lhs, const Image< TEqualityComparable, VImageDimension > &rhs)
 
template<typename TEqualityComparable >
std::enable_if_t< std::is_same_v< TEqualityComparable, TPixel >, bool > operator== (const Image< TEqualityComparable, VImageDimension > &lhs, const Image< TEqualityComparable, VImageDimension > &rhs)
 

Additional Inherited Members

- Static Public Attributes inherited from itk::ImageBase< 2 >
static constexpr ImageDimensionType ImageDimension
 
- Protected Attributes inherited from itk::ImageBase< 2 >
SpacingType m_Spacing
 
PointType m_Origin
 
DirectionType m_Direction
 
DirectionType m_InverseDirection
 
DirectionType m_IndexToPhysicalPoint
 
DirectionType m_PhysicalPointToIndex
 
- Protected Attributes inherited from itk::LightObject
std::atomic< int > m_ReferenceCount {}
 

Member Typedef Documentation

◆ AccessorFunctorType

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::AccessorFunctorType = DefaultPixelAccessorFunctor<Self>

Definition at line 124 of file itkImage.h.

◆ AccessorType

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::AccessorType = DefaultPixelAccessor<PixelType>

Accessor type that convert data between internal and external representations.

Definition at line 123 of file itkImage.h.

◆ ConstPointer

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::ConstPointer = SmartPointer<const Self>

Definition at line 97 of file itkImage.h.

◆ ConstWeakPointer

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::ConstWeakPointer = WeakPointer<const Self>

Definition at line 98 of file itkImage.h.

◆ InternalPixelType

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::InternalPixelType = TPixel

Internal Pixel representation. Used to maintain a uniform API with Image Adaptors and allow to keep a particular internal representation of data while showing a different external representation.

Definition at line 117 of file itkImage.h.

◆ IOPixelType

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::IOPixelType = PixelType

Definition at line 119 of file itkImage.h.

◆ NeighborhoodAccessorFunctorType

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::NeighborhoodAccessorFunctorType = NeighborhoodAccessorFunctor<Self>

Typedef for the functor used to access a neighborhood of pixel pointers.

Definition at line 128 of file itkImage.h.

◆ PixelContainer

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::PixelContainer = ImportImageContainer<SizeValueType, PixelType>

Container used to store pixels in the image.

Definition at line 145 of file itkImage.h.

◆ PixelContainerConstPointer

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::PixelContainerConstPointer = typename PixelContainer::ConstPointer

Definition at line 165 of file itkImage.h.

◆ PixelContainerPointer

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::PixelContainerPointer = typename PixelContainer::Pointer

A pointer to the pixel container.

Definition at line 164 of file itkImage.h.

◆ PixelType

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::PixelType = TPixel

Pixel type alias support. Used to declare pixel type in filters or other operations.

Definition at line 108 of file itkImage.h.

◆ Pointer

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::Pointer = SmartPointer<Self>

Definition at line 96 of file itkImage.h.

◆ RebindImageType

template<typename TPixel , unsigned int VImageDimension = 2>
template<typename UPixelType , unsigned int VUImageDimension = VImageDimension>
using itk::Image< TPixel, VImageDimension >::RebindImageType = itk::Image<UPixelType, VUImageDimension>

Definition at line 184 of file itkImage.h.

◆ Self

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::Self = Image

Standard class type aliases

Definition at line 94 of file itkImage.h.

◆ Superclass

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::Superclass = ImageBase<VImageDimension>

Definition at line 95 of file itkImage.h.

◆ ValueType

template<typename TPixel , unsigned int VImageDimension = 2>
using itk::Image< TPixel, VImageDimension >::ValueType = TPixel

Typedef alias for PixelType

Definition at line 111 of file itkImage.h.

Constructor & Destructor Documentation

◆ Image()

template<typename TPixel , unsigned int VImageDimension = 2>
itk::Image< TPixel, VImageDimension >::Image ( )
protecteddefault

◆ ~Image()

template<typename TPixel , unsigned int VImageDimension = 2>
itk::Image< TPixel, VImageDimension >::~Image ( )
overrideprotecteddefault

Member Function Documentation

◆ Allocate()

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::Allocate ( bool  initializePixels = false)
overridevirtual

Allocate the image memory. The size of the image must already be set, e.g. by calling SetRegions().

Reimplemented from itk::ImageBase< 2 >.

◆ ComputeIndexToPhysicalPointMatrices()

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::ComputeIndexToPhysicalPointMatrices ( )
overrideprotectedvirtual

Compute helper matrices used to transform Index coordinates to PhysicalPoint coordinates and back. This method is virtual and will be overloaded in derived classes in order to provide backward compatibility behavior in classes that did not used to take image orientation into account.

Reimplemented from itk::ImageBase< 2 >.

◆ FillBuffer()

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::FillBuffer ( const TPixel &  value)

Fill the image buffer with a value. Be sure to call Allocate() first.

◆ GetBufferPointer() [1/2]

template<typename TPixel , unsigned int VImageDimension = 2>
virtual TPixel * itk::Image< TPixel, VImageDimension >::GetBufferPointer ( )
inlinevirtual

Return a pointer to the beginning of the buffer. This is used by the image iterator class.

Reimplemented in itk::GPUImage< TPixel, VImageDimension >.

Definition at line 251 of file itkImage.h.

◆ GetBufferPointer() [2/2]

template<typename TPixel , unsigned int VImageDimension = 2>
virtual const TPixel * itk::Image< TPixel, VImageDimension >::GetBufferPointer ( ) const
inlinevirtual

Return a pointer to the beginning of the buffer. This is used by the image iterator class.

Reimplemented in itk::GPUImage< TPixel, VImageDimension >.

Definition at line 256 of file itkImage.h.

◆ GetNameOfClass()

template<typename TPixel , unsigned int VImageDimension = 2>
const char * itk::Image< TPixel, VImageDimension >::GetNameOfClass ( ) const
overridevirtual
See also
LightObject::GetNameOfClass()

Reimplemented from itk::DataObject.

◆ GetNeighborhoodAccessor() [1/2]

template<typename TPixel , unsigned int VImageDimension = 2>
NeighborhoodAccessorFunctorType itk::Image< TPixel, VImageDimension >::GetNeighborhoodAccessor ( )
inline

Return the NeighborhoodAccessor functor

Definition at line 309 of file itkImage.h.

◆ GetNeighborhoodAccessor() [2/2]

template<typename TPixel , unsigned int VImageDimension = 2>
const NeighborhoodAccessorFunctorType itk::Image< TPixel, VImageDimension >::GetNeighborhoodAccessor ( ) const
inline

Return the NeighborhoodAccessor functor

Definition at line 316 of file itkImage.h.

◆ GetNumberOfComponentsPerPixel()

template<typename TPixel , unsigned int VImageDimension = 2>
unsigned int itk::Image< TPixel, VImageDimension >::GetNumberOfComponentsPerPixel ( ) const
overridevirtual

INTERNAL This method is used internally by filters to copy meta-data from the output to the input. Users should not have a need to use this method.

Filters that override the ProcessObject's GenerateOutputInformation() should generally have the following line if they want to propagate meta- data for both Image and VectorImage

outputImage->SetNumberOfComponentsPerPixel(
inputImage->GetNumberOfComponentsPerPixel() )
See also
ImageBase, VectorImage

Returns/Sets the number of components in the image. Note that in the ImageBase implementation, this always returns 1. Image returns the

returned from NumericTraits for the pixel type, and VectorImage

returns the vector length set by the user.

Reimplemented from itk::ImageBase< 2 >.

◆ GetPixel() [1/2]

template<typename TPixel , unsigned int VImageDimension = 2>
TPixel & itk::Image< TPixel, VImageDimension >::GetPixel ( const IndexType index)
inline

Get a reference to a pixel (e.g. for editing).

For efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 230 of file itkImage.h.

◆ GetPixel() [2/2]

template<typename TPixel , unsigned int VImageDimension = 2>
const TPixel & itk::Image< TPixel, VImageDimension >::GetPixel ( const IndexType index) const
inline

Get a pixel (read only version).

For efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 219 of file itkImage.h.

◆ GetPixelAccessor() [1/2]

template<typename TPixel , unsigned int VImageDimension = 2>
AccessorType itk::Image< TPixel, VImageDimension >::GetPixelAccessor ( )
inline

Return the Pixel Accessor object

Definition at line 295 of file itkImage.h.

◆ GetPixelAccessor() [2/2]

template<typename TPixel , unsigned int VImageDimension = 2>
const AccessorType itk::Image< TPixel, VImageDimension >::GetPixelAccessor ( ) const
inline

Return the Pixel Accesor object

Definition at line 302 of file itkImage.h.

◆ GetPixelContainer() [1/2]

template<typename TPixel , unsigned int VImageDimension = 2>
PixelContainer * itk::Image< TPixel, VImageDimension >::GetPixelContainer ( )
inline

Return a pointer to the container.

Definition at line 264 of file itkImage.h.

◆ GetPixelContainer() [2/2]

template<typename TPixel , unsigned int VImageDimension = 2>
const PixelContainer * itk::Image< TPixel, VImageDimension >::GetPixelContainer ( ) const
inline

Definition at line 270 of file itkImage.h.

◆ Graft() [1/2]

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::Graft ( const DataObject image)
overrideprotectedvirtual

Graft the data and information from one image to another. This is a convenience method to setup a second image with all the meta information of another image and use the same pixel container. Note that this method is different than just using two SmartPointers to the same image since separate DataObjects are still maintained. This method is similar to ImageSource::GraftOutput(). The implementation in ImageBase simply calls CopyInformation() and copies the region ivars. Subclasses of ImageBase are responsible for copying the pixel container.

Reimplemented from itk::ImageBase< 2 >.

◆ Graft() [2/2]

template<typename TPixel , unsigned int VImageDimension = 2>
virtual void itk::Image< TPixel, VImageDimension >::Graft ( const Self image)
virtual

Graft the data and information from one image to another. This is a convenience method to setup a second image with all the meta information of another image and use the same pixel container. Note that this method is different than just using two SmartPointers to the same image since separate DataObjects are still maintained. This method is similar to ImageSource::GraftOutput(). The implementation in ImageBase simply calls CopyInformation() and copies the region ivars. The implementation here refers to the superclass' implementation and then copies over the pixel container.

Reimplemented from itk::ImageBase< 2 >.

Reimplemented in itk::GPUImage< TPixel, VImageDimension >, and itk::GPUImage< TPixel, VImageDimension >.

◆ Initialize()

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::Initialize ( )
overridevirtual

Restore the data object to its initial state. This means releasing memory.

Reimplemented from itk::DataObject.

◆ New()

template<typename TPixel , unsigned int VImageDimension = 2>
static Pointer itk::Image< TPixel, VImageDimension >::New ( )
static

◆ operator[]() [1/2]

template<typename TPixel , unsigned int VImageDimension = 2>
TPixel & itk::Image< TPixel, VImageDimension >::operator[] ( const IndexType index)
inline

Access a pixel. This version can be an lvalue.

For efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 240 of file itkImage.h.

◆ operator[]() [2/2]

template<typename TPixel , unsigned int VImageDimension = 2>
const TPixel & itk::Image< TPixel, VImageDimension >::operator[] ( const IndexType index) const
inline

Access a pixel. This version can only be an rvalue.

For efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 246 of file itkImage.h.

◆ PrintSelf()

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual

Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from itk::DataObject.

◆ SetPixel()

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::SetPixel ( const IndexType index,
const TPixel &  value 
)
inline

Set a pixel value.

Allocate() needs to have been called first – for efficiency, this function does not check that the image has actually been allocated yet.

Definition at line 208 of file itkImage.h.

◆ SetPixelContainer()

template<typename TPixel , unsigned int VImageDimension = 2>
void itk::Image< TPixel, VImageDimension >::SetPixelContainer ( PixelContainer container)

Set the container to use. Note that this does not cause the DataObject to be modified.

Friends And Related Function Documentation

◆ operator!=

template<typename TPixel , unsigned int VImageDimension = 2>
template<typename TEqualityComparable >
std::enable_if_t< std::is_same_v< TEqualityComparable, TPixel >, bool > operator!= ( const Image< TEqualityComparable, VImageDimension > &  lhs,
const Image< TEqualityComparable, VImageDimension > &  rhs 
)
friend

Returns (image1 != image2).

Definition at line 371 of file itkImage.h.

◆ operator==

template<typename TPixel , unsigned int VImageDimension = 2>
template<typename TEqualityComparable >
std::enable_if_t< std::is_same_v< TEqualityComparable, TPixel >, bool > operator== ( const Image< TEqualityComparable, VImageDimension > &  lhs,
const Image< TEqualityComparable, VImageDimension > &  rhs 
)
friend

Returns (image1 == image2).

Note
operator== and operator!= are defined as function templates (rather than as non-templates), just to allow template instantiation of itk::Image for non-EqualityComparable pixel types.

Definition at line 331 of file itkImage.h.

Member Data Documentation

◆ m_Buffer

template<typename TPixel , unsigned int VImageDimension = 2>
PixelContainerPointer itk::Image< TPixel, VImageDimension >::m_Buffer { PixelContainer::New() }
private

Memory for the current buffer.

Definition at line 397 of file itkImage.h.


The documentation for this class was generated from the following file: