SimpleITK  
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
itk::simple::Image Class Reference

The Image class for SimpleITK. More...

#include <sitkImage.h>

+ Collaboration diagram for itk::simple::Image:

Detailed Description

The Image class for SimpleITK.

This Image class can represent 2D, 3D, and 4D images. The pixel types may be a scalar, a multi-component vector or a run-length-encoded (RLE) "label". The dimension, pixel type and size is specified at construction.

A fundamental concept of ITK images is that they occupy physical space where the image is defined by an origin, spacing, and direction cosine matrix. The attributes are taken into consideration when doing most operations on an image. A meta-data dictionary is also associated with the image, which may contain additional fields from reading but these attributes are not propagated by image filters.

The SimpleITK Image provides a single facade interface to several ITK image types. Internally, the SimpleITK Image maintains a pointer to the ITK image class, and performs reference counting and lazy copying. This means that deep copying of an image including it's buffer is delayed until the image is modified. This removes the need to use pointers to SimpleITK Image class, as copying and returning by value do not unnecessarily duplicate the data.

See also
itk::Image itk::VectorImage itk::LabelMap itk::ImageBase
Examples
BufferImportExport.cxx, CppCMake/Source/sitk_example.cxx, CppInPlace/CppInPlace.cxx, DemonsRegistration1/DemonsRegistration1.cxx, DemonsRegistration2/DemonsRegistration2.cxx, DicomSeriesReader/DicomSeriesReader.cxx, FastMarchingSegmentation/FastMarchingSegmentation.cxx, FilterProgressReporting/FilterProgressReporting.cxx, HelloWorld/HelloWorld.cxx, ImageIOSelection/ImageIOSelection.cxx, ImageRegistrationMethod1/ImageRegistrationMethod1.cxx, ImageRegistrationMethod2/ImageRegistrationMethod2.cxx, ImageRegistrationMethodBSpline1/ImageRegistrationMethodBSpline1.cxx, ImageRegistrationMethodBSpline3/ImageRegistrationMethodBSpline3.cxx, ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx, ITKIntegration/ITKIntegration.cxx, N4BiasFieldCorrection/N4BiasFieldCorrection.cxx, Segmentation/ConnectedThresholdImageFilter.cxx, Segmentation/NeighborhoodConnectedImageFilter.cxx, SimpleGaussian/SimpleGaussian.cxx, SimpleGaussianFunctional.cxx, SimpleIO/SimpleIO.cxx, and SliceBySliceDecorator/SliceBySliceDecorator.cxx.

Definition at line 76 of file sitkImage.h.

Public Types

using Self = Image
 

Public Member Functions

void CopyInformation (const Image &srcImage)
 Copy common meta-data from an image to this one. More...
 
bool EraseMetaData (const std::string &key)
 Remove an entry from the meta-data dictionary. More...
 
std::vector< double > EvaluateAtContinuousIndex (const std::vector< double > &index, InterpolatorEnum interp=sitkLinear) const
 Interpolate pixel value at a continuous index. More...
 
std::vector< double > EvaluateAtPhysicalPoint (const std::vector< double > &point, InterpolatorEnum interp=sitkLinear) const
 
unsigned int GetDepth () const
 
unsigned int GetDimension () const
 
unsigned int GetHeight () const
 
std::string GetMetaData (const std::string &key) const
 Get the value of a meta-data dictionary entry as a string. More...
 
std::vector< std::string > GetMetaDataKeys () const
 get a vector of keys in from the meta-data dictionary More...
 
unsigned int GetNumberOfComponentsPerPixel () const
 Get the number of components for each pixel. More...
 
uint64_t GetNumberOfPixels () const
 Get the number of pixels in the image. More...
 
PixelIDValueEnum GetPixelID () const
 
std::string GetPixelIDTypeAsString () const
 
PixelIDValueType GetPixelIDValue () const
 
std::vector< unsigned int > GetSize () const
 
unsigned int GetSizeOfPixelComponent () const
 Get the number of bytes per component of a pixel. More...
 
unsigned int GetWidth () const
 
bool HasMetaDataKey (const std::string &key) const
 Query the meta-data dictionary for the existence of a key. More...
 
 Image ()
 Default constructor, creates an image of size 0. More...
 
 Image (const Image &img)
 
 Image (Image &&img) noexcept
 Move constructor and assignment. More...
 
bool IsCongruentImageGeometry (const Image &otherImage, double coordinateTolerance, double directionTolerance) const
 
bool IsSameImageGeometryAs (const Image &otherImage, double=DefaultImageCoordinateTolerance, double=DefaultImageDirectionTolerance) const
 
bool IsUnique () const
 Returns true if no other SimpleITK Image object refers to the same internal data structure. More...
 
void MakeUnique ()
 Performs actually coping if needed to make object unique. More...
 
Imageoperator= (const Image &img)
 
Imageoperator= (Image &&img) noexcept
 
Image ProxyForInPlaceOperation ()
 Advanced method not commonly needed. More...
 
void SetMetaData (const std::string &key, const std::string &value)
 Set an entry in the meta-data dictionary. More...
 
Image ToScalarImage (bool inPlace=true)
 Convert a image of vector pixel type to a scalar image with N+1 dimensions. More...
 
std::string ToString () const
 
Image ToVectorImage (bool inPlace=true)
 Convert the first dimension to the components for image with vector pixel type. More...
 
