ITK  4.0.0
Insight Segmentation and Registration Toolkit
itkImageAdaptor.h
Go to the documentation of this file.
00001 /*=========================================================================
00002  *
00003  *  Copyright Insight Software Consortium
00004  *
00005  *  Licensed under the Apache License, Version 2.0 (the "License");
00006  *  you may not use this file except in compliance with the License.
00007  *  You may obtain a copy of the License at
00008  *
00009  *         http://www.apache.org/licenses/LICENSE-2.0.txt
00010  *
00011  *  Unless required by applicable law or agreed to in writing, software
00012  *  distributed under the License is distributed on an "AS IS" BASIS,
00013  *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
00014  *  See the License for the specific language governing permissions and
00015  *  limitations under the License.
00016  *
00017  *=========================================================================*/
00018 #ifndef __itkImageAdaptor_h
00019 #define __itkImageAdaptor_h
00020 
00021 #include "itkImage.h"
00022 
00023 namespace itk
00024 {
00051 template< class TImage, class TAccessor >
00052 class ITK_EXPORT ImageAdaptor:public ImageBase< ::itk::GetImageDimension< TImage >::ImageDimension >
00053 {
00054 public:
00055 
00060   itkStaticConstMacro(ImageDimension, unsigned int, TImage::ImageDimension);
00061 
00063   typedef ImageAdaptor                                        Self;
00064   typedef ImageBase< itkGetStaticConstMacro(ImageDimension) > Superclass;
00065   typedef SmartPointer< Self >                                Pointer;
00066   typedef SmartPointer< const Self >                          ConstPointer;
00067   typedef WeakPointer< const Self >                           ConstWeakPointer;
00068 
00070   itkTypeMacro(ImageAdaptor, ImageBase);
00071 
00073   typedef TImage InternalImageType;
00074 
00076   itkNewMacro(Self);
00077 
00080   typedef typename TAccessor::ExternalType PixelType;
00081 
00084   typedef typename TAccessor::InternalType InternalPixelType;
00085 
00086   typedef PixelType IOPixelType;
00087 
00090   typedef   TAccessor AccessorType;
00091 
00094   typedef DefaultPixelAccessorFunctor< Self > AccessorFunctorType;
00095 
00097   typedef typename Superclass::IndexType     IndexType;
00098   typedef typename IndexType::IndexValueType IndexValueType;
00099 
00101   typedef typename Superclass::SizeType    SizeType;
00102   typedef typename SizeType::SizeValueType SizeValueType;
00103 
00105   typedef typename Superclass::OffsetType      OffsetType;
00106   typedef typename OffsetType::OffsetValueType OffsetValueType;
00107 
00110   typedef typename Superclass::RegionType RegionType;
00111 
00114   typedef typename Superclass::SpacingType SpacingType;
00115 
00118   typedef typename Superclass::PointType PointType;
00119 
00123   typedef typename Superclass::DirectionType DirectionType;
00124 
00131   virtual void SetLargestPossibleRegion(const RegionType & region);
00132 
00136   virtual void SetBufferedRegion(const RegionType & region);
00137 
00141   virtual void SetRequestedRegion(const RegionType & region);
00142 
00147   virtual void SetRequestedRegion(const DataObject *data);
00148 
00155   virtual const RegionType & GetRequestedRegion() const;
00156 
00165   virtual const RegionType & GetLargestPossibleRegion() const;
00166 
00172   virtual const RegionType & GetBufferedRegion() const;
00173 
00175   inline void Allocate()
00176   {
00177     m_Image->Allocate();
00178   }
00179 
00182   virtual void Initialize();
00183 
00185   void SetPixel(const IndexType & index, const PixelType & value)
00186   { m_PixelAccessor.Set(m_Image->GetPixel(index), value); }
00187 
00189   PixelType GetPixel(const IndexType & index) const
00190   { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
00191 
00193   PixelType operator[](const IndexType & index) const
00194   { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
00195 
00197   const OffsetValueType * GetOffsetTable() const;
00198 
00200   IndexType ComputeIndex(OffsetValueType offset) const;
00201 
00204   typedef typename TImage::PixelContainer             PixelContainer;
00205   typedef typename TImage::PixelContainerPointer      PixelContainerPointer;
00206   typedef typename TImage::PixelContainerConstPointer PixelContainerConstPointer;
00207 
00209   PixelContainerPointer GetPixelContainer()
00210   { return m_Image->GetPixelContainer(); }
00211 
00212   const PixelContainer * GetPixelContainer() const
00213   { return m_Image->GetPixelContainer(); }
00214 
00217   void SetPixelContainer(PixelContainer *container);
00218 
00229   virtual void Graft(const DataObject *data);
00230 
00232   typedef InternalPixelType *InternalPixelPointerType;
00233 
00236   InternalPixelType * GetBufferPointer();
00237 
00238   const InternalPixelType * GetBufferPointer() const;
00239 
00241   virtual void SetSpacing(const SpacingType & values);
00242 
00243   virtual void SetSpacing(const double *values /*[ImageDimension]*/);
00244 
00245   virtual void SetSpacing(const float *values /*[ImageDimension]*/);
00246 
00250   virtual const SpacingType & GetSpacing() const;
00251 
00255   virtual const PointType & GetOrigin() const;
00256 
00258   virtual void SetOrigin(const PointType values);
00259 
00260   virtual void SetOrigin(const double *values /*[ImageDimension]*/);
00261 
00262   virtual void SetOrigin(const float *values /*[ImageDimension]*/);
00263 
00265   virtual void SetDirection(const DirectionType direction);
00266 
00270   virtual const DirectionType & GetDirection() const;
00271 
00273   virtual void SetImage(TImage *);
00274 
00276   virtual void Modified() const;
00277 
00279   virtual unsigned long GetMTime() const;
00280 
00282   AccessorType & GetPixelAccessor(void)
00283   { return m_PixelAccessor; }
00284 
00286   const AccessorType & GetPixelAccessor(void) const
00287   { return m_PixelAccessor; }
00288 
00290   void SetPixelAccessor(const AccessorType & accessor)
00291   { m_PixelAccessor = accessor; }
00292 
00294   virtual void Update();
00295 
00296   virtual void CopyInformation(const DataObject *data);
00297 
00300   virtual void UpdateOutputInformation();
00301 
00302   virtual void SetRequestedRegionToLargestPossibleRegion();
00303 
00304   virtual void PropagateRequestedRegion()
00305   throw ( InvalidRequestedRegionError );
00306 
00307   virtual void UpdateOutputData();
00308 
00309   virtual bool VerifyRequestedRegion();
00310 
00315   template< class TCoordRep >
00316   bool TransformPhysicalPointToContinuousIndex(
00317     const Point< TCoordRep,
00318                  itkGetStaticConstMacro(ImageDimension) > & point,
00319     ContinuousIndex< TCoordRep,
00320                      itkGetStaticConstMacro(ImageDimension) > & index) const
00321   {
00322     return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
00323   }
00324 
00329   template< class TCoordRep >
00330   bool TransformPhysicalPointToIndex(
00331     const Point< TCoordRep,
00332                  itkGetStaticConstMacro(ImageDimension) > & point,
00333     IndexType & index) const
00334   {
00335     return m_Image->TransformPhysicalPointToIndex(point, index);
00336   }
00337 
00342   template< class TCoordRep >
00343   void TransformContinuousIndexToPhysicalPoint(
00344     const ContinuousIndex< TCoordRep,
00345                            itkGetStaticConstMacro(ImageDimension) > & index,
00346     Point< TCoordRep,
00347            itkGetStaticConstMacro(ImageDimension) > & point) const
00348   {
00349     m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
00350   }
00351 
00357   template< class TCoordRep >
00358   void TransformIndexToPhysicalPoint(
00359     const IndexType & index,
00360     Point< TCoordRep,
00361            itkGetStaticConstMacro(ImageDimension) > & point) const
00362   {
00363     m_Image->TransformIndexToPhysicalPoint(index, point);
00364   }
00365 
00366   template< class TCoordRep >
00367   void TransformLocalVectorToPhysicalVector(
00368     const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
00369     FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
00370   {
00371     m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
00372   }
00373 
00374   template< class TCoordRep >
00375   void TransformPhysicalVectorToLocalVector(
00376     const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
00377     FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
00378   {
00379     m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
00380   }
00381 
00382 protected:
00383 
00384   ImageAdaptor();
00385   virtual ~ImageAdaptor();
00386   void PrintSelf(std::ostream & os, Indent indent) const;
00387 
00388 private:
00389 
00390   ImageAdaptor(const Self &);   //purposely not implemented
00391   void operator=(const Self &); //purposely not implemented
00392 
00393   // Adapted image, most of the calls to ImageAdaptor
00394   // will be delegated to this image
00395   typename TImage::Pointer m_Image;
00396 
00397   // Data accessor object,
00398   // it converts the presentation of a pixel
00399   AccessorType m_PixelAccessor;
00400 };
00401 } // end namespace itk
00402 
00403 #ifndef ITK_MANUAL_INSTANTIATION
00404 #include "itkImageAdaptor.hxx"
00405 #endif
00406 
00407 #endif
00408