DGtal  1.1.0
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes
DGtal::ImageContainerByITKImage< TDomain, TValue > Class Template Reference

Aim: implements a model of CImageContainer using a ITK Image. More...

#include <DGtal/images/ImageContainerByITKImage.h>

Collaboration diagram for DGtal::ImageContainerByITKImage< TDomain, TValue >:
[legend]

Public Types

typedef TValue Value
 
typedef TDomain Domain
 
typedef ImageContainerByITKImage< TDomain, TValue > Self
 
typedef Domain::Point Point
 
typedef Domain::Vector Vector
 
typedef Domain::Dimension Dimension
 
typedef Domain::Integer Integer
 
typedef Domain::Size Size
 
typedef Point Vertex
 
typedef PointVector< dimension, double > RealPoint
 
typedef RealPoint PhysicalPoint
 
typedef itk::Image< TValue, dimensionITKImage
 
typedef ITKImage::Pointer ITKImagePointer
 
typedef itk::ImageRegionConstIterator< ITKImageConstIterator
 
typedef itk::ImageRegionIterator< ITKImageIterator
 
typedef DefaultConstImageRange< SelfConstRange
 
typedef DefaultImageRange< SelfRange
 

Public Member Functions

 BOOST_CONCEPT_ASSERT ((concepts::CLabel< TValue >))
 
 BOOST_CONCEPT_ASSERT ((concepts::CDomain< TDomain >))
 
 ImageContainerByITKImage (const Domain &aDomain)
 
 ImageContainerByITKImage (const ITKImagePointer &aRef)
 
 ImageContainerByITKImage (const ImageContainerByITKImage &other)
 
ImageContainerByITKImageoperator= (const ImageContainerByITKImage &other)
 
 ~ImageContainerByITKImage ()
 
void updateDomain (const Point &inputDomainShift=Point())
 
ConstRange constRange () const
 
Range range ()
 
Value operator() (const Point &domainPoint) const
 
Value operator() (const ConstIterator &it) const
 
Value operator() (const Iterator &it) const
 
void setValue (const Point &domainPoint, const Value &aValue)
 
void setValue (Iterator &it, const Value &V)
 
const Domaindomain () const
 
ITKImagePointer getITKImagePointer () const
 
const PointgetDomainShift () const
 
void selfDisplay (std::ostream &out) const
 
bool isValid () const
 
ConstIterator begin () const
 
Iterator begin ()
 
const ConstIterator end () const
 
Iterator end ()
 
Point getDomainPointFromIndex (const Point &indexPoint) const
 
Point getIndexFromDomainPoint (const Point &domainPoint) const
 
Point getDomainPointFromItkIndex (const typename ITKImage::IndexType &itkIndexPoint) const
 
ITKImage::IndexType getItkIndexFromDomainPoint (const Point &domainPoint) const
 
PhysicalPoint getPhysicalPointFromDomainPoint (const Point &domainPoint) const
 
Point getDomainPointFromPhysicalPoint (const PhysicalPoint &physicalPoint) const
 
PhysicalPoint getLowerBoundAsPhysicalPoint () const
 
PhysicalPoint getUpperBoundAsPhysicalPoint () const
 

Static Public Attributes

static const Domain::Dimension dimension = Domain::dimension
 

Protected Member Functions

 ImageContainerByITKImage ()
 

Private Attributes

ITKImagePointer myITKImagePointer
 
Domain myDomain
 
Point myDomainShift = Point()
 

Detailed Description

template<typename TDomain, typename TValue>
class DGtal::ImageContainerByITKImage< TDomain, TValue >

Aim: implements a model of CImageContainer using a ITK Image.

Description of template class 'ImageContainerByITKImage'

Using this container, you can switch from DGtal alogrithms to ITK processing pipeline. The Ownership of the underlying ITK image is shared between the wrapper and the ITK pipeline. If the ITK image region is modified, one should manually update the domain of the wrapper. This is done by calling the updateDomain() method.

See also
testITKImage.cpp

Definition at line 93 of file ImageContainerByITKImage.h.

Member Typedef Documentation

◆ ConstIterator

template<typename TDomain , typename TValue >
typedef itk::ImageRegionConstIterator< ITKImage > DGtal::ImageContainerByITKImage< TDomain, TValue >::ConstIterator

Definition at line 119 of file ImageContainerByITKImage.h.

◆ ConstRange

template<typename TDomain , typename TValue >
typedef DefaultConstImageRange<Self> DGtal::ImageContainerByITKImage< TDomain, TValue >::ConstRange

