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: 2002/09/11 19:57:12 $
00007   Version:   $Revision: 1.44 $
00008 
00009   Copyright (c) 2002 Insight 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 "itkImageRegion.h"
00029 
00030 namespace itk
00031 {
00032 
00039 template <typename TImage>
00040 struct GetImageDimension
00041 {
00042   itkStaticConstMacro(ImageDimension, unsigned int, TImage::ImageDimension);
00043 }; 
00044   
00045 
00073 template<unsigned int VImageDimension=2>
00074 class ITK_EXPORT ImageBase : public DataObject
00075 {
00076 public:
00078   typedef ImageBase           Self;
00079   typedef DataObject  Superclass;
00080   typedef SmartPointer<Self>  Pointer;
00081   typedef SmartPointer<const Self>  ConstPointer;
00082   
00084   itkNewMacro(Self);
00085 
00087   itkTypeMacro(ImageBase, DataObject);
00088 
00093   itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension );
00094 
00096   typedef Index<VImageDimension>  IndexType;
00097   typedef typename IndexType::IndexValueType  IndexValueType;
00098   
00101   typedef Offset<VImageDimension>  OffsetType;
00102   typedef typename OffsetType::OffsetValueType OffsetValueType;
00103 
00105   typedef Size<VImageDimension>  SizeType;
00106   typedef typename SizeType::SizeValueType SizeValueType;
00107     
00109   typedef ImageRegion<VImageDimension>  RegionType;
00110 
00112   void Initialize();
00113 
00115   static unsigned int GetImageDimension() 
00116     { return VImageDimension; }
00117 
00122   virtual const double* GetSpacing() const;
00123 
00128   virtual const double * GetOrigin() const;
00129 
00136   virtual void SetLargestPossibleRegion(const RegionType &region);
00137 
00144   virtual const RegionType& GetLargestPossibleRegion() const
00145     { return m_LargestPossibleRegion;};
00146 
00150   virtual void SetBufferedRegion(const RegionType &region);
00151 
00155   virtual const RegionType& GetBufferedRegion() const
00156   { return m_BufferedRegion;};
00157   
00162   virtual void SetRequestedRegion(const RegionType &region);
00163 
00168   virtual void SetRequestedRegion(DataObject *data);
00169 
00174   virtual const RegionType& GetRequestedRegion() const
00175   { return m_RequestedRegion;};
00176 
00186   const OffsetValueType *GetOffsetTable() const { return m_OffsetTable; };
00187   
00190   OffsetValueType ComputeOffset(const IndexType &ind) const
00191   {
00192     // need to add bounds checking for the region/buffer?
00193     OffsetValueType offset=0;
00194     const IndexType &bufferedRegionIndex = m_BufferedRegion.GetIndex();
00195   
00196     // data is arranged as [][][][slice][row][col]
00197     // with Index[0] = col, Index[1] = row, Index[2] = slice
00198     for (int i=VImageDimension-1; i > 0; i--)
00199       {
00200       offset += (ind[i] - bufferedRegionIndex[i])*m_OffsetTable[i];
00201       }
00202     offset += (ind[0] - bufferedRegionIndex[0]);
00203 
00204     return offset;
00205   }
00206 
00209   IndexType ComputeIndex(OffsetValueType offset) const
00210   {
00211     IndexType index;
00212     const IndexType &bufferedRegionIndex = m_BufferedRegion.GetIndex();
00213     
00214     for (int i=VImageDimension-1; i > 0; i--)
00215       {
00216       index[i] = static_cast<IndexValueType>(offset / m_OffsetTable[i]);
00217       offset -= (index[i] * m_OffsetTable[i]);
00218       index[i] += bufferedRegionIndex[i];
00219       }
00220     index[0] = bufferedRegionIndex[0] + static_cast<IndexValueType>(offset);
00221 
00222     return index;
00223   }
00224 
00234   virtual void CopyInformation(const DataObject *data);
00235 
00243   virtual void UpdateOutputInformation();
00244 
00248   virtual void SetRequestedRegionToLargestPossibleRegion();
00249 
00259   virtual bool RequestedRegionIsOutsideOfTheBufferedRegion();
00260 
00269   virtual bool VerifyRequestedRegion();
00270   
00271 protected:
00272   ImageBase();
00273   ~ImageBase();
00274   virtual void PrintSelf(std::ostream& os, Indent indent) const;
00275 
00280   void ComputeOffsetTable();
00281 
00282 protected:
00284   double                m_Spacing[VImageDimension];
00285   double                m_Origin[VImageDimension];
00286 private:
00287   ImageBase(const Self&); //purposely not implemented
00288   void operator=(const Self&); //purposely not implemented
00289 
00290   OffsetValueType  m_OffsetTable[VImageDimension+1];
00291 
00292   RegionType          m_LargestPossibleRegion;
00293   RegionType          m_RequestedRegion;
00294   RegionType          m_BufferedRegion;
00295 };
00296 
00297 
00298 } // end namespace itk
00299 
00300 #ifndef ITK_MANUAL_INSTANTIATION
00301 #include "itkImageBase.txx"
00302 #endif
00303 
00304 #endif
00305 

Generated at Fri May 21 01:14:52 2004 for ITK by doxygen 1.2.15 written by Dimitri van Heesch, © 1997-2000