std::vector< double > TransformContinuousIndexToPhysicalPoint (const std::vector< double > &index) const
 
std::vector< double > TransformIndexToPhysicalPoint (const std::vector< int64_t > &index) const
 
std::vector< double > TransformPhysicalPointToContinuousIndex (const std::vector< double > &point) const
 
std::vector< int64_t > TransformPhysicalPointToIndex (const std::vector< double > &point) const
 
virtual ~Image ()
 
 Image (unsigned int width, unsigned int height, PixelIDValueEnum valueEnum)
 Constructors for 2D, 3D an optionally 4D images where pixel type and number of components can be specified. More...
 
 Image (unsigned int width, unsigned int height, unsigned int depth, PixelIDValueEnum valueEnum)
 Constructors for 2D, 3D an optionally 4D images where pixel type and number of components can be specified. More...
 
 Image (const std::vector< unsigned int > &size, PixelIDValueEnum valueEnum, unsigned int numberOfComponents=0)
 Constructors for 2D, 3D an optionally 4D images where pixel type and number of components can be specified. More...
 
template<typename TImageType >
 Image (itk::SmartPointer< TImageType > image)
 Construct an SimpleITK Image from an pointer to an ITK image. More...
 
template<typename TImageType >
 Image (TImageType *image)
 Construct an SimpleITK Image from an pointer to an ITK image. More...
 
itk::DataObjectGetITKBase ()
 
const itk::DataObjectGetITKBase () const
 
std::vector< double > GetOrigin () const
 
void SetOrigin (const std::vector< double > &origin)
 
std::vector< double > GetSpacing () const
 
void SetSpacing (const std::vector< double > &spacing)
 
std::vector< double > GetDirection () const
 Set/Get the Direction. More...
 
void SetDirection (const std::vector< double > &direction)
 Set/Get the Direction. More...
 