Definition at line 122 of file ImageContainerByITKImage.h.

◆ Dimension

template<typename TDomain , typename TValue >
typedef Domain::Dimension DGtal::ImageContainerByITKImage< TDomain, TValue >::Dimension

Definition at line 110 of file ImageContainerByITKImage.h.

◆ Domain

template<typename TDomain , typename TValue >
typedef TDomain DGtal::ImageContainerByITKImage< TDomain, TValue >::Domain

Definition at line 102 of file ImageContainerByITKImage.h.

◆ Integer

template<typename TDomain , typename TValue >
typedef Domain::Integer DGtal::ImageContainerByITKImage< TDomain, TValue >::Integer

Definition at line 111 of file ImageContainerByITKImage.h.

◆ Iterator

template<typename TDomain , typename TValue >
typedef itk::ImageRegionIterator< ITKImage > DGtal::ImageContainerByITKImage< TDomain, TValue >::Iterator

Definition at line 120 of file ImageContainerByITKImage.h.

◆ ITKImage

template<typename TDomain , typename TValue >
typedef itk::Image< TValue, dimension> DGtal::ImageContainerByITKImage< TDomain, TValue >::ITKImage

Definition at line 117 of file ImageContainerByITKImage.h.

◆ ITKImagePointer

template<typename TDomain , typename TValue >
typedef ITKImage::Pointer DGtal::ImageContainerByITKImage< TDomain, TValue >::ITKImagePointer

Definition at line 118 of file ImageContainerByITKImage.h.

◆ PhysicalPoint

template<typename TDomain , typename TValue >
typedef RealPoint DGtal::ImageContainerByITKImage< TDomain, TValue >::PhysicalPoint

Definition at line 115 of file ImageContainerByITKImage.h.

◆ Point

template<typename TDomain , typename TValue >
typedef Domain::Point DGtal::ImageContainerByITKImage< TDomain, TValue >::Point

Definition at line 108 of file ImageContainerByITKImage.h.

◆ Range

template<typename TDomain , typename TValue >
typedef DefaultImageRange<Self> DGtal::ImageContainerByITKImage< TDomain, TValue >::Range

Definition at line 123 of file ImageContainerByITKImage.h.

◆ RealPoint

template<typename TDomain , typename TValue >
typedef PointVector<dimension, double> DGtal::ImageContainerByITKImage< TDomain, TValue >::RealPoint

Definition at line 114 of file ImageContainerByITKImage.h.

◆ Self

template<typename TDomain , typename TValue >
typedef ImageContainerByITKImage<TDomain, TValue> DGtal::ImageContainerByITKImage< TDomain, TValue >::Self

Definition at line 103 of file ImageContainerByITKImage.h.

◆ Size

template<typename TDomain , typename TValue >
typedef Domain::Size DGtal::ImageContainerByITKImage< TDomain, TValue >::Size

Definition at line 112 of file ImageContainerByITKImage.h.

◆ Value

template<typename TDomain , typename TValue >
typedef TValue DGtal::ImageContainerByITKImage< TDomain, TValue >::Value

Definition at line 101 of file ImageContainerByITKImage.h.

◆ Vector

template<typename TDomain , typename TValue >
typedef Domain::Vector DGtal::ImageContainerByITKImage< TDomain, TValue >::Vector

Definition at line 109 of file ImageContainerByITKImage.h.

◆ Vertex

template<typename TDomain , typename TValue >
typedef Point DGtal::ImageContainerByITKImage< TDomain, TValue >::Vertex

Definition at line 113 of file ImageContainerByITKImage.h.

Constructor & Destructor Documentation

◆ ImageContainerByITKImage() [1/4]

template<typename TDomain , typename TValue >
DGtal::ImageContainerByITKImage< TDomain, TValue >::ImageContainerByITKImage ( const Domain aDomain)

Constructor.

Parameters
aDomainthe image domain.

◆ ImageContainerByITKImage() [2/4]

template<typename TDomain , typename TValue >
DGtal::ImageContainerByITKImage< TDomain, TValue >::ImageContainerByITKImage ( const ITKImagePointer aRef)

Constructor.

Parameters
aRefa reference to an ITKImage

◆ ImageContainerByITKImage() [3/4]

template<typename TDomain , typename TValue >
DGtal::ImageContainerByITKImage< TDomain, TValue >::ImageContainerByITKImage ( const ImageContainerByITKImage< TDomain, TValue > &  other)

Copy constructor

Parameters
otherthe object to copy.

