ITK  4.0.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 
00171   void Allocate();
00172 
00176   void SetRegions(RegionType region)
00177   {
00178     this->SetLargestPossibleRegion(region);
00179     this->SetBufferedRegion(region);
00180     this->SetRequestedRegion(region);
00181   }
00183 
00184   void SetRegions(SizeType size)
00185   {
00186     RegionType region; region.SetSize(size);
00187 
00188     this->SetLargestPossibleRegion(region);
00189     this->SetBufferedRegion(region);
00190     this->SetRequestedRegion(region);
00191   }
00192 
00195   virtual void Initialize();
00196 
00199   void FillBuffer(const TPixel & value);
00200 
00206   void SetPixel(const IndexType & index, const TPixel & value)
00207   {
00208     OffsetValueType offset = this->ComputeOffset(index);
00209     ( *m_Buffer )[offset] = value;
00210   }
00211 
00216   const TPixel & GetPixel(const IndexType & index) const
00217   {
00218     OffsetValueType offset = this->ComputeOffset(index);
00219     return ( ( *m_Buffer )[offset] );
00220   }
00221 
00226   TPixel & GetPixel(const IndexType & index)
00227   {
00228     OffsetValueType offset = this->ComputeOffset(index);
00229     return ( ( *m_Buffer )[offset] );
00230   }
00231 
00236   TPixel & operator[](const IndexType & index)
00237   { return this->GetPixel(index); }
00238 
00243   const TPixel & operator[](const IndexType & index) const
00244   { return this->GetPixel(index); }
00245 
00248   virtual TPixel * GetBufferPointer()
00249   { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00250   virtual const TPixel * GetBufferPointer() const
00251   { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00253 
00255   PixelContainer * GetPixelContainer()
00256   { return m_Buffer.GetPointer(); }
00257 
00258   const PixelContainer * GetPixelContainer() const
00259   { return m_Buffer.GetPointer(); }
00260 
00263   void SetPixelContainer(PixelContainer *container);
00264 
00275   virtual void Graft(const DataObject *data);
00276 
00278   AccessorType GetPixelAccessor(void)
00279   { return AccessorType(); }
00280 
00282   const AccessorType GetPixelAccessor(void) const
00283   { return AccessorType(); }
00284 
00286   NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
00287   { return NeighborhoodAccessorFunctorType(); }
00288 
00290   const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
00291   { return NeighborhoodAccessorFunctorType(); }
00292 
00293   virtual unsigned int GetNumberOfComponentsPerPixel() const;
00294 
00295 protected:
00296   Image();
00297   void PrintSelf(std::ostream & os, Indent indent) const;
00298 
00299   virtual ~Image() {}
00300 
00306   virtual void ComputeIndexToPhysicalPointMatrices();
00307 
00308 private:
00309   Image(const Self &);          //purposely not implemented
00310   void operator=(const Self &); //purposely not implemented
00311 
00313   PixelContainerPointer m_Buffer;
00314 };
00315 } // end namespace itk
00316 
00317 // Define instantiation macro for this template.
00318 #define ITK_TEMPLATE_Image(_, EXPORT, TypeX, TypeY)     \
00319   namespace itk                                         \
00320   {                                                     \
00321   _( 2 ( class EXPORT Image< ITK_TEMPLATE_2 TypeX > ) ) \
00322   namespace Templates                                   \
00323   {                                                     \
00324   typedef Image< ITK_TEMPLATE_2 TypeX > Image##TypeY; \
00325   }                                                     \
00326   }
00327 
00328 #if ITK_TEMPLATE_EXPLICIT
00329 //template <class TPixel, unsigned int VImageDimension> const unsigned int
00330 // itk::Image<TPixel,VImageDimension>::ImageDimension;
00331 #include "Templates/itkImage+-.h"
00332 #endif
00333 
00334 #if ITK_TEMPLATE_TXX
00335 #include "itkImage.hxx"
00336 #endif
00337 
00338 #endif
00339