Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

itkImageAdaptor.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkImageAdaptor.h,v $
00005   Language:  C++
00006   Date:      $Date: 2008-10-19 12:33:32 $
00007   Version:   $Revision: 1.70 $
00008 
00009   Copyright (c) Insight Software Consortium. All rights reserved.
00010   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
00011 
00012      This software is distributed WITHOUT ANY WARRANTY; without even 
00013      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
00014      PURPOSE.  See the above copyright notices for more information.
00015 
00016 =========================================================================*/
00017 #ifndef __itkImageAdaptor_h
00018 #define __itkImageAdaptor_h
00019 
00020 #include "itkImage.h"
00021 #include "itkDefaultPixelAccessorFunctor.h"
00022 
00023 namespace itk
00024 {
00025 
00047 template <class TImage, class TAccessor >
00048 class ITK_EXPORT ImageAdaptor : public ImageBase< ::itk::GetImageDimension<TImage>::ImageDimension>
00049 {
00050 public:
00051 
00056   itkStaticConstMacro(ImageDimension, unsigned int, TImage::ImageDimension);
00057 
00059   typedef ImageAdaptor                                      Self;
00060   typedef ImageBase<itkGetStaticConstMacro(ImageDimension)> Superclass;
00061   typedef SmartPointer<Self>                                Pointer;
00062   typedef SmartPointer<const Self>                          ConstPointer;
00063   typedef WeakPointer<const Self>                           ConstWeakPointer;
00064 
00066   itkTypeMacro(ImageAdaptor, ImageBase);
00067 
00069   typedef TImage InternalImageType;
00070 
00072   itkNewMacro(Self);  
00073 
00076   typedef typename TAccessor::ExternalType PixelType;
00077 
00080   typedef typename TAccessor::InternalType InternalPixelType;
00081 
00082   typedef PixelType  IOPixelType;
00083 
00086   typedef   TAccessor   AccessorType;
00087 
00090   typedef DefaultPixelAccessorFunctor< Self > AccessorFunctorType;
00091 
00093   typedef typename Superclass::IndexType     IndexType;
00094   typedef typename IndexType::IndexValueType IndexValueType;
00095 
00097   typedef typename Superclass::SizeType    SizeType;
00098   typedef typename SizeType::SizeValueType SizeValueType;
00099 
00101   typedef typename Superclass::OffsetType      OffsetType;
00102   typedef typename OffsetType::OffsetValueType OffsetValueType;
00103 
00106   typedef typename Superclass::RegionType RegionType;
00107 
00110   typedef typename Superclass::SpacingType SpacingType;
00111 
00114   typedef typename Superclass::PointType PointType;
00115 
00119   typedef typename Superclass::DirectionType DirectionType;
00120 
00127   virtual void SetLargestPossibleRegion(const RegionType &region);
00128 
00132   virtual void SetBufferedRegion(const RegionType &region);
00133 
00137   virtual void SetRequestedRegion(const RegionType &region);
00138 
00143   virtual void SetRequestedRegion(DataObject *data);
00144 
00151   virtual const RegionType & GetRequestedRegion() const;
00152 
00161   virtual const RegionType& GetLargestPossibleRegion() const;
00162 
00168   virtual const RegionType& GetBufferedRegion() const;
00169 
00171   inline void Allocate()
00172     {
00173     m_Image->Allocate();
00174     }
00175 
00176 
00179   virtual void Initialize();
00180 
00182   void SetPixel(const IndexType &index, const PixelType & value)
00183     { m_PixelAccessor.Set( m_Image->GetPixel(index), value ); }
00184 
00186   PixelType GetPixel(const IndexType &index) const
00187     { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
00188 
00190   PixelType operator[](const IndexType &index) const
00191     { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
00192 
00194   const OffsetValueType *GetOffsetTable() const;
00195 
00197   IndexType ComputeIndex(OffsetValueType offset) const;
00198 
00201   typedef typename TImage::PixelContainer             PixelContainer;
00202   typedef typename TImage::PixelContainerPointer      PixelContainerPointer;
00203   typedef typename TImage::PixelContainerConstPointer PixelContainerConstPointer;
00204 
00206   PixelContainerPointer GetPixelContainer()
00207     { return m_Image->GetPixelContainer(); }
00208 
00209   const PixelContainer* GetPixelContainer() const
00210     { return m_Image->GetPixelContainer(); }
00211 
00214   void SetPixelContainer( PixelContainer *container );
00215 
00226   virtual void Graft(const DataObject *data);
00227 
00229   typedef InternalPixelType * InternalPixelPointerType;
00230 
00233   InternalPixelType *GetBufferPointer();
00234   const InternalPixelType *GetBufferPointer() const;
00236 
00238   virtual void SetSpacing( const SpacingType &values );
00239   virtual void SetSpacing( const double* values /*[ImageDimension]*/ );
00240   virtual void SetSpacing( const float* values /*[ImageDimension]*/ );
00242 
00246   virtual const SpacingType& GetSpacing() const;
00247 
00251   virtual const PointType& GetOrigin() const;
00252 
00254   virtual void SetOrigin( const PointType values);
00255   virtual void SetOrigin( const double* values /*[ImageDimension]*/ );
00256   virtual void SetOrigin( const float* values /*[ImageDimension]*/ );
00258 
00260   virtual void SetDirection( const DirectionType direction );
00261 
00265   virtual const DirectionType& GetDirection() const;
00266 
00268   virtual void SetImage( TImage * );
00269 
00271   virtual void Modified() const;
00272 
00274   virtual unsigned long GetMTime() const;
00275 
00277   AccessorType & GetPixelAccessor( void ) 
00278     { return m_PixelAccessor; }
00279 
00281   const AccessorType & GetPixelAccessor( void ) const
00282     { return m_PixelAccessor; }
00283 
00285   void SetPixelAccessor( const AccessorType & accessor ) 
00286     { m_PixelAccessor = accessor; }
00287 
00289   virtual void Update();
00290   virtual void CopyInformation(const DataObject *data);
00292 
00295   virtual void UpdateOutputInformation();
00296   virtual void SetRequestedRegionToLargestPossibleRegion();
00297   virtual void PropagateRequestedRegion() throw (InvalidRequestedRegionError);
00298   virtual void UpdateOutputData();
00299   virtual bool VerifyRequestedRegion();
00301 
00306   template<class TCoordRep>
00307   bool TransformPhysicalPointToContinuousIndex(
00308               const Point<TCoordRep,
00309               itkGetStaticConstMacro(ImageDimension)>& point,
00310               ContinuousIndex<TCoordRep,
00311               itkGetStaticConstMacro(ImageDimension)>& index   ) const
00312     {
00313     return m_Image->TransformPhysicalPointToContinuousIndex( point, index );
00314     }
00315 
00320   template<class TCoordRep>
00321   bool TransformPhysicalPointToIndex(
00322             const Point<TCoordRep,
00323             itkGetStaticConstMacro(ImageDimension)>& point,
00324             IndexType & index                                ) const
00325     {
00326     return m_Image->TransformPhysicalPointToIndex( point, index );
00327     }
00328 
00333   template<class TCoordRep>
00334   void TransformContinuousIndexToPhysicalPoint(
00335             const ContinuousIndex<TCoordRep,
00336             itkGetStaticConstMacro(ImageDimension)>& index,
00337             Point<TCoordRep,
00338             itkGetStaticConstMacro(ImageDimension)>& point        ) const
00339     {
00340     m_Image->TransformContinuousIndexToPhysicalPoint( index, point );
00341     }
00342 
00348   template<class TCoordRep>
00349   void TransformIndexToPhysicalPoint(
00350                       const IndexType & index,
00351                       Point<TCoordRep,
00352                       itkGetStaticConstMacro(ImageDimension)>& point ) const
00353     {
00354     m_Image->TransformIndexToPhysicalPoint( index, point );
00355     }
00356 
00357   template<class TCoordRep>
00358   void TransformLocalVectorToPhysicalVector(
00359     const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
00360           FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient ) const
00361     {
00362     m_Image->TransformLocalVectorToPhysicalVector( inputGradient, outputGradient );
00363     }
00364 
00365 protected:
00366 
00367   ImageAdaptor();
00368   virtual ~ImageAdaptor();
00369   void PrintSelf(std::ostream& os, Indent indent) const;
00370   
00371 private:
00372 
00373   ImageAdaptor(const Self&); //purposely not implemented
00374   void operator=(const Self&); //purposely not implemented
00375   
00376   // Adapted image, most of the calls to ImageAdaptor
00377   // will be delegated to this image
00378   typename TImage::Pointer   m_Image;
00379 
00380   // Data accessor object, 
00381   // it converts the presentation of a pixel
00382   AccessorType               m_PixelAccessor;
00383   
00384 
00385 };
00386   
00387 } // end namespace itk
00388 
00389 #ifndef ITK_MANUAL_INSTANTIATION
00390 #include "itkImageAdaptor.txx"
00391 #endif
00392 
00393 #endif
00394 

Generated at Mon Jul 12 2010 18:35:43 for ITK by doxygen 1.7.1 written by Dimitri van Heesch, © 1997-2000