ITK  4.6.0
Insight Segmentation and Registration Toolkit
itkImageAdaptor.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef __itkImageAdaptor_h
19 #define __itkImageAdaptor_h
20 
21 #include "itkImage.h"
22 
23 namespace itk
24 {
25 
26 template <typename TPixelType, unsigned int VImageDimension > class VectorImage;
27 
54 template< typename TImage, typename TAccessor >
55 class ImageAdaptor:public ImageBase< TImage::ImageDimension >
56 {
57 public:
58 
63  itkStaticConstMacro(ImageDimension, unsigned int, TImage::ImageDimension);
64 
66  typedef ImageAdaptor Self;
71 
73  itkTypeMacro(ImageAdaptor, ImageBase);
74 
76  typedef TImage InternalImageType;
77 
79  itkNewMacro(Self);
80 
83  typedef typename TAccessor::ExternalType PixelType;
84 
87  typedef typename TAccessor::InternalType InternalPixelType;
88 
90 
93  typedef TAccessor AccessorType;
94 
97  typedef typename InternalImageType::AccessorFunctorType::template Rebind< Self >::Type AccessorFunctorType;
98 
102 
104  typedef typename Superclass::SizeType SizeType;
106 
110 
114 
118 
122 
127 
128 
134  template <typename UPixelType, unsigned int UImageDimension = TImage::ImageDimension>
135  struct Rebind
136  {
138  };
139 
140 
147  virtual void SetLargestPossibleRegion(const RegionType & region);
148 
152  virtual void SetBufferedRegion(const RegionType & region);
153 
157  virtual void SetRequestedRegion(const RegionType & region);
158 
163  virtual void SetRequestedRegion(const DataObject *data);
164 
171  virtual const RegionType & GetRequestedRegion() const;
172 
181  virtual const RegionType & GetLargestPossibleRegion() const;
182 
188  virtual const RegionType & GetBufferedRegion() const;
189 
191  virtual void Allocate(bool initialize = false) ITK_OVERRIDE;
192 
195  virtual void Initialize();
196 
198  void SetPixel(const IndexType & index, const PixelType & value)
199  { m_PixelAccessor.Set(m_Image->GetPixel(index), value); }
200 
202  PixelType GetPixel(const IndexType & index) const
203  { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
204 
206  PixelType operator[](const IndexType & index) const
207  { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
208 
210  const OffsetValueType * GetOffsetTable() const;
211 
213  IndexType ComputeIndex(OffsetValueType offset) const;
214 
217  typedef typename TImage::PixelContainer PixelContainer;
218  typedef typename TImage::PixelContainerPointer PixelContainerPointer;
219  typedef typename TImage::PixelContainerConstPointer PixelContainerConstPointer;
220 
223  { return m_Image->GetPixelContainer(); }
224 
226  { return m_Image->GetPixelContainer(); }
227 
230  void SetPixelContainer(PixelContainer *container);
231 
242  virtual void Graft(const DataObject *data);
243 
246 
250 
251  const InternalPixelType * GetBufferPointer() const;
252 
254  virtual void SetSpacing(const SpacingType & values);
255 
256  virtual void SetSpacing(const double *values /*[ImageDimension]*/);
257 
258  virtual void SetSpacing(const float *values /*[ImageDimension]*/);
259 
263  virtual const SpacingType & GetSpacing() const;
264 
268  virtual const PointType & GetOrigin() const;
269 
271  virtual void SetOrigin(const PointType values);
272 
273  virtual void SetOrigin(const double *values /*[ImageDimension]*/);
274 
275  virtual void SetOrigin(const float *values /*[ImageDimension]*/);
276 
278  virtual void SetDirection(const DirectionType & direction);
279 
283  virtual const DirectionType & GetDirection() const;
284 
286  virtual void SetImage(TImage *);
287 
289  virtual void Modified() const;
290 
292  virtual ModifiedTimeType GetMTime() const;
293 
296  { return m_PixelAccessor; }
297 
299  const AccessorType & GetPixelAccessor(void) const
300  { return m_PixelAccessor; }
301 
303  void SetPixelAccessor(const AccessorType & accessor)
304  { m_PixelAccessor = accessor; }
305 
307  virtual void Update();
308 
309  virtual void CopyInformation(const DataObject *data);
310 
313  virtual void UpdateOutputInformation();
314 
316 
317  virtual void PropagateRequestedRegion();
318 
319  virtual void UpdateOutputData();
320 
321  virtual bool VerifyRequestedRegion();
322 
327  template< typename TCoordRep >
329  const Point< TCoordRep,
330  itkGetStaticConstMacro(ImageDimension) > & point,
331  ContinuousIndex< TCoordRep,
332  itkGetStaticConstMacro(ImageDimension) > & index) const
333  {
334  return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
335  }
336 
341  template< typename TCoordRep >
343  const Point< TCoordRep,
344  itkGetStaticConstMacro(ImageDimension) > & point,
345  IndexType & index) const
346  {
347  return m_Image->TransformPhysicalPointToIndex(point, index);
348  }
349 
354  template< typename TCoordRep >
356  const ContinuousIndex< TCoordRep,
357  itkGetStaticConstMacro(ImageDimension) > & index,
358  Point< TCoordRep,
359  itkGetStaticConstMacro(ImageDimension) > & point) const
360  {
361  m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
362  }
363 
369  template< typename TCoordRep >
371  const IndexType & index,
372  Point< TCoordRep,
373  itkGetStaticConstMacro(ImageDimension) > & point) const
374  {
375  m_Image->TransformIndexToPhysicalPoint(index, point);
376  }
377 
378  template< typename TCoordRep >
380  const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
381  FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
382  {
383  m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
384  }
385 
386  template< typename TCoordRep >
388  const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
389  FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
390  {
391  m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
392  }
393 
394 protected:
395 
396  ImageAdaptor();
397  virtual ~ImageAdaptor();
398  void PrintSelf(std::ostream & os, Indent indent) const;
399 
400 private:
401 
402  ImageAdaptor(const Self &); //purposely not implemented
403  void operator=(const Self &); //purposely not implemented
404 
405  // a specialized method to update PixelAccessors for VectorImages,
406  // to have the correct vector length of the image.
407  template< typename TPixelType >
408  void UpdateAccessor( typename ::itk::VectorImage< TPixelType, ImageDimension > * itkNotUsed( dummy ) )
409  {
410  this->m_PixelAccessor.SetVectorLength( this->m_Image->GetNumberOfComponentsPerPixel() );
411  }
412 
413  // The other image types don't expect an accessor which needs any updates
414  template< typename T > void UpdateAccessor( T *itkNotUsed( dummy ) ) { }
415 
416  // Adapted image, most of the calls to ImageAdaptor
417  // will be delegated to this image
418  typename TImage::Pointer m_Image;
419 
420  // Data accessor object,
421  // it converts the presentation of a pixel
423 };
424 } // end namespace itk
425 
426 #ifndef ITK_MANUAL_INSTANTIATION
427 #include "itkImageAdaptor.hxx"
428 #endif
429 
430 #endif
SizeType::SizeValueType SizeValueType
void UpdateAccessor(T *)
Superclass::RegionType RegionType
virtual void SetDirection(const DirectionType &direction)
Index< VImageDimension > IndexType
Definition: itkImageBase.h:134
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &point) const
itk::SizeValueType SizeValueType
Definition: itkSize.h:60
const PixelContainer * GetPixelContainer() const
InternalPixelType * GetBufferPointer()
Superclass::SpacingType SpacingType
TImage::PixelContainer PixelContainer
virtual void UpdateOutputData()
signed long OffsetValueType
Definition: itkIntTypes.h:154
virtual void Allocate(bool initialize=false) ITK_OVERRIDE
virtual const RegionType & GetBufferedRegion() const
IndexType::IndexValueType IndexValueType
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &index, Point< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &point) const
virtual const SpacingType & GetSpacing() const
ImageBase< itkGetStaticConstMacro(ImageDimension) > Superclass
unsigned long ModifiedTimeType
Definition: itkIntTypes.h:164
InternalPixelType * InternalPixelPointerType
Image< UPixelType, UImageDimension > Type
AccessorType m_PixelAccessor
virtual void SetLargestPossibleRegion(const RegionType &region)
ImageAdaptor Self
Implements a weak reference to an object.
void UpdateAccessor(typename::itk::VectorImage< TPixelType, ImageDimension > *)
PixelType GetPixel(const IndexType &index) const
Size< VImageDimension > SizeType
Definition: itkImageBase.h:143
virtual bool VerifyRequestedRegion()
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:159
virtual void UpdateOutputInformation()
SmartPointer< Self > Pointer
virtual const PointType & GetOrigin() const
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Definition: itkImageBase.h:164
static const unsigned int ImageDimension
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
InternalImageType::AccessorFunctorType::template Rebind< Self >::Type AccessorFunctorType
const AccessorType & GetPixelAccessor(void) const
virtual void SetImage(TImage *)
PixelContainerPointer GetPixelContainer()
virtual void SetRequestedRegionToLargestPossibleRegion()
WeakPointer< const Self > ConstWeakPointer
virtual void Update()
virtual void SetBufferedRegion(const RegionType &region)
void SetPixelAccessor(const AccessorType &accessor)
TImage::PixelContainerConstPointer PixelContainerConstPointer
::itk::IndexValueType IndexValueType
Definition: itkIndex.h:79
TAccessor::InternalType InternalPixelType
TAccessor::ExternalType PixelType
virtual void SetSpacing(const SpacingType &values)
bool TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &point, ContinuousIndex< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &index) const
Get the continuous index from a physical point.
AccessorType & GetPixelAccessor(void)
Superclass::DirectionType DirectionType
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:139
virtual void SetOrigin(const PointType values)
Superclass::SizeType SizeType
Vector< SpacingValueType, VImageDimension > SpacingType
Definition: itkImageBase.h:154
void operator=(const Self &)
virtual void Initialize()
IndexType ComputeIndex(OffsetValueType offset) const
virtual const RegionType & GetLargestPossibleRegion() const
virtual void Modified() const
void SetPixel(const IndexType &index, const PixelType &value)
PixelType operator[](const IndexType &index) const
virtual void CopyInformation(const DataObject *data)
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &inputGradient, FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &outputGradient) const
virtual void Graft(const DataObject *data)
Base class for templated image classes.
Definition: itkImageBase.h:112
A templated class holding a point in n-Dimensional image space.
virtual ~ImageAdaptor()
virtual const DirectionType & GetDirection() const
void SetPixelContainer(PixelContainer *container)
Control indentation during Print() invocation.
Definition: itkIndent.h:49
OffsetType::OffsetValueType OffsetValueType
bool TransformPhysicalPointToIndex(const Point< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &point, IndexType &index) const
Superclass::IndexType IndexType
virtual const RegionType & GetRequestedRegion() const
Superclass::PointType PointType
const OffsetValueType * GetOffsetTable() const
TImage::PixelContainerPointer PixelContainerPointer
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:147
virtual void PropagateRequestedRegion()
Give access to partial aspects of voxels from an Image.
itk::OffsetValueType OffsetValueType
Definition: itkOffset.h:69
Superclass::OffsetType OffsetType
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &inputGradient, FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > &outputGradient) const
SmartPointer< const Self > ConstPointer
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
TAccessor AccessorType
virtual void SetRequestedRegion(const RegionType &region)
virtual ModifiedTimeType GetMTime() const
TImage::Pointer m_Image
void PrintSelf(std::ostream &os, Indent indent) const