ITK  4.2.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 {
00079 template< class TPixel, unsigned int VImageDimension = 2 >
00080 class ITK_EXPORT Image:public ImageBase< VImageDimension >
00081 {
00082 public:
00084   typedef Image                        Self;
00085   typedef ImageBase< VImageDimension > Superclass;
00086   typedef SmartPointer< Self >         Pointer;
00087   typedef SmartPointer< const Self >   ConstPointer;
00088   typedef WeakPointer< const Self >    ConstWeakPointer;
00089 
00091   itkNewMacro(Self);
00092 
00094   itkTypeMacro(Image, ImageBase);
00095 
00098   typedef TPixel PixelType;
00099 
00101   typedef TPixel ValueType;
00102 
00107   typedef TPixel InternalPixelType;
00108 
00109   typedef PixelType IOPixelType;
00110 
00113   typedef DefaultPixelAccessor< PixelType >   AccessorType;
00114   typedef DefaultPixelAccessorFunctor< Self > AccessorFunctorType;
00115 
00118   typedef NeighborhoodAccessorFunctor< Self > NeighborhoodAccessorFunctorType;
00119 
00124   itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
00125 
00127   typedef typename Superclass::IndexType      IndexType;
00128   typedef typename Superclass::IndexValueType IndexValueType;
00129 
00131   typedef typename Superclass::OffsetType OffsetType;
00132 
00134   typedef typename Superclass::SizeType      SizeType;
00135   typedef typename Superclass::SizeValueType SizeValueType;
00136 
00138   typedef ImportImageContainer< SizeValueType, PixelType > PixelContainer;
00139 
00141   typedef typename Superclass::DirectionType DirectionType;
00142 
00145   typedef typename Superclass::RegionType RegionType;
00146 
00149   typedef typename Superclass::SpacingType      SpacingType;
00150   typedef typename Superclass::SpacingValueType SpacingValueType;
00151 
00154   typedef typename Superclass::PointType PointType;
00155 
00157   typedef typename PixelContainer::Pointer      PixelContainerPointer;
00158   typedef typename PixelContainer::ConstPointer PixelContainerConstPointer;
00159 
00161   typedef typename Superclass::OffsetValueType OffsetValueType;
00162 
00168   template <class UPixelType, unsigned int UImageDimension = VImageDimension>
00169   struct Rebind
00170     {
00171       typedef itk::Image<UPixelType, UImageDimension>  Type;
00172     };
00173 
00174 
00177   void Allocate();
00178 
00181   virtual void Initialize();
00182 
00185   void FillBuffer(const TPixel & value);
00186 
00192   void SetPixel(const IndexType & index, const TPixel & value)
00193   {
00194     OffsetValueType offset = this->ComputeOffset(index);
00195     ( *m_Buffer )[offset] = value;
00196   }
00197 
00202   const TPixel & GetPixel(const IndexType & index) const
00203   {
00204     OffsetValueType offset = this->ComputeOffset(index);
00205     return ( ( *m_Buffer )[offset] );
00206   }
00207 
00212   TPixel & GetPixel(const IndexType & index)
00213   {
00214     OffsetValueType offset = this->ComputeOffset(index);
00215     return ( ( *m_Buffer )[offset] );
00216   }
00217 
00222   TPixel & operator[](const IndexType & index)
00223   { return this->GetPixel(index); }
00224 
00229   const TPixel & operator[](const IndexType & index) const
00230   { return this->GetPixel(index); }
00231 
00234   virtual TPixel * GetBufferPointer()
00235   { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00236   virtual const TPixel * GetBufferPointer() const
00237   { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00239 
00241   PixelContainer * GetPixelContainer()
00242   { return m_Buffer.GetPointer(); }
00243 
00244   const PixelContainer * GetPixelContainer() const
00245   { return m_Buffer.GetPointer(); }
00246 
00249   void SetPixelContainer(PixelContainer *container);
00250 
00261   virtual void Graft(const DataObject *data);
00262 
00264   AccessorType GetPixelAccessor(void)
00265   { return AccessorType(); }
00266 
00268   const AccessorType GetPixelAccessor(void) const
00269   { return AccessorType(); }
00270 
00272   NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
00273   { return NeighborhoodAccessorFunctorType(); }
00274 
00276   const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
00277   { return NeighborhoodAccessorFunctorType(); }
00278 
00279   virtual unsigned int GetNumberOfComponentsPerPixel() const;
00280 
00281 protected:
00282   Image();
00283   void PrintSelf(std::ostream & os, Indent indent) const;
00284 
00285   virtual ~Image() {}
00286 
00292   virtual void ComputeIndexToPhysicalPointMatrices();
00293 
00294 private:
00295   Image(const Self &);          //purposely not implemented
00296   void operator=(const Self &); //purposely not implemented
00297 
00299   PixelContainerPointer m_Buffer;
00300 };
00301 } // end namespace itk
00302 
00303 // Define instantiation macro for this template.
00304 #define ITK_TEMPLATE_Image(_, EXPORT, TypeX, TypeY)     \
00305   namespace itk                                         \
00306   {                                                     \
00307   _( 2 ( class EXPORT Image< ITK_TEMPLATE_2 TypeX > ) ) \
00308   namespace Templates                                   \
00309   {                                                     \
00310   typedef Image< ITK_TEMPLATE_2 TypeX > Image##TypeY; \
00311   }                                                     \
00312   }
00313 
00314 #if ITK_TEMPLATE_EXPLICIT
00315 //template <class TPixel, unsigned int VImageDimension> const unsigned int
00316 // itk::Image<TPixel,VImageDimension>::ImageDimension;
00317 #include "Templates/itkImage+-.h"
00318 #endif
00319 
00320 #if ITK_TEMPLATE_TXX
00321 #include "itkImage.hxx"
00322 #endif
00323 
00324 #endif
00325