◆ ~ImageContainerByITKImage()

template<typename TDomain , typename TValue >
DGtal::ImageContainerByITKImage< TDomain, TValue >::~ImageContainerByITKImage ( )

Destructor.

◆ ImageContainerByITKImage() [4/4]

template<typename TDomain , typename TValue >
DGtal::ImageContainerByITKImage< TDomain, TValue >::ImageContainerByITKImage ( )
protected

Constructor. Forbidden by default (protected to avoid g++ warnings).

Member Function Documentation

◆ begin() [1/2]

template<typename TDomain , typename TValue >
Iterator DGtal::ImageContainerByITKImage< TDomain, TValue >::begin ( )
inline

begin() const iterator.

Definition at line 291 of file ImageContainerByITKImage.h.

292  {
293  Iterator iter = Iterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
294  iter.GoToBegin();
295  return iter;
296  }

References DGtal::ImageContainerByITKImage< TDomain, TValue >::myITKImagePointer.

◆ begin() [2/2]

template<typename TDomain , typename TValue >
ConstIterator DGtal::ImageContainerByITKImage< TDomain, TValue >::begin ( ) const
inline

begin() const iterator.

Definition at line 279 of file ImageContainerByITKImage.h.

280  {
281  ConstIterator iter = ConstIterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
282  iter.GoToBegin();
283  return iter;
284  }

References DGtal::ImageContainerByITKImage< TDomain, TValue >::myITKImagePointer.

◆ BOOST_CONCEPT_ASSERT() [1/2]

template<typename TDomain , typename TValue >
DGtal::ImageContainerByITKImage< TDomain, TValue >::BOOST_CONCEPT_ASSERT ( (concepts::CDomain< TDomain >)  )

◆ BOOST_CONCEPT_ASSERT() [2/2]

template<typename TDomain , typename TValue >
DGtal::ImageContainerByITKImage< TDomain, TValue >::BOOST_CONCEPT_ASSERT ( (concepts::CLabel< TValue >)  )

◆ constRange()

template<typename TDomain , typename TValue >
ConstRange DGtal::ImageContainerByITKImage< TDomain, TValue >::constRange ( ) const
inline
Returns
the range providing begin and end iterators to scan the values of image.

Definition at line 179 of file ImageContainerByITKImage.h.

180  {
181  return ConstRange(*this);
182  }

◆ domain()

template<typename TDomain , typename TValue >
const Domain& DGtal::ImageContainerByITKImage< TDomain, TValue >::domain ( ) const
inline
Returns
the domain associated to the image.

Definition at line 239 of file ImageContainerByITKImage.h.

240  {
241  return myDomain;
242  }

References DGtal::ImageContainerByITKImage< TDomain, TValue >::myDomain.

◆ end() [1/2]

template<typename TDomain , typename TValue >
Iterator DGtal::ImageContainerByITKImage< TDomain, TValue >::end ( )
inline

end() iterator.

Definition at line 315 of file ImageContainerByITKImage.h.

316  {
317  Iterator iter = Iterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
318  iter.GoToEnd();
319  return iter;
320  }

References DGtal::ImageContainerByITKImage< TDomain, TValue >::myITKImagePointer.

◆ end() [2/2]

template<typename TDomain , typename TValue >
const ConstIterator DGtal::ImageContainerByITKImage< TDomain, TValue >::end ( ) const
inline

end() const iterator.

Definition at line 303 of file ImageContainerByITKImage.h.

304  {
305  ConstIterator iter = ConstIterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
306  iter.GoToEnd();
307  return iter;
308  }

References DGtal::ImageContainerByITKImage< TDomain, TValue >::myITKImagePointer.

◆ getDomainPointFromIndex()

template<typename TDomain , typename TValue >
Point DGtal::ImageContainerByITKImage< TDomain, TValue >::getDomainPointFromIndex ( const Point indexPoint) const
inline

IndexPoint refers to a valid ITKImage::IndexPoint DomainPoint refers to a Point between lowerBound and upperBound of myDomain. They are different only when myDomainShift is different than Zero.

Parameters
indexPointa point holding valid index coordinates of the ITK image.
Returns
domainPoint a point between lowerBound and upperBound

◆ getDomainPointFromItkIndex()

template<typename TDomain , typename TValue >
Point DGtal::ImageContainerByITKImage< TDomain, TValue >::getDomainPointFromItkIndex ( const typename ITKImage::IndexType &  itkIndexPoint) const
inline

