ITK  4.1.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 
00125 
00131   template <class UPixelType, unsigned int UImageDimension = ::itk::GetImageDimension< TImage >::ImageDimension>
00132   struct Rebind
00133     {
00134       typedef itk::Image<UPixelType, UImageDimension>  Type;
00135     };
00136 
00137 
00144   virtual void SetLargestPossibleRegion(const RegionType & region);
00145 
00149   virtual void SetBufferedRegion(const RegionType & region);
00150 
00154   virtual void SetRequestedRegion(const RegionType & region);
00155 
00160   virtual void SetRequestedRegion(const DataObject *data);
00161 
00168   virtual const RegionType & GetRequestedRegion() const;
00169 
00178   virtual const RegionType & GetLargestPossibleRegion() const;
00179 
00185   virtual const RegionType & GetBufferedRegion() const;
00186 
00188   inline void Allocate()
00189   {
00190     m_Image->Allocate();
00191   }
00192 
00195   virtual void Initialize();
00196 
00198   void SetPixel(const IndexType & index, const PixelType & value)
00199   { m_PixelAccessor.Set(m_Image->GetPixel(index), value); }
00200 
00202   PixelType GetPixel(const IndexType & index) const
00203   { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
00204 
00206   PixelType operator[](const IndexType & index) const
00207   { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
00208 
00210   const OffsetValueType * GetOffsetTable() const;
00211 
00213   IndexType ComputeIndex(OffsetValueType offset) const;
00214 
00217   typedef typename TImage::PixelContainer             PixelContainer;
00218   typedef typename TImage::PixelContainerPointer      PixelContainerPointer;
00219   typedef typename TImage::PixelContainerConstPointer PixelContainerConstPointer;
00220 
00222   PixelContainerPointer GetPixelContainer()
00223   { return m_Image->GetPixelContainer(); }
00224 
00225   const PixelContainer * GetPixelContainer() const
00226   { return m_Image->GetPixelContainer(); }
00227 
00230   void SetPixelContainer(PixelContainer *container);
00231 
00242   virtual void Graft(const DataObject *data);
00243 
00245   typedef InternalPixelType *InternalPixelPointerType;
00246 
00249   InternalPixelType * GetBufferPointer();
00250 
00251   const InternalPixelType * GetBufferPointer() const;
00252 
00254   virtual void SetSpacing(const SpacingType & values);
00255 
00256   virtual void SetSpacing(const double *values /*[ImageDimension]*/);
00257 
00258   virtual void SetSpacing(const float *values /*[ImageDimension]*/);
00259 
00263   virtual const SpacingType & GetSpacing() const;
00264 
00268   virtual const PointType & GetOrigin() const;
00269 
00271   virtual void SetOrigin(const PointType values);
00272 
00273   virtual void SetOrigin(const double *values /*[ImageDimension]*/);
00274 
00275   virtual void SetOrigin(const float *values /*[ImageDimension]*/);
00276 
00278   virtual void SetDirection(const DirectionType direction);
00279 
00283   virtual const DirectionType & GetDirection() const;
00284 
00286   virtual void SetImage(TImage *);
00287 
00289   virtual void Modified() const;
00290 
00292   virtual unsigned long GetMTime() const;
00293 
00295   AccessorType & GetPixelAccessor(void)
00296   { return m_PixelAccessor; }
00297 
00299   const AccessorType & GetPixelAccessor(void) const
00300   { return m_PixelAccessor; }
00301 
00303   void SetPixelAccessor(const AccessorType & accessor)
00304   { m_PixelAccessor = accessor; }
00305 
00307   virtual void Update();
00308 
00309   virtual void CopyInformation(const DataObject *data);
00310 
00313   virtual void UpdateOutputInformation();
00314 
00315   virtual void SetRequestedRegionToLargestPossibleRegion();
00316 
00317   virtual void PropagateRequestedRegion()
00318   throw ( InvalidRequestedRegionError );
00319 
00320   virtual void UpdateOutputData();
00321 
00322   virtual bool VerifyRequestedRegion();
00323 
00328   template< class TCoordRep >
00329   bool TransformPhysicalPointToContinuousIndex(
00330     const Point< TCoordRep,
00331                  itkGetStaticConstMacro(ImageDimension) > & point,
00332     ContinuousIndex< TCoordRep,
00333                      itkGetStaticConstMacro(ImageDimension) > & index) const
00334   {
00335     return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
00336   }
00337 
00342   template< class TCoordRep >
00343   bool TransformPhysicalPointToIndex(
00344     const Point< TCoordRep,
00345                  itkGetStaticConstMacro(ImageDimension) > & point,
00346     IndexType & index) const
00347   {
00348     return m_Image->TransformPhysicalPointToIndex(point, index);
00349   }
00350 
00355   template< class TCoordRep >
00356   void TransformContinuousIndexToPhysicalPoint(
00357     const ContinuousIndex< TCoordRep,
00358                            itkGetStaticConstMacro(ImageDimension) > & index,
00359     Point< TCoordRep,
00360            itkGetStaticConstMacro(ImageDimension) > & point) const
00361   {
00362     m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
00363   }
00364 
00370   template< class TCoordRep >
00371   void TransformIndexToPhysicalPoint(
00372     const IndexType & index,
00373     Point< TCoordRep,
00374            itkGetStaticConstMacro(ImageDimension) > & point) const
00375   {
00376     m_Image->TransformIndexToPhysicalPoint(index, point);
00377   }
00378 
00379   template< class TCoordRep >
00380   void TransformLocalVectorToPhysicalVector(
00381     const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
00382     FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
00383   {
00384     m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
00385   }
00386 
00387   template< class TCoordRep >
00388   void TransformPhysicalVectorToLocalVector(
00389     const FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & inputGradient,
00390     FixedArray< TCoordRep, itkGetStaticConstMacro(ImageDimension) > & outputGradient) const
00391   {
00392     m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
00393   }
00394 
00395 protected:
00396 
00397   ImageAdaptor();
00398   virtual ~ImageAdaptor();
00399   void PrintSelf(std::ostream & os, Indent indent) const;
00400 
00401 private:
00402 
00403   ImageAdaptor(const Self &);   //purposely not implemented
00404   void operator=(const Self &); //purposely not implemented
00405 
00406   // Adapted image, most of the calls to ImageAdaptor
00407   // will be delegated to this image
00408   typename TImage::Pointer m_Image;
00409 
00410   // Data accessor object,
00411   // it converts the presentation of a pixel
00412   AccessorType m_PixelAccessor;
00413 };
00414 } // end namespace itk
00415 
00416 #ifndef ITK_MANUAL_INSTANTIATION
00417 #include "itkImageAdaptor.hxx"
00418 #endif
00419 
00420 #endif
00421