int8_t GetPixelAsInt8 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
uint8_t GetPixelAsUInt8 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
int16_t GetPixelAsInt16 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
uint16_t GetPixelAsUInt16 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
int32_t GetPixelAsInt32 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
uint32_t GetPixelAsUInt32 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
int64_t GetPixelAsInt64 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
uint64_t GetPixelAsUInt64 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
float GetPixelAsFloat (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
double GetPixelAsDouble (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< int8_t > GetPixelAsVectorInt8 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< uint8_t > GetPixelAsVectorUInt8 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< int16_t > GetPixelAsVectorInt16 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< uint16_t > GetPixelAsVectorUInt16 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< int32_t > GetPixelAsVectorInt32 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< uint32_t > GetPixelAsVectorUInt32 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< int64_t > GetPixelAsVectorInt64 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< uint64_t > GetPixelAsVectorUInt64 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< float > GetPixelAsVectorFloat32 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::vector< double > GetPixelAsVectorFloat64 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::complex< float > GetPixelAsComplexFloat32 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
std::complex< double > GetPixelAsComplexFloat64 (const std::vector< uint32_t > &idx) const
 Get the value of a pixel. More...
 
void SetPixelAsInt8 (const std::vector< uint32_t > &idx, int8_t v)
 Set the value of a pixel. More...
 
void SetPixelAsUInt8 (const std::vector< uint32_t > &idx, uint8_t v)
 Set the value of a pixel. More...
 
void SetPixelAsInt16 (const std::vector< uint32_t > &idx, int16_t v)
 Set the value of a pixel. More...
 
void SetPixelAsUInt16 (const std::vector< uint32_t > &idx, uint16_t v)
 Set the value of a pixel. More...
 
void SetPixelAsInt32 (const std::vector< uint32_t > &idx, int32_t v)
 Set the value of a pixel. More...
 
void SetPixelAsUInt32 (const std::vector< uint32_t > &idx, uint32_t v)
 Set the value of a pixel. More...
 
void SetPixelAsInt64 (const std::vector< uint32_t > &idx, int64_t v)
 Set the value of a pixel. More...
 
void SetPixelAsUInt64 (const std::vector< uint32_t > &idx, uint64_t v)
 Set the value of a pixel. More...
 
void SetPixelAsFloat (const std::vector< uint32_t > &idx, float v)
 Set the value of a pixel. More...
 
void SetPixelAsDouble (const std::vector< uint32_t > &idx, double v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorInt8 (const std::vector< uint32_t > &idx, const std::vector< int8_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorUInt8 (const std::vector< uint32_t > &idx, const std::vector< uint8_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorInt16 (const std::vector< uint32_t > &idx, const std::vector< int16_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorUInt16 (const std::vector< uint32_t > &idx, const std::vector< uint16_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorInt32 (const std::vector< uint32_t > &idx, const std::vector< int32_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorUInt32 (const std::vector< uint32_t > &idx, const std::vector< uint32_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorInt64 (const std::vector< uint32_t > &idx, const std::vector< int64_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorUInt64 (const std::vector< uint32_t > &idx, const std::vector< uint64_t > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorFloat32 (const std::vector< uint32_t > &idx, const std::vector< float > &v)
 Set the value of a pixel. More...
 
void SetPixelAsVectorFloat64 (const std::vector< uint32_t > &idx, const std::vector< double > &v)
 Set the value of a pixel. More...
 
void SetPixelAsComplexFloat32 (const std::vector< uint32_t > &idx, const std::complex< float > v)
 Set the value of a pixel. More...
 
void SetPixelAsComplexFloat64 (const std::vector< uint32_t > &idx, const std::complex< double > v)
 Set the value of a pixel. More...
 
int8_t * GetBufferAsInt8 ()
 Get a pointer to the image buffer. More...
 
uint8_t * GetBufferAsUInt8 ()
 Get a pointer to the image buffer. More...
 
int16_t * GetBufferAsInt16 ()
 Get a pointer to the image buffer. More...
 
uint16_t * GetBufferAsUInt16 ()
 Get a pointer to the image buffer. More...
 
int32_t * GetBufferAsInt32 ()
 Get a pointer to the image buffer. More...
 
uint32_t * GetBufferAsUInt32 ()
 Get a pointer to the image buffer. More...
 
int64_t * GetBufferAsInt64 ()
 Get a pointer to the image buffer. More...
 
uint64_t * GetBufferAsUInt64 ()
 Get a pointer to the image buffer. More...
 
float * GetBufferAsFloat ()
 Get a pointer to the image buffer. More...
 
double * GetBufferAsDouble ()
 Get a pointer to the image buffer. More...
 
void * GetBufferAsVoid ()
 Get a pointer to the image buffer. More...
 
const int8_t * GetBufferAsInt8 () const
 Get a pointer to the image buffer. More...
 
const uint8_t * GetBufferAsUInt8 () const
 Get a pointer to the image buffer. More...
 
const int16_t * GetBufferAsInt16 () const
 Get a pointer to the image buffer. More...
 
const uint16_t * GetBufferAsUInt16 () const
 Get a pointer to the image buffer. More...
 
const int32_t * GetBufferAsInt32 () const
 Get a pointer to the image buffer. More...
 
const uint32_t * GetBufferAsUInt32 () const
 Get a pointer to the image buffer. More...
 
const int64_t * GetBufferAsInt64 () const
 Get a pointer to the image buffer. More...
 
const uint64_t * GetBufferAsUInt64 () const
 Get a pointer to the image buffer. More...
 
const float * GetBufferAsFloat () const
 Get a pointer to the image buffer. More...
 
const double * GetBufferAsDouble () const
 Get a pointer to the image buffer. More...
 
const void * GetBufferAsVoid () const
 Get a pointer to the image buffer. More...
 

Static Public Attributes

static constexpr double DefaultImageCoordinateTolerance = 1e-6
 
static constexpr double DefaultImageDirectionTolerance = 1e-6
 

Protected Member Functions

void Allocate (const std::vector< unsigned int > &size, PixelIDValueEnum valueEnum, unsigned int numberOfComponents)
 Methods called by the constructor to allocate and initialize an image. More...
 
template<class TImageType >
void AllocateInternal (const std::vector< unsigned int > &size, unsigned int numberOfComponents)
 Allocates images of different types. More...
 
template<class TImageType >
Image ToVectorInternal (bool inPlace)
 Internal methods for converting images between vectors and scalars. More...
 
template<class TImageType >
Image ToScalarInternal (bool inPlace)
 Internal methods for converting images between vectors and scalars. More...
 

Private Member Functions

template<typename TImageType >
PimpleImageBaseDispatchedInternalInitialization (itk::DataObject *image)
 
 Image (std::unique_ptr< PimpleImageBase > pimpleImage)
 
void InternalInitialization (PixelIDValueType type, unsigned int dimension, itk::DataObject *image)
 

Private Attributes

std::unique_ptr< PimpleImageBasem_PimpleImage
 

Friends

struct AllocateMemberFunctionAddressor
 
struct DispatchedInternalInitialiationAddressor
 
struct ToScalarAddressor
 
struct ToVectorAddressor
 

Member Typedef Documentation

◆ Self

Definition at line 79 of file sitkImage.h.

Constructor & Destructor Documentation

◆ ~Image()

virtual itk::simple::Image::~Image ( )
virtual

◆ Image() [1/9]

itk::simple::Image::Image ( )

Default constructor, creates an image of size 0.

◆ Image() [2/9]

itk::simple::Image::Image ( const Image img)

◆ Image() [3/9]

itk::simple::Image::Image ( Image &&  img)
noexcept

Move constructor and assignment.

Parameters
imgAfter the operation img is valid only for destructing and assignment; all other operations have undefined behavior.

◆ Image() [4/9]

itk::simple::Image::Image ( unsigned int  width,
unsigned int  height,
PixelIDValueEnum  valueEnum 
)

Constructors for 2D, 3D an optionally 4D images where pixel type and number of components can be specified.

If the pixel type is a scalar or a label pixel type, then the number of components must be specified as 0 or 1.

If the pixel type is a vector pixel type, then the number of components defaults to the image dimension, unless the numberOfComponents is explicitly specified.

Unlike the standard convention for Dimensional Vectors the size parameter must be the exact dimension requesting. That is, it must be of length 2 of a 2D image, 3 for a 3D image and 4 for a 4D image.

◆ Image() [5/9]

itk::simple::Image::Image ( unsigned int  width,
unsigned int  height,
unsigned int  depth,
PixelIDValueEnum  valueEnum 
)

Constructors for 2D, 3D an optionally 4D images where pixel type and number of components can be specified.

If the pixel type is a scalar or a label pixel type, then the number of components must be specified as 0 or 1.

If the pixel type is a vector pixel type, then the number of components defaults to the image dimension, unless the numberOfComponents is explicitly specified.

Unlike the standard convention for Dimensional Vectors the size parameter must be the exact dimension requesting. That is, it must be of length 2 of a 2D image, 3 for a 3D image and 4 for a 4D image.

◆ Image() [6/9]

itk::simple::Image::Image ( const std::vector< unsigned int > &  size,
PixelIDValueEnum  valueEnum,
unsigned int  numberOfComponents = 0 
)

Constructors for 2D, 3D an optionally 4D images where pixel type and number of components can be specified.

If the pixel type is a scalar or a label pixel type, then the number of components must be specified as 0 or 1.

If the pixel type is a vector pixel type, then the number of components defaults to the image dimension, unless the numberOfComponents is explicitly specified.

Unlike the standard convention for Dimensional Vectors the size parameter must be the exact dimension requesting. That is, it must be of length 2 of a 2D image, 3 for a 3D image and 4 for a 4D image.

◆ Image() [7/9]

template<typename TImageType >
itk::simple::Image::Image ( itk::SmartPointer< TImageType >  image)
inlineexplicit

Construct an SimpleITK Image from an pointer to an ITK image.

The SimpleITK image will add a reference to the underlying the ITK image and hold a pointer to the image. If the image is manipulated directly from the ITK interface, SimpleITK may be unaware of it, and may cause complication related to aliasing and SimpleITK copy on write policy.

If simpleITK does not support the image type, a compile-time error or assertion will fail.

The ITK image must be fully buffered, and must have a zero starting index for the Buffered/Largest regions.

Definition at line 165 of file sitkImage.h.

◆ Image() [8/9]

template<typename TImageType >
itk::simple::Image::Image ( TImageType *  image)
inlineexplicit

Construct an SimpleITK Image from an pointer to an ITK image.

The SimpleITK image will add a reference to the underlying the ITK image and hold a pointer to the image. If the image is manipulated directly from the ITK interface, SimpleITK may be unaware of it, and may cause complication related to aliasing and SimpleITK copy on write policy.

If simpleITK does not support the image type, a compile-time error or assertion will fail.

The ITK image must be fully buffered, and must have a zero starting index for the Buffered/Largest regions.

Definition at line 170 of file sitkImage.h.

References SITK_MAX_DIMENSION, and itk::simple::sitkUnknown.

◆ Image() [9/9]

itk::simple::Image::Image ( std::unique_ptr< PimpleImageBase pimpleImage)
private

Member Function Documentation

◆ Allocate()

void itk::simple::Image::Allocate ( const std::vector< unsigned int > &  size,
PixelIDValueEnum  valueEnum,
unsigned int  numberOfComponents 
)
protected

Methods called by the constructor to allocate and initialize an image.

This method internally utilizes the member function factory to dispatch to methods instantiated on the image of the pixel ID

◆ AllocateInternal()

template<class TImageType >
void itk::simple::Image::AllocateInternal ( const std::vector< unsigned int > &  size,
unsigned int  numberOfComponents 
)
protected

Allocates images of different types.

◆ CopyInformation()

void itk::simple::Image::CopyInformation ( const Image srcImage)

Copy common meta-data from an image to this one.

Copies the Origin, Spacing, and Direction from the source image to this image. The meta-data dictionary is not copied.

It is required for the source Image's dimension and size to match, this image's attributes, otherwise an exception will be generated.

Examples
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

◆ DispatchedInternalInitialization()

template<typename TImageType >
PimpleImageBase* itk::simple::Image::DispatchedInternalInitialization ( itk::DataObject image)
private

◆ EraseMetaData()

bool itk::simple::Image::EraseMetaData ( const std::string &  key)

Remove an entry from the meta-data dictionary.

Returns true, when the value exists in the dictionary and is removed, false otherwise.

◆ EvaluateAtContinuousIndex()

std::vector<double> itk::simple::Image::EvaluateAtContinuousIndex ( const std::vector< double > &  index,
InterpolatorEnum  interp = sitkLinear 
) const

Interpolate pixel value at a continuous index.

This method is not supported for Label pixel types.

The valid range of continuous index is [-0.5, size-0.5] for each dimension. An exception is thrown if index is out of bounds.

Parameters
indexThe continuous index must be at least the length of the image dimension.
interpThe interpolation type to use, only sitkNearest and sitkLinear are supported for Vector and Complex pixel types.
Returns
All supported pixel types are returned as an array, where complex numbers are returned with the real followed by the complex component.

◆ EvaluateAtPhysicalPoint()

std::vector<double> itk::simple::Image::EvaluateAtPhysicalPoint ( const std::vector< double > &  point,
InterpolatorEnum  interp = sitkLinear 
) const

Interpolate pixel value at a physical point.

This method is not supported for Label pixel types.

An exception is thrown if the point is out of the defined region for the image.

Parameters
pointThe physical point at which the interpolation is computed.
interpThe interpolation type to use, only sitkNearest and sitkLinear are supported for Vector and Complex pixel types.
Returns
All supported pixel types are returned as an array, where complex numbers are returned with the real followed by the complex component.

◆ GetBufferAsDouble() [1/2]

double* itk::simple::Image::GetBufferAsDouble ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsDouble() [2/2]

const double* itk::simple::Image::GetBufferAsDouble ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsFloat() [1/2]

float* itk::simple::Image::GetBufferAsFloat ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsFloat() [2/2]

const float* itk::simple::Image::GetBufferAsFloat ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt16() [1/2]

int16_t* itk::simple::Image::GetBufferAsInt16 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt16() [2/2]

const int16_t* itk::simple::Image::GetBufferAsInt16 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt32() [1/2]

int32_t* itk::simple::Image::GetBufferAsInt32 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt32() [2/2]

const int32_t* itk::simple::Image::GetBufferAsInt32 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt64() [1/2]

int64_t* itk::simple::Image::GetBufferAsInt64 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt64() [2/2]

const int64_t* itk::simple::Image::GetBufferAsInt64 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt8() [1/2]

int8_t* itk::simple::Image::GetBufferAsInt8 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsInt8() [2/2]

const int8_t* itk::simple::Image::GetBufferAsInt8 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsUInt16() [1/2]

uint16_t* itk::simple::Image::GetBufferAsUInt16 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsUInt16() [2/2]

const uint16_t* itk::simple::Image::GetBufferAsUInt16 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsUInt32() [1/2]

uint32_t* itk::simple::Image::GetBufferAsUInt32 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsUInt32() [2/2]

const uint32_t* itk::simple::Image::GetBufferAsUInt32 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsUInt64() [1/2]

uint64_t* itk::simple::Image::GetBufferAsUInt64 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsUInt64() [2/2]

const uint64_t* itk::simple::Image::GetBufferAsUInt64 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsUInt8() [1/2]

uint8_t* itk::simple::Image::GetBufferAsUInt8 ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue
Examples
BufferImportExport.cxx.

◆ GetBufferAsUInt8() [2/2]

const uint8_t* itk::simple::Image::GetBufferAsUInt8 ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsVoid() [1/2]

void* itk::simple::Image::GetBufferAsVoid ( )

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetBufferAsVoid() [2/2]

const void* itk::simple::Image::GetBufferAsVoid ( ) const

Get a pointer to the image buffer.

Warning
this is dangerous

The size of the buffer is the number of components*Xsize*Ysize and then Zsize of a 3D image. The buffer should be accessed as a 1-D array. For example a 3D image buffer should be accessed:

uint8_t *buffer = img->GetBufferAsUInt8();
buffer[c + numComponents*(x+xSize*(y+ySize*z))]

The pointer to the buffer is not referenced counted. Additionally, while this image is made unique before returning the pointer, additional copying and usage may introduce unexpected aliasing of the image's buffer.

Vector and Complex pixel types are both accessed via the appropriate component type method.

The correct method for the current pixel type of the image must be called or else an exception will be generated. For vector pixel types the type of the component of the vector must be called.

See also
Image::GetPixelIDValue

◆ GetDepth()

unsigned int itk::simple::Image::GetDepth ( ) const

Get the number of pixels the Image is in the third dimension or 0 if the Image is only 2D

◆ GetDimension()

unsigned int itk::simple::Image::GetDimension ( ) const

◆ GetDirection()

std::vector<double> itk::simple::Image::GetDirection ( ) const

Set/Get the Direction.

Internally, the Direction is represented by a matrix 2x2 for a 2D and 3x3 for a 3D image. The matrix is passed as a 1D array in row-major form.

◆ GetHeight()

unsigned int itk::simple::Image::GetHeight ( ) const

Get the number of pixels the Image is in the second dimension

◆ GetITKBase() [1/2]

itk::DataObject* itk::simple::Image::GetITKBase ( )

Get access to internal ITK data object.

The return value should immediately be assigned to as itk::SmartPointer.

In many cases the value may need to be dynamically casted to the actual image type. The GetPixelIDValue() method should return an PixelID which identifies the image type which the DataObject points to.

If this object has been moved, then nullptr is returned.

Examples
ITKIntegration/ITKIntegration.cxx.

Referenced by itk::simple::ProcessObject::CastImageToITK().

◆ GetITKBase() [2/2]

const itk::DataObject* itk::simple::Image::GetITKBase ( ) const

Get access to internal ITK data object.

The return value should immediately be assigned to as itk::SmartPointer.

In many cases the value may need to be dynamically casted to the actual image type. The GetPixelIDValue() method should return an PixelID which identifies the image type which the DataObject points to.

If this object has been moved, then nullptr is returned.

◆ GetMetaData()

std::string itk::simple::Image::GetMetaData ( const std::string &  key) const

Get the value of a meta-data dictionary entry as a string.

If the key is not in the dictionary then an exception is thrown.

string types in the dictionary are returned as their native strings. Other types are printed to string before returning.

◆ GetMetaDataKeys()

std::vector<std::string> itk::simple::Image::GetMetaDataKeys ( ) const

get a vector of keys in from the meta-data dictionary

Returns a vector of keys to the key/value entries in the image's meta-data dictionary. Iterate through with these keys to get the values.

◆ GetNumberOfComponentsPerPixel()

unsigned int itk::simple::Image::GetNumberOfComponentsPerPixel ( ) const

Get the number of components for each pixel.

For images with scalar or complex pixel types this method returns one. For images with a vector pixel type the method returns the number of vector components per pixel.

◆ GetNumberOfPixels()

uint64_t itk::simple::Image::GetNumberOfPixels ( ) const

Get the number of pixels in the image.

To calculate the total number of values stored continuously for the image's buffer, the NumberOfPixels should be multiplied by NumberOfComponentsPerPixel in order to account for multiple component images.

◆ GetOrigin()

std::vector<double> itk::simple::Image::GetOrigin ( ) const

Get/Set the Origin in physical space

◆ GetPixelAsComplexFloat32()

std::complex<float> itk::simple::Image::GetPixelAsComplexFloat32 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsComplexFloat64()

std::complex<double> itk::simple::Image::GetPixelAsComplexFloat64 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsDouble()

double itk::simple::Image::GetPixelAsDouble ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsFloat()

float itk::simple::Image::GetPixelAsFloat ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsInt16()

int16_t itk::simple::Image::GetPixelAsInt16 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsInt32()

int32_t itk::simple::Image::GetPixelAsInt32 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsInt64()

int64_t itk::simple::Image::GetPixelAsInt64 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsInt8()

int8_t itk::simple::Image::GetPixelAsInt8 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsUInt16()

uint16_t itk::simple::Image::GetPixelAsUInt16 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsUInt32()

uint32_t itk::simple::Image::GetPixelAsUInt32 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsUInt64()

uint64_t itk::simple::Image::GetPixelAsUInt64 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsUInt8()

uint8_t itk::simple::Image::GetPixelAsUInt8 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorFloat32()

std::vector<float> itk::simple::Image::GetPixelAsVectorFloat32 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorFloat64()

std::vector<double> itk::simple::Image::GetPixelAsVectorFloat64 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorInt16()

std::vector<int16_t> itk::simple::Image::GetPixelAsVectorInt16 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorInt32()

std::vector<int32_t> itk::simple::Image::GetPixelAsVectorInt32 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorInt64()

std::vector<int64_t> itk::simple::Image::GetPixelAsVectorInt64 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorInt8()

std::vector<int8_t> itk::simple::Image::GetPixelAsVectorInt8 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorUInt16()

std::vector<uint16_t> itk::simple::Image::GetPixelAsVectorUInt16 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorUInt32()

std::vector<uint32_t> itk::simple::Image::GetPixelAsVectorUInt32 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorUInt64()

std::vector<uint64_t> itk::simple::Image::GetPixelAsVectorUInt64 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelAsVectorUInt8()

std::vector<uint8_t> itk::simple::Image::GetPixelAsVectorUInt8 ( const std::vector< uint32_t > &  idx) const

Get the value of a pixel.

Returns the value of a pixel for the given index. The index follows standard SimpleITK conventions for it's length. The correct method must be called for the underlying Image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
See also
Image::GetPixelIDValue

◆ GetPixelID()

PixelIDValueEnum itk::simple::Image::GetPixelID ( ) const

Get the pixel type

The pixel type is set at construction type and can not be manually changed, unless by assignment. The value may be -1 or "Unknown".

Examples
BufferImportExport.cxx, DemonsRegistration2/DemonsRegistration2.cxx, FilterProgressReporting/FilterProgressReporting.cxx, SimpleGaussian/SimpleGaussian.cxx, and SimpleGaussianFunctional.cxx.

◆ GetPixelIDTypeAsString()

std::string itk::simple::Image::GetPixelIDTypeAsString ( ) const

Return the pixel type as a human readable string value.

Referenced by itk::simple::ProcessObject::CastImageToITK().

◆ GetPixelIDValue()

PixelIDValueType itk::simple::Image::GetPixelIDValue ( ) const

◆ GetSize()

std::vector<unsigned int> itk::simple::Image::GetSize ( ) const

Get the number of pixels the Image is in each dimension as a std::vector. The size of the vector is equal to the number of dimensions for the image.

Examples
DicomSeriesReader/DicomSeriesReader.cxx, ImageIOSelection/ImageIOSelection.cxx, ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx, and SliceBySliceDecorator/SliceBySliceDecorator.cxx.

◆ GetSizeOfPixelComponent()

unsigned int itk::simple::Image::GetSizeOfPixelComponent ( ) const

Get the number of bytes per component of a pixel.

Returns the sizeof the pixel component type.

◆ GetSpacing()

std::vector<double> itk::simple::Image::GetSpacing ( ) const

Get/Set the Spacing of the Image as an std::vector .

The spacing describes the physical size of each pixel. The length of the vector is equal to the dimension of the Image.

◆ GetWidth()

unsigned int itk::simple::Image::GetWidth ( ) const

Get the number of pixels the Image is in the first dimension

◆ HasMetaDataKey()

bool itk::simple::Image::HasMetaDataKey ( const std::string &  key) const

Query the meta-data dictionary for the existence of a key.

◆ InternalInitialization()

void itk::simple::Image::InternalInitialization ( PixelIDValueType  type,
unsigned int  dimension,
itk::DataObject image 
)
private

Method called by certain constructors to convert ITK images into SimpleITK ones.

This is the single method which needs to be explicitly instantiated to separate the internal ITK and Pimple image from the external SimpleITK interface. Template parameters have been chosen carefully to flexibly enable this.

◆ IsCongruentImageGeometry()

bool itk::simple::Image::IsCongruentImageGeometry ( const Image otherImage,
double  coordinateTolerance,
double  directionTolerance 
) const

Checks whether the images' pixels at the same index occupy the same physical space.

Compares the origin, spacing, and direction for equality within provided tolerances. There is no check for matching regions in between the images.

If The dimensions of the images do not match, false is returned.

◆ IsSameImageGeometryAs()

bool itk::simple::Image::IsSameImageGeometryAs ( const Image otherImage,
double  = DefaultImageCoordinateTolerance,
double  = DefaultImageDirectionTolerance 
) const

Check whether the images have the same grid in physical space.

Compares largest possible regions for equality, and the origin, spacing, and direction cosines for equality within provided tolerances.

If the dimensions of the images do not match, false is returned.

◆ IsUnique()

bool itk::simple::Image::IsUnique ( ) const

Returns true if no other SimpleITK Image object refers to the same internal data structure.

◆ MakeUnique()

void itk::simple::Image::MakeUnique ( )

Performs actually coping if needed to make object unique.

The Image class by default performs lazy coping and assignment. This method make sure that coping actually happens to the itk::Image pointed to is only pointed to by this object.

◆ operator=() [1/2]

Image& itk::simple::Image::operator= ( const Image img)

◆ operator=() [2/2]

Image& itk::simple::Image::operator= ( Image &&  img)
noexcept

◆ ProxyForInPlaceOperation()

Image itk::simple::Image::ProxyForInPlaceOperation ( )

Advanced method not commonly needed.

This method is designed to support implementations "in-place" object behavior for methods which operate on r-value references. The returned image is a new image which has a low level pointer to this object's image buffer, without the SimpleITK or ITK reference counting. This is implemented by setting the new ITK Image's buffer to the same as this objects without ownership.

Warning
This method bypasses the SimpleITK reference counting, and the reference needs to be manually maintained in the scope. The resulting object is designed only to be a temporary.

In the following example this method is used instead of an std::move call when the filter's first argument takes an r-value reference. The img object will container the results of the filter execution, and the img image buffer will be preserved in case of exceptions, and the meta-data will remain in the img object.

filter.Execute( img.ProxyForInPlaceOperation() );

The meta-data dictionary is not copied to the returned proxy image.

Referenced by itk::simple::operator%=(), itk::simple::operator&=(), itk::simple::operator*=(), itk::simple::operator+=(), itk::simple::operator-=(), itk::simple::operator/=(), itk::simple::operator^=(), and itk::simple::operator|=().

◆ SetDirection()

void itk::simple::Image::SetDirection ( const std::vector< double > &  direction)

Set/Get the Direction.

Internally, the Direction is represented by a matrix 2x2 for a 2D and 3x3 for a 3D image. The matrix is passed as a 1D array in row-major form.

◆ SetMetaData()

void itk::simple::Image::SetMetaData ( const std::string &  key,
const std::string &  value 
)

Set an entry in the meta-data dictionary.

Replaces or creates an entry in the image's meta-data dictionary.

◆ SetOrigin()

void itk::simple::Image::SetOrigin ( const std::vector< double > &  origin)

Get/Set the Origin in physical space

◆ SetPixelAsComplexFloat32()

void itk::simple::Image::SetPixelAsComplexFloat32 ( const std::vector< uint32_t > &  idx,
const std::complex< float >  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsComplexFloat64()

void itk::simple::Image::SetPixelAsComplexFloat64 ( const std::vector< uint32_t > &  idx,
const std::complex< double >  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsDouble()

void itk::simple::Image::SetPixelAsDouble ( const std::vector< uint32_t > &  idx,
double  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsFloat()

void itk::simple::Image::SetPixelAsFloat ( const std::vector< uint32_t > &  idx,
float  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsInt16()

void itk::simple::Image::SetPixelAsInt16 ( const std::vector< uint32_t > &  idx,
int16_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsInt32()

void itk::simple::Image::SetPixelAsInt32 ( const std::vector< uint32_t > &  idx,
int32_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsInt64()

void itk::simple::Image::SetPixelAsInt64 ( const std::vector< uint32_t > &  idx,
int64_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsInt8()

void itk::simple::Image::SetPixelAsInt8 ( const std::vector< uint32_t > &  idx,
int8_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsUInt16()

void itk::simple::Image::SetPixelAsUInt16 ( const std::vector< uint32_t > &  idx,
uint16_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsUInt32()

void itk::simple::Image::SetPixelAsUInt32 ( const std::vector< uint32_t > &  idx,
uint32_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsUInt64()

void itk::simple::Image::SetPixelAsUInt64 ( const std::vector< uint32_t > &  idx,
uint64_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsUInt8()

void itk::simple::Image::SetPixelAsUInt8 ( const std::vector< uint32_t > &  idx,
uint8_t  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorFloat32()

void itk::simple::Image::SetPixelAsVectorFloat32 ( const std::vector< uint32_t > &  idx,
const std::vector< float > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorFloat64()

void itk::simple::Image::SetPixelAsVectorFloat64 ( const std::vector< uint32_t > &  idx,
const std::vector< double > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorInt16()

void itk::simple::Image::SetPixelAsVectorInt16 ( const std::vector< uint32_t > &  idx,
const std::vector< int16_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorInt32()

void itk::simple::Image::SetPixelAsVectorInt32 ( const std::vector< uint32_t > &  idx,
const std::vector< int32_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorInt64()

void itk::simple::Image::SetPixelAsVectorInt64 ( const std::vector< uint32_t > &  idx,
const std::vector< int64_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorInt8()

void itk::simple::Image::SetPixelAsVectorInt8 ( const std::vector< uint32_t > &  idx,
const std::vector< int8_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorUInt16()

void itk::simple::Image::SetPixelAsVectorUInt16 ( const std::vector< uint32_t > &  idx,
const std::vector< uint16_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorUInt32()

void itk::simple::Image::SetPixelAsVectorUInt32 ( const std::vector< uint32_t > &  idx,
const std::vector< uint32_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorUInt64()

void itk::simple::Image::SetPixelAsVectorUInt64 ( const std::vector< uint32_t > &  idx,
const std::vector< uint64_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetPixelAsVectorUInt8()

void itk::simple::Image::SetPixelAsVectorUInt8 ( const std::vector< uint32_t > &  idx,
const std::vector< uint8_t > &  v 
)

Set the value of a pixel.

Set the value of a pixel at the provided index. The index follows standard SimpleITK conventions for it's length. The correct method must be called which should match the underlying image type, otherwise an exception will be thrown.

Parameters
idxthe zero based index into the image. It's length must be at least the value of GetDimension(), additional elements will be ignored. Boundary checking is performed on idx, if it is out of bounds an exception will be thrown.
vvalue to set the pixel to
See also
Image::GetPixelIDValue

◆ SetSpacing()

void itk::simple::Image::SetSpacing ( const std::vector< double > &  spacing)

Get/Set the Spacing of the Image as an std::vector .

The spacing describes the physical size of each pixel. The length of the vector is equal to the dimension of the Image.

◆ ToScalarImage()

Image itk::simple::Image::ToScalarImage ( bool  inPlace = true)

Convert a image of vector pixel type to a scalar image with N+1 dimensions.

This method will convert a vector image to a scalar image with the size of the first dimension equal to the number of components. If the image is already a scalar image then the image is returned.

For the additional dimension the origin is set to zero, the spacing to one, and the new components of the direction cosine to the identity matrix.

An exception is thrown if the image is has SITK_MAX_DIMENSION dimensions or if the pixel type is a label or complex pixel type.

Parameters
inPlaceIf true then the image is made unique and converted in place updating this image, otherwise a copy of the image is made and returned.
See also
ToVectorImage

◆ ToScalarInternal()

template<class TImageType >
Image itk::simple::Image::ToScalarInternal ( bool  inPlace)
protected

Internal methods for converting images between vectors and scalars.

◆ ToString()

std::string itk::simple::Image::ToString ( ) const

◆ ToVectorImage()

Image itk::simple::Image::ToVectorImage ( bool  inPlace = true)

Convert the first dimension to the components for image with vector pixel type.

This method will convert a scalar image to a vector image with the number of components equal to the size of the first dimension. If the image is already a vector image then the image is returned.

The components of the direction cosine matrix for the first dimension must be the identity matrix, or else an exception is thrown.

An exception is thrown if the image is 2D or if the pixel type is a label or complex pixel type.

Parameters
inPlaceIf true then the image is made unique and converted in place updating this image, otherwise a copy of the image is made and returned.
See also
ToScalarImage

◆ ToVectorInternal()

template<class TImageType >
Image itk::simple::Image::ToVectorInternal ( bool  inPlace)
protected

Internal methods for converting images between vectors and scalars.

◆ TransformContinuousIndexToPhysicalPoint()

std::vector<double> itk::simple::Image::TransformContinuousIndexToPhysicalPoint ( const std::vector< double > &  index) const

Transform continuous index to physical point

◆ TransformIndexToPhysicalPoint()

std::vector<double> itk::simple::Image::TransformIndexToPhysicalPoint ( const std::vector< int64_t > &  index) const

Transform index to physical point

◆ TransformPhysicalPointToContinuousIndex()

std::vector<double> itk::simple::Image::TransformPhysicalPointToContinuousIndex ( const std::vector< double > &  point) const

Transform physical point to continuous index

◆ TransformPhysicalPointToIndex()

std::vector<int64_t> itk::simple::Image::TransformPhysicalPointToIndex ( const std::vector< double > &  point) const

Transform physical point to index

Friends And Related Function Documentation

◆ AllocateMemberFunctionAddressor

friend struct AllocateMemberFunctionAddressor
friend

Definition at line 752 of file sitkImage.h.

◆ DispatchedInternalInitialiationAddressor

friend struct DispatchedInternalInitialiationAddressor
friend

Definition at line 751 of file sitkImage.h.

◆ ToScalarAddressor

friend struct ToScalarAddressor
friend

Definition at line 754 of file sitkImage.h.

◆ ToVectorAddressor

friend struct ToVectorAddressor
friend

Definition at line 753 of file sitkImage.h.

Member Data Documentation

◆ DefaultImageCoordinateTolerance

constexpr double itk::simple::Image::DefaultImageCoordinateTolerance = 1e-6
staticconstexpr

Definition at line 694 of file sitkImage.h.

◆ DefaultImageDirectionTolerance

constexpr double itk::simple::Image::DefaultImageDirectionTolerance = 1e-6
staticconstexpr

Definition at line 695 of file sitkImage.h.

◆ m_PimpleImage

std::unique_ptr<PimpleImageBase> itk::simple::Image::m_PimpleImage
private

Definition at line 757 of file sitkImage.h.


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