DGtal 1.4.0
Loading...
Searching...
No Matches
ImageContainerByITKImage.h
1
17#pragma once
18
33#if defined(ImageContainerByITKImage_RECURSES)
34#error Recursive header files inclusion detected in ImageContainerByITKImage.h
35#else // defined(ImageContainerByITKImage_RECURSES)
37#define ImageContainerByITKImage_RECURSES
38
39#if !defined ImageContainerByITKImage_h
41#define ImageContainerByITKImage_h
42
44// Inclusions
45#include "DGtal/base/Common.h"
46#include "DGtal/base/CLabel.h"
47#include "DGtal/kernel/domains/CDomain.h"
48#include "DGtal/kernel/domains/CDomain.h"
49#include "DGtal/kernel/PointVector.h"
50#include "DGtal/images/DefaultConstImageRange.h"
51#include "DGtal/images/DefaultImageRange.h"
52
53#if defined(__GNUG__)
54#pragma GCC diagnostic push
55#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
56#pragma GCC diagnostic ignored "-Wpedantic"
57#endif
58#if defined(__clang__)
59#pragma clang diagnostic push
60#pragma clang diagnostic ignored "-Wdocumentation"
61#endif
62#include <itkImage.h>
63#include <itkImageRegionConstIterator.h>
64#include <itkImageRegionIterator.h>
65#include <iostream>
66#if defined(__clang__)
67#pragma clang diagnostic pop
68#endif
69#if defined(__GNUG__)
70#pragma GCC diagnostic pop
71#endif
73
74namespace DGtal
75{
76
78 // template class ImageContainerByITKImage
92 template <typename TDomain, typename TValue>
94 {
95 // ----------------------- Standard services ------------------------------
96 public:
97
100
101 typedef TValue Value;
102 typedef TDomain Domain;
104
105 // static constants
106 static const typename Domain::Dimension dimension = Domain::dimension;
107
108 typedef typename Domain::Point Point;
109 typedef typename Domain::Vector Vector;
110 typedef typename Domain::Dimension Dimension;
111 typedef typename Domain::Integer Integer;
112 typedef typename Domain::Size Size;
113 typedef Point Vertex;
116
117 typedef typename itk::Image< TValue, dimension> ITKImage;
118 typedef typename itk::ImageBase<dimension>::SpacingValueType ITKSpacingValueType;
120 typedef typename ITKImage::Pointer ITKImagePointer;
121 typedef typename ITKImage::PixelContainer Container;
122 typedef typename itk::ImageRegionConstIterator< ITKImage > ConstIterator;
123 typedef typename itk::ImageRegionIterator< ITKImage > Iterator;
124
127
134
141
149
157
162
163 // ----------------------- Interface --------------------------------------
164 public:
165
176 void updateDomain(const Point & inputDomainShift = Point());
177
183 {
184 return ConstRange(*this);
185 }
186
192 {
193 return Range(*this);
194 }
195
200 const Container & container() const
201 {
202 return *(myITKImagePointer->GetPixelContainer());
203 }
209 {
210 return *(myITKImagePointer->GetPixelContainer());
211 }
212
219 Value operator()(const Point &domainPoint) const;
220
228
235 Value operator()(const Iterator &it) const;
236
243 void setValue(const Point &domainPoint, const Value &aValue);
244
251 void setValue(Iterator &it, const Value &V);
252
253 // ------------------------- methods ------------------------------
254
255
259 const Domain& domain() const
260 {
261 return myDomain;
262 }
263
267 inline
269 {
270 return myITKImagePointer;
271 }
272
273 inline
274 const Point & getDomainShift() const
275 {
276 return myDomainShift;
277 }
278
279 // ------------------------- stream ------------------------------
280
285 void selfDisplay ( std::ostream & out ) const;
286
291 bool isValid() const;
292
293 // ------------------------- Iterators ------------------------------
298 inline
300 {
301 ConstIterator iter = ConstIterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
302 iter.GoToBegin();
303 return iter;
304 }
305
310 inline
312 {
313 Iterator iter = Iterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
314 iter.GoToBegin();
315 return iter;
316 }
317
322 inline
323 const ConstIterator end() const
324 {
325 ConstIterator iter = ConstIterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
326 iter.GoToEnd();
327 return iter;
328 }
329
334 inline
336 {
337 Iterator iter = Iterator(myITKImagePointer, myITKImagePointer->GetLargestPossibleRegion());
338 iter.GoToEnd();
339 return iter;
340 }
341
342 // ----------------- DomainPoint to/from IndexPoint interface -------------
343
352 inline Point getDomainPointFromIndex(const Point &indexPoint) const;
353
354 inline Point getIndexFromDomainPoint(const Point &domainPoint) const;
355
363 inline Point getDomainPointFromItkIndex(const typename ITKImage::IndexType &itkIndexPoint) const;
364
365 inline typename ITKImage::IndexType getItkIndexFromDomainPoint(const Point &domainPoint) const;
366
367
368 // ------------------------- PhysicalPoint interface ----------------------
369
381 inline PhysicalPoint getPhysicalPointFromDomainPoint(const Point &domainPoint) const;
382
383 inline Point getDomainPointFromPhysicalPoint(const PhysicalPoint &physicalPoint) const;
384
406 inline void setImageSpacing(const ImageSpacing& s) const;
407
408 // ------------------------- Private Datas --------------------------------
409 private:
410
411 // ------------------------- Hidden services ------------------------------
412 protected:
413
419
420 // ------------------------- Internals ------------------------------------
421 private:
422
424 Domain myDomain; // cached from myITKImagePointer region. updated when calling update().
425
443 }; // end of class ImageContainerByITKImage
444
451 template <typename T, typename TV>
452 std::ostream&
453 operator<< ( std::ostream & out, const ImageContainerByITKImage<T, TV> & object );
454
455}
457// Includes inline functions.
458#include "DGtal/images/ImageContainerByITKImage.ih"
459
460//
462
463#endif // !defined ImageContainerByITKImage_h
464
465#undef ImageContainerByITKImage_RECURSES
466#endif // else defined(ImageContainerByITKImage_RECURSES)
Aim: model of CConstBidirectionalRangeFromPoint that adapts the domain of an image in order to iterat...
Aim: model of CConstBidirectionalRangeFromPoint and CBidirectionalRangeWithWritableIteratorFromPoint ...
Space::Dimension Dimension
Aim: implements a model of CImageContainer using a ITK Image.
BOOST_CONCEPT_ASSERT((concepts::CDomain< TDomain >))
Point getDomainPointFromIndex(const Point &indexPoint) const
void updateDomain(const Point &inputDomainShift=Point())
Value operator()(const ConstIterator &it) const
ImageContainerByITKImage(const ITKImagePointer &aRef)
DefaultConstImageRange< Self > ConstRange
BOOST_CONCEPT_ASSERT((concepts::CLabel< TValue >))
itk::ImageRegionIterator< ITKImage > Iterator
Value operator()(const Iterator &it) const
PointVector< dimension, double > RealPoint
itk::Image< TValue, dimension > ITKImage
itk::ImageBase< dimension >::SpacingValueType ITKSpacingValueType
static const Domain::Dimension dimension
void selfDisplay(std::ostream &out) const
ImageContainerByITKImage(const Domain &aDomain)
ImageContainerByITKImage & operator=(const ImageContainerByITKImage &other)
ImageContainerByITKImage(const ImageContainerByITKImage &other)
ImageSpacing getImageSpacing() const
Point getDomainPointFromItkIndex(const typename ITKImage::IndexType &itkIndexPoint) const
ITKImage::IndexType getItkIndexFromDomainPoint(const Point &domainPoint) const
void setValue(const Point &domainPoint, const Value &aValue)
ImageContainerByITKImage< TDomain, TValue > Self
Value operator()(const Point &domainPoint) const
void setValue(Iterator &it, const Value &V)
PhysicalPoint getLowerBoundAsPhysicalPoint() const
void setImageSpacing(const ImageSpacing &s) const
Point getIndexFromDomainPoint(const Point &domainPoint) const
PhysicalPoint getPhysicalPointFromDomainPoint(const Point &domainPoint) const
PhysicalPoint getUpperBoundAsPhysicalPoint() const
Point getDomainPointFromPhysicalPoint(const PhysicalPoint &physicalPoint) const
itk::ImageRegionConstIterator< ITKImage > ConstIterator
DGtal is the top-level namespace which contains all DGtal functions and types.
std::ostream & operator<<(std::ostream &out, const ClosedIntegerHalfPlane< TSpace > &object)
Aim: This concept represents a digital domain, i.e. a non mutable subset of points of the given digit...
Definition CDomain.h:130
Aim: Define the concept of DGtal labels. Models of CLabel can be default-constructible,...
Definition CLabel.h:93