ITK  4.2.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< class TImage, class TAccessor >
55 class ITK_EXPORT ImageAdaptor:public ImageBase< ::itk::GetImageDimension< 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 
100  typedef typename Superclass::IndexType IndexType;
102 
104  typedef typename Superclass::SizeType SizeType;
106 
108  typedef typename Superclass::OffsetType OffsetType;
110 
113  typedef typename Superclass::RegionType RegionType;
114 
117  typedef typename Superclass::SpacingType SpacingType;
118 
121  typedef typename Superclass::PointType PointType;
122 
126  typedef typename Superclass::DirectionType DirectionType;
127 
128 
134  template <class UPixelType, unsigned int UImageDimension = ::itk::GetImageDimension< 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  inline void Allocate()
192  {
193  m_Image->Allocate();
194  }
195 
198  virtual void Initialize();
199 
201  void SetPixel(const IndexType & index, const PixelType & value)
202  { m_PixelAccessor.Set(m_Image->GetPixel(index), value); }
203 
205  PixelType GetPixel(const IndexType & index) const
206  { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
207 
209  PixelType operator[](const IndexType & index) const
210  { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
211 
213  const OffsetValueType * GetOffsetTable() const;
214 
216  IndexType ComputeIndex(OffsetValueType offset) const;
217 
220  typedef typename TImage::PixelContainer PixelContainer;
221  typedef typename TImage::PixelContainerPointer PixelContainerPointer;
222  typedef typename TImage::PixelContainerConstPointer PixelContainerConstPointer;
223 
225  PixelContainerPointer GetPixelContainer()
226  { return m_Image->GetPixelContainer(); }
227 
228  const PixelContainer * GetPixelContainer() const
229  { return m_Image->GetPixelContainer(); }
230 
233  void SetPixelContainer(PixelContainer *container);
234 
245  virtual void Graft(const DataObject *data);
246 
249 
252  InternalPixelType * GetBufferPointer();
253 
254  const InternalPixelType * GetBufferPointer() const;
255 
257  virtual void SetSpacing(const SpacingType & values);
258 
259  virtual void SetSpacing(const double *values /*[ImageDimension]*/);
260 
261  virtual void SetSpacing(const float *values /*[ImageDimension]*/);
262 
266  virtual const SpacingType & GetSpacing() const;
267 
271  virtual const PointType & GetOrigin() const;
272 
274  virtual void SetOrigin(const PointType values);
275 
276  virtual void SetOrigin(const double *values /*[ImageDimension]*/);
277 
278  virtual void SetOrigin(const float *values /*[ImageDimension]*/);
279 
281  virtual void SetDirection(const DirectionType direction);
282 
286  virtual const DirectionType & GetDirection() const;
287 
289  virtual void SetImage(TImage *);
290 
292  virtual void Modified() const;
293 
295  virtual unsigned long GetMTime() const;
296 
298  AccessorType & GetPixelAccessor(void)
299  { return m_PixelAccessor; }
300 
302  const AccessorType & GetPixelAccessor(void) const
303  { return m_PixelAccessor; }
304 
306  void SetPixelAccessor(const AccessorType & accessor)
307  { m_PixelAccessor = accessor; }
308 
310  virtual void Update();
311 
312  virtual void CopyInformation(const DataObject *data);
313 
316  virtual void UpdateOutputInformation();
317 
318  virtual void SetRequestedRegionToLargestPossibleRegion();
319 
320  virtual void PropagateRequestedRegion()
321  throw ( InvalidRequestedRegionError );
322 
323  virtual void UpdateOutputData();
324 
325  virtual bool VerifyRequestedRegion();
326 
331  template< class TCoordRep >
332  bool TransformPhysicalPointToContinuousIndex(
333  const Point< TCoordRep,
334  itkGetStaticConstMacro(ImageDimension) > & point,
335  ContinuousIndex< TCoordRep,
336  itkGetStaticConstMacro(ImageDimension) > & index) const
337  {
338  return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
339  }
340 
345  template< class TCoordRep >
346  bool TransformPhysicalPointToIndex(
347  const Point< TCoordRep,
348  itkGetStaticConstMacro(ImageDimension) > & point,
349  IndexType & index) const
350  {
351  return m_Image->TransformPhysicalPointToIndex(point, index);
352  }
353 
358  template< class TCoordRep >
359  void TransformContinuousIndexToPhysicalPoint(
360  const ContinuousIndex< TCoordRep,
361  itkGetStaticConstMacro(ImageDimension) > & index,
362  Point< TCoordRep,
363  itkGetStaticConstMacro(ImageDimension) > & point) const
364  {
365  m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
366  }
367 
373  template< class TCoordRep >
374  void TransformIndexToPhysicalPoint(
375  const IndexType & index,
376  Point< TCoordRep,
377  itkGetStaticConstMacro(ImageDimension) > & point) const
378  {
379  m_Image->TransformIndexToPhysicalPoint(index, point);
380  }
381 
382  template< class TCoordRep >
383  void TransformLocalVectorToPhysicalVector(
384  const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
385  FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
386  {
387  m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
388  }
389 
390  template< class TCoordRep >
391  void TransformPhysicalVectorToLocalVector(
392  const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
393  FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
394  {
395  m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
396  }
397 
398 protected:
399 
400  ImageAdaptor();
401  virtual ~ImageAdaptor();
402  void PrintSelf(std::ostream & os, Indent indent) const;
403 
404 private:
405 
406  ImageAdaptor(const Self &); //purposely not implemented
407  void operator=(const Self &); //purposely not implemented
408 
409  // a specialized method to update PixelAccessors for VectorImages,
410  // to have the correct vector length of the image.
411  template< class TPixelType >
412  void UpdateAccessor( typename itk::VectorImage< TPixelType, ImageDimension > * itkNotUsed( dummy ) )
413  {
414  this->m_PixelAccessor.SetVectorLength( this->m_Image->GetNumberOfComponentsPerPixel() );
415  }
416 
417  // The other image types don't expect an accessor which needs any updates
418  template< class T > void UpdateAccessor( T *itkNotUsed( dummy ) ) { }
419 
420  // Adapted image, most of the calls to ImageAdaptor
421  // will be delegated to this image
422  typename TImage::Pointer m_Image;
423 
424  // Data accessor object,
425  // it converts the presentation of a pixel
427 };
428 } // end namespace itk
429 
430 #ifndef ITK_MANUAL_INSTANTIATION
431 #include "itkImageAdaptor.hxx"
432 #endif
433 
434 #endif
435