ITK  4.1.0
Insight Segmentation and Registration Toolkit
itkImage.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 __itkImage_h
00019 #define __itkImage_h
00020 
00021 #include "itkImageRegion.h"
00022 #include "itkImportImageContainer.h"
00023 #include "itkDefaultPixelAccessor.h"
00024 #include "itkDefaultPixelAccessorFunctor.h"
00025 #include "itkPoint.h"
00026 #include "itkFixedArray.h"
00027 #include "itkWeakPointer.h"
00028 #include "itkNeighborhoodAccessorFunctor.h"
00029 
00030 namespace itk
00031 {
00085 template< class TPixel, unsigned int VImageDimension = 2 >
00086 class ITK_EXPORT Image:public ImageBase< VImageDimension >
00087 {
00088 public:
00090   typedef Image                        Self;
00091   typedef ImageBase< VImageDimension > Superclass;
00092   typedef SmartPointer< Self >         Pointer;
00093   typedef SmartPointer< const Self >   ConstPointer;
00094   typedef WeakPointer< const Self >    ConstWeakPointer;
00095 
00097   itkNewMacro(Self);
00098 
00100   itkTypeMacro(Image, ImageBase);
00101 
00104   typedef TPixel PixelType;
00105 
00107   typedef TPixel ValueType;
00108 
00113   typedef TPixel InternalPixelType;
00114 
00115   typedef PixelType IOPixelType;
00116 
00119   typedef DefaultPixelAccessor< PixelType >   AccessorType;
00120   typedef DefaultPixelAccessorFunctor< Self > AccessorFunctorType;
00121 
00124   typedef NeighborhoodAccessorFunctor< Self > NeighborhoodAccessorFunctorType;
00125 
00130   itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
00131 
00133   typedef typename Superclass::IndexType      IndexType;
00134   typedef typename Superclass::IndexValueType IndexValueType;
00135 
00137   typedef typename Superclass::OffsetType OffsetType;
00138 
00140   typedef typename Superclass::SizeType      SizeType;
00141   typedef typename Superclass::SizeValueType SizeValueType;
00142 
00144   typedef ImportImageContainer< SizeValueType, PixelType > PixelContainer;
00145 
00147   typedef typename Superclass::DirectionType DirectionType;
00148 
00151   typedef typename Superclass::RegionType RegionType;
00152 
00155   typedef typename Superclass::SpacingType      SpacingType;
00156   typedef typename Superclass::SpacingValueType SpacingValueType;
00157 
00160   typedef typename Superclass::PointType PointType;
00161 
00163   typedef typename PixelContainer::Pointer      PixelContainerPointer;
00164   typedef typename PixelContainer::ConstPointer PixelContainerConstPointer;
00165 
00167   typedef typename Superclass::OffsetValueType OffsetValueType;
00168 
00174   template <class UPixelType, unsigned int UImageDimension = VImageDimension>
00175   struct Rebind
00176     {
00177       typedef itk::Image<UPixelType, UImageDimension>  Type;
00178     };
00179 
00180 
00183   void Allocate();
00184 
00188   void SetRegions(RegionType region)
00189   {
00190     this->SetLargestPossibleRegion(region);
00191     this->SetBufferedRegion(region);
00192     this->SetRequestedRegion(region);
00193   }
00195 
00196   void SetRegions(SizeType size)
00197   {
00198     RegionType region; region.SetSize(size);
00199 
00200     this->SetLargestPossibleRegion(region);
00201     this->SetBufferedRegion(region);
00202     this->SetRequestedRegion(region);
00203   }
00204 
00207   virtual void Initialize();
00208 
00211   void FillBuffer(const TPixel & value);
00212 
00218   void SetPixel(const IndexType & index, const TPixel & value)
00219   {
00220     OffsetValueType offset = this->ComputeOffset(index);
00221     ( *m_Buffer )[offset] = value;
00222   }
00223 
00228   const TPixel & GetPixel(const IndexType & index) const
00229   {
00230     OffsetValueType offset = this->ComputeOffset(index);
00231     return ( ( *m_Buffer )[offset] );
00232   }
00233 
00238   TPixel & GetPixel(const IndexType & index)
00239   {
00240     OffsetValueType offset = this->ComputeOffset(index);
00241     return ( ( *m_Buffer )[offset] );
00242   }
00243 
00248   TPixel & operator[](const IndexType & index)
00249   { return this->GetPixel(index); }
00250 
00255   const TPixel & operator[](const IndexType & index) const
00256   { return this->GetPixel(index); }
00257 
00260   virtual TPixel * GetBufferPointer()
00261   { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00262   virtual const TPixel * GetBufferPointer() const
00263   { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00265 
00267   PixelContainer * GetPixelContainer()
00268   { return m_Buffer.GetPointer(); }
00269 
00270   const PixelContainer * GetPixelContainer() const
00271   { return m_Buffer.GetPointer(); }
00272 
00275   void SetPixelContainer(PixelContainer *container);
00276 
00287   virtual void Graft(const DataObject *data);
00288 
00290   AccessorType GetPixelAccessor(void)
00291   { return AccessorType(); }
00292 
00294   const AccessorType GetPixelAccessor(void) const
00295   { return AccessorType(); }
00296 
00298   NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
00299   { return NeighborhoodAccessorFunctorType(); }
00300 
00302   const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
00303   { return NeighborhoodAccessorFunctorType(); }
00304 
00305   virtual unsigned int GetNumberOfComponentsPerPixel() const;
00306 
00307 protected:
00308   Image();
00309   void PrintSelf(std::ostream & os, Indent indent) const;
00310 
00311   virtual ~Image() {}
00312 
00318   virtual void ComputeIndexToPhysicalPointMatrices();
00319 
00320 private:
00321   Image(const Self &);          //purposely not implemented
00322   void operator=(const Self &); //purposely not implemented
00323 
00325   PixelContainerPointer m_Buffer;
00326 };
00327 } // end namespace itk
00328 
00329 // Define instantiation macro for this template.
00330 #define ITK_TEMPLATE_Image(_, EXPORT, TypeX, TypeY)     \
00331   namespace itk                                         \
00332   {                                                     \
00333   _( 2 ( class EXPORT Image< ITK_TEMPLATE_2 TypeX > ) ) \
00334   namespace Templates                                   \
00335   {                                                     \
00336   typedef Image< ITK_TEMPLATE_2 TypeX > Image##TypeY; \
00337   }                                                     \
00338   }
00339 
00340 #if ITK_TEMPLATE_EXPLICIT
00341 //template <class TPixel, unsigned int VImageDimension> const unsigned int
00342 // itk::Image<TPixel,VImageDimension>::ImageDimension;
00343 #include "Templates/itkImage+-.h"
00344 #endif
00345 
00346 #if ITK_TEMPLATE_TXX
00347 #include "itkImage.hxx"
00348 #endif
00349 
00350 #endif
00351