The same as getDomainPointFromItkIndex and getIndexFromDomainPoint but using ITK types.

Parameters
itkIndexPointan IndexType of ITK.
Returns
domainPoint a point between lowerBound and upperBound

◆ getDomainPointFromPhysicalPoint()

template<typename TDomain , typename TValue >
Point DGtal::ImageContainerByITKImage< TDomain, TValue >::getDomainPointFromPhysicalPoint ( const PhysicalPoint physicalPoint) const
inline

◆ getDomainShift()

template<typename TDomain , typename TValue >
const Point& DGtal::ImageContainerByITKImage< TDomain, TValue >::getDomainShift ( ) const
inline

Definition at line 254 of file ImageContainerByITKImage.h.

255  {
256  return myDomainShift;
257  }

References DGtal::ImageContainerByITKImage< TDomain, TValue >::myDomainShift.

◆ getIndexFromDomainPoint()

template<typename TDomain , typename TValue >
Point DGtal::ImageContainerByITKImage< TDomain, TValue >::getIndexFromDomainPoint ( const Point domainPoint) const
inline

◆ getITKImagePointer()

template<typename TDomain , typename TValue >
ITKImagePointer DGtal::ImageContainerByITKImage< TDomain, TValue >::getITKImagePointer ( ) const
inline

Returns a copy of the itkImage smartPointer

Definition at line 248 of file ImageContainerByITKImage.h.

249  {
250  return myITKImagePointer;
251  }

References DGtal::ImageContainerByITKImage< TDomain, TValue >::myITKImagePointer.

◆ getItkIndexFromDomainPoint()

template<typename TDomain , typename TValue >
ITKImage::IndexType DGtal::ImageContainerByITKImage< TDomain, TValue >::getItkIndexFromDomainPoint ( const Point domainPoint) const
inline

◆ getLowerBoundAsPhysicalPoint()

template<typename TDomain , typename TValue >
PhysicalPoint DGtal::ImageContainerByITKImage< TDomain, TValue >::getLowerBoundAsPhysicalPoint ( ) const
inline

Returns the lower and upper bounds as physical points. Note that the location of pixels in ITK are in the center of the square. But in DGtal a lowerBound is in the south-west corner of that square. And an upperBound is the north-east corner of the last pixel.

Returns
physical point of the location of the start index of the region.

◆ getPhysicalPointFromDomainPoint()

template<typename TDomain , typename TValue >
PhysicalPoint DGtal::ImageContainerByITKImage< TDomain, TValue >::getPhysicalPointFromDomainPoint ( const Point domainPoint) const
inline

Get PhysicalPoints from a domain point in DGtal and viceversa.

Remember that GetOrigin() in ITK is the physical location of the index {0,0...}. Not the location of the start index of the region.

Parameters
domainPointa point holding a point in the DGtal domain. It will be converted to a valid index of the ITK image, taking into account the value of myDomainShift.
Returns
physical point of the index.

◆ getUpperBoundAsPhysicalPoint()

template<typename TDomain , typename TValue >
PhysicalPoint DGtal::ImageContainerByITKImage< TDomain, TValue >::getUpperBoundAsPhysicalPoint ( ) const
inline

◆ isValid()

template<typename TDomain , typename TValue >
bool DGtal::ImageContainerByITKImage< TDomain, TValue >::isValid ( ) const

Checks the validity/consistency of the object.

Returns
'true' if the object is valid, 'false' otherwise.

◆ operator()() [1/3]

template<typename TDomain , typename TValue >
Value DGtal::ImageContainerByITKImage< TDomain, TValue >::operator() ( const ConstIterator it) const

Get the value of an image at a given position.

Parameters
itposition in the image.
Returns
the value of the point pointed by the iterator.

◆ operator()() [2/3]

template<typename TDomain , typename TValue >
Value DGtal::ImageContainerByITKImage< TDomain, TValue >::operator() ( const Iterator it) const

Get the value of an image at a given position.

Parameters
itposition in the image.
Returns
the value of the point pointed by the iterator.

◆ operator()() [3/3]

template<typename TDomain , typename TValue >
Value DGtal::ImageContainerByITKImage< TDomain, TValue >::operator() ( const Point domainPoint) const

Get the value of an image at a given position.

Parameters
domainPointposition in the image.
Returns
the value at indexPoint.

◆ operator=()

template<typename TDomain , typename TValue >
ImageContainerByITKImage& DGtal::ImageContainerByITKImage< TDomain, TValue >::operator= ( const ImageContainerByITKImage< TDomain, TValue > &  other)

Assignment.

Parameters
otherthe object to copy.
Returns
a reference on 'this'.

◆ range()

template<typename TDomain , typename TValue >
Range DGtal::ImageContainerByITKImage< TDomain, TValue >::range ( )
inline
Returns
the range providing begin and end iterators to scan the values of image.

Definition at line 188 of file ImageContainerByITKImage.h.

189  {
190  return Range(*this);
191  }

◆ selfDisplay()

template<typename TDomain , typename TValue >
void DGtal::ImageContainerByITKImage< TDomain, TValue >::selfDisplay ( std::ostream &  out) const

Writes/Displays the object on an output stream.

Parameters
outthe output stream where the object is written.

◆ setValue() [1/2]

template<typename TDomain , typename TValue >
void DGtal::ImageContainerByITKImage< TDomain, TValue >::setValue ( const Point domainPoint,
const Value aValue 
)

Set a value on an Image at domainPoint.

Parameters
domainPointlocation of the point to associate with aValue.
aValuethe value.

◆ setValue() [2/2]

template<typename TDomain , typename TValue >
void DGtal::ImageContainerByITKImage< TDomain, TValue >::setValue ( Iterator it,
const Value V 
)

Set a value on an Image at a given position

Parameters
itlocation of the point (Iterator) to associate with aValue.
Vthe value.

◆ updateDomain()

template<typename TDomain , typename TValue >
void DGtal::ImageContainerByITKImage< TDomain, TValue >::updateDomain ( const Point inputDomainShift = Point())

update internal domain cache. should be called after modifying underlying ITK image or to set myDomainShift.

Parameters
inputDomainShiftapplies a domainShift to the lowerBound and upperBound of myDomain.
See also
myDomainShift

Field Documentation

◆ dimension

template<typename TDomain , typename TValue >
const Domain::Dimension DGtal::ImageContainerByITKImage< TDomain, TValue >::dimension = Domain::dimension
static

Definition at line 106 of file ImageContainerByITKImage.h.

◆ myDomain

template<typename TDomain , typename TValue >
Domain DGtal::ImageContainerByITKImage< TDomain, TValue >::myDomain
private

◆ myDomainShift

template<typename TDomain , typename TValue >
Point DGtal::ImageContainerByITKImage< TDomain, TValue >::myDomainShift = Point()
private

Apply a shift in the lower and upper bounds. Useful to represent multiple images, or images with physical information, i.e with non-default Origin (see ITK for Origin, Spacing and Direction metadata)

Set it using the function updateDomain(inputDomainShift)

Please note that this is a workaround, DGtal cannot fully work in the Physical Domain. The spacing in DGtal is fixed to 1. This allows to work to represent multiple images with different Origins, but they need to have the same Spacing and Direction for the visualization to be meaningful.

Default to zero.

When is not zero, use getIndexFromDomainPoint and getDomainPointFromIndex to switch between points in the domain and indices in the itk image.

Definition at line 412 of file ImageContainerByITKImage.h.

Referenced by DGtal::ImageContainerByITKImage< TDomain, TValue >::getDomainShift().

◆ myITKImagePointer

template<typename TDomain , typename TValue >
ITKImagePointer DGtal::ImageContainerByITKImage< TDomain, TValue >::myITKImagePointer
private

The documentation for this class was generated from the following file:
ConstIterator
MyDigitalSurface::ConstIterator ConstIterator
Definition: greedy-plane-segmentation-ex2.cpp:93
DGtal::ImageContainerByITKImage::myITKImagePointer
ITKImagePointer myITKImagePointer
Definition: ImageContainerByITKImage.h:393
DGtal::ImageContainerByITKImage::myDomain
Domain myDomain
Definition: ImageContainerByITKImage.h:394
DGtal::ImageContainerByITKImage::myDomainShift
Point myDomainShift
Definition: ImageContainerByITKImage.h:412
DGtal::ImageContainerByITKImage::Range
DefaultImageRange< Self > Range
Definition: ImageContainerByITKImage.h:123
DGtal::ImageContainerByITKImage::ConstIterator
itk::ImageRegionConstIterator< ITKImage > ConstIterator
Definition: ImageContainerByITKImage.h:119
DGtal::ImageContainerByITKImage::ConstRange
DefaultConstImageRange< Self > ConstRange
Definition: ImageContainerByITKImage.h:122
DGtal::ImageContainerByITKImage::Iterator
itk::ImageRegionIterator< ITKImage > Iterator
Definition: ImageContainerByITKImage.h:120