Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

itkImageBase.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkImageBase.h,v $
00005   Language:  C++
00006   Date:      $Date: 2003/12/10 15:21:33 $
00007   Version:   $Revision: 1.50 $
00008 
00009   Copyright (c) Insight Software Consortium. All rights reserved.
00010   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
00011 
00012   Portions of this code are covered under the VTK copyright.
00013   See VTKCopyright.txt or http://www.kitware.com/VTKCopyright.htm for details.
00014 
00015      This software is distributed WITHOUT ANY WARRANTY; without even 
00016      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
00017      PURPOSE.  See the above copyright notices for more information.
00018 
00019 =========================================================================*/
00020 #ifndef __itkImageBase_h
00021 #define __itkImageBase_h
00022 
00023 #include "itkDataObject.h"
00024 #include "itkProcessObject.h"
00025 #include "itkIndex.h"
00026 #include "itkOffset.h"
00027 #include "itkSize.h"
00028 #include "itkFixedArray.h"
00029 #include "itkPoint.h"
00030 #include "itkImageRegion.h"
00031 
00032 namespace itk
00033 {
00034 
00041 template <typename TImage>
00042 struct GetImageDimension
00043 {
00044   itkStaticConstMacro(ImageDimension, unsigned int, TImage::ImageDimension);
00045 }; 
00046   
00047 
00075 template<unsigned int VImageDimension=2>
00076 class ITK_EXPORT ImageBase : public DataObject
00077 {
00078 public:
00080   typedef ImageBase           Self;
00081   typedef DataObject  Superclass;
00082   typedef SmartPointer<Self>  Pointer;
00083   typedef SmartPointer<const Self>  ConstPointer;
00084   
00086   itkNewMacro(Self);
00087 
00089   itkTypeMacro(ImageBase, DataObject);
00090 
00095   itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension );
00096 
00098   typedef Index<VImageDimension>  IndexType;
00099   typedef typename IndexType::IndexValueType  IndexValueType;
00100   
00103   typedef Offset<VImageDimension>  OffsetType;
00104   typedef typename OffsetType::OffsetValueType OffsetValueType;
00105 
00107   typedef Size<VImageDimension>  SizeType;
00108   typedef typename SizeType::SizeValueType SizeValueType;
00109     
00111   typedef ImageRegion<VImageDimension>  RegionType;
00112 
00115   typedef Vector<double, VImageDimension> SpacingType;
00116 
00119   typedef Point<double, VImageDimension> PointType;
00120 
00122   void Initialize();
00123 
00125   static unsigned int GetImageDimension() 
00126     { return VImageDimension; }
00127 
00132   itkGetConstReferenceMacro(Spacing, SpacingType);
00133 
00138   itkGetConstReferenceMacro(Origin, PointType);
00139 
00146   virtual void SetLargestPossibleRegion(const RegionType &region);
00147 
00154   virtual const RegionType& GetLargestPossibleRegion() const
00155     { return m_LargestPossibleRegion;};
00156 
00160   virtual void SetBufferedRegion(const RegionType &region);
00161 
00165   virtual const RegionType& GetBufferedRegion() const
00166   { return m_BufferedRegion;};
00167   
00175   virtual void SetRequestedRegion(const RegionType &region);
00176 
00184   virtual void SetRequestedRegion(DataObject *data);
00185 
00190   virtual const RegionType& GetRequestedRegion() const
00191   { return m_RequestedRegion;};
00192 
00202   const OffsetValueType *GetOffsetTable() const { return m_OffsetTable; };
00203   
00206   OffsetValueType ComputeOffset(const IndexType &ind) const
00207   {
00208     // need to add bounds checking for the region/buffer?
00209     OffsetValueType offset=0;
00210     const IndexType &bufferedRegionIndex = m_BufferedRegion.GetIndex();
00211   
00212     // data is arranged as [][][][slice][row][col]
00213     // with Index[0] = col, Index[1] = row, Index[2] = slice
00214     for (int i=VImageDimension-1; i > 0; i--)
00215       {
00216       offset += (ind[i] - bufferedRegionIndex[i])*m_OffsetTable[i];
00217       }
00218     offset += (ind[0] - bufferedRegionIndex[0]);
00219 
00220     return offset;
00221   }
00222 
00225   IndexType ComputeIndex(OffsetValueType offset) const
00226   {
00227     IndexType index;
00228     const IndexType &bufferedRegionIndex = m_BufferedRegion.GetIndex();
00229     
00230     for (int i=VImageDimension-1; i > 0; i--)
00231       {
00232       index[i] = static_cast<IndexValueType>(offset / m_OffsetTable[i]);
00233       offset -= (index[i] * m_OffsetTable[i]);
00234       index[i] += bufferedRegionIndex[i];
00235       }
00236     index[0] = bufferedRegionIndex[0] + static_cast<IndexValueType>(offset);
00237 
00238     return index;
00239   }
00240 
00250   virtual void CopyInformation(const DataObject *data);
00251 
00259   virtual void UpdateOutputInformation();
00260 
00264   virtual void SetRequestedRegionToLargestPossibleRegion();
00265 
00275   virtual bool RequestedRegionIsOutsideOfTheBufferedRegion();
00276 
00285   virtual bool VerifyRequestedRegion();
00286   
00287 protected:
00288   ImageBase();
00289   ~ImageBase();
00290   virtual void PrintSelf(std::ostream& os, Indent indent) const;
00291 
00296   void ComputeOffsetTable();
00297 
00298 protected:
00301   SpacingType  m_Spacing;
00302   PointType   m_Origin;
00303 private:
00304   ImageBase(const Self&); //purposely not implemented
00305   void operator=(const Self&); //purposely not implemented
00306 
00307   OffsetValueType  m_OffsetTable[VImageDimension+1];
00308 
00309   RegionType          m_LargestPossibleRegion;
00310   RegionType          m_RequestedRegion;
00311   RegionType          m_BufferedRegion;
00312 };
00313 
00314 #ifdef ITK_EXPLICIT_INSTANTIATION
00315    extern template class ImageBase<2>;
00316    extern template class ImageBase<3>;
00317 #endif
00318 
00319 } // end namespace itk
00320 
00321 #ifndef ITK_MANUAL_INSTANTIATION
00322 #include "itkImageBase.txx"
00323 #endif
00324 
00325 #endif
00326 

Generated at Sun Jan 25 13:19:36 2004 for ITK by doxygen 1.3.3 written by Dimitri van Heesch, © 1997-2000