00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #ifndef __itkImage_h
00018 #define __itkImage_h
00019
00020 #include "itkImageBase.h"
00021 #include "itkImageRegion.h"
00022 #include "itkImportImageContainer.h"
00023 #include "itkDefaultPixelAccessor.h"
00024 #include "itkPoint.h"
00025 #include "itkContinuousIndex.h"
00026 #include "itkTransform.h"
00027
00028 namespace itk
00029 {
00030
00073 template <class TPixel, unsigned int VImageDimension=2>
00074 class ITK_EXPORT Image : public ImageBase<VImageDimension>
00075 {
00076 public:
00078 typedef Image Self;
00079 typedef ImageBase<VImageDimension> Superclass;
00080 typedef SmartPointer<Self> Pointer;
00081 typedef SmartPointer<const Self> ConstPointer;
00082
00084 itkNewMacro(Self);
00085
00087 itkTypeMacro(Image, ImageBase);
00088
00091 typedef TPixel PixelType;
00092
00094 typedef TPixel ValueType ;
00095
00100 typedef TPixel InternalPixelType;
00101
00104 typedef DefaultPixelAccessor< PixelType > AccessorType;
00105
00110 itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
00111
00113 typedef ImportImageContainer<unsigned long, PixelType> PixelContainer;
00114
00116 typedef Index<VImageDimension> IndexType;
00117
00119 typedef Offset<VImageDimension> OffsetType;
00120
00122 typedef Size<VImageDimension> SizeType;
00123
00125 typedef ImageRegion<VImageDimension> RegionType;
00126
00128 typedef typename PixelContainer::Pointer PixelContainerPointer;
00129
00136 typedef Transform< double, VImageDimension, VImageDimension > TransformType;
00137 typedef typename TransformType::Pointer TransformPointer;
00138
00139
00142 void Allocate();
00143
00147 void SetRegions(RegionType region)
00148 {
00149 this->SetLargestPossibleRegion(region);
00150 this->SetBufferedRegion(region);
00151 this->SetRequestedRegion(region);
00152 };
00153
00154 void SetRegions(SizeType size)
00155 {
00156 RegionType region; region.SetSize(size);
00157 this->SetLargestPossibleRegion(region);
00158 this->SetBufferedRegion(region);
00159 this->SetRequestedRegion(region);
00160 };
00161
00164 virtual void Initialize();
00165
00167 void FillBuffer (const TPixel& value);
00168
00173 void SetPixel(const IndexType &index, const TPixel& value)
00174 {
00175 typename Superclass::OffsetValueType offset = this->ComputeOffset(index);
00176 (*m_Buffer)[offset] = value;
00177 }
00178
00183 const TPixel& GetPixel(const IndexType &index) const
00184 {
00185 typename Superclass::OffsetValueType offset = this->ComputeOffset(index);
00186 return ( (*m_Buffer)[offset] );
00187 }
00188
00193 TPixel& GetPixel(const IndexType &index)
00194 {
00195 typename Superclass::OffsetValueType offset = this->ComputeOffset(index);
00196 return ( (*m_Buffer)[offset] );
00197 }
00198
00203 TPixel & operator[](const IndexType &index)
00204 { return this->GetPixel(index); }
00205
00210 const TPixel& operator[](const IndexType &index) const
00211 { return this->GetPixel(index); }
00212
00215 TPixel *GetBufferPointer()
00216 { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00217 const TPixel *GetBufferPointer() const
00218 { return m_Buffer ? m_Buffer->GetBufferPointer() : 0; }
00219
00221 PixelContainer* GetPixelContainer()
00222 { return m_Buffer.GetPointer(); }
00223
00226 void SetPixelContainer( PixelContainer *container );
00227
00229 AccessorType GetPixelAccessor( void )
00230 { return AccessorType(); }
00231
00233 const AccessorType GetPixelAccessor( void ) const
00234 { return AccessorType(); }
00235
00240 virtual void SetSpacing( const double spacing[VImageDimension] );
00241 virtual void SetSpacing( const float spacing[VImageDimension] );
00242
00247 virtual void SetOrigin( const double origin[VImageDimension] );
00248 virtual void SetOrigin( const float origin[VImageDimension] );
00249
00255 TransformPointer GetIndexToPhysicalTransform(void) const;
00256
00262 TransformPointer GetPhysicalToIndexTransform(void) const;
00263
00265 itkSetObjectMacro(IndexToPhysicalTransform, TransformType );
00266 itkSetObjectMacro(PhysicalToIndexTransform, TransformType );
00267
00274 virtual void RebuildTransforms(void) throw ( std::exception );
00275
00286 template<class TCoordRep>
00287 bool TransformPhysicalPointToContinuousIndex(
00288 const Point<TCoordRep, VImageDimension>& point,
00289 ContinuousIndex<TCoordRep, VImageDimension>& index ) const
00290 {
00291
00292 if ( !m_PhysicalToIndexTransform )
00293 {
00294 itkExceptionMacro("The Image lacks a PhysicalToIndexTransform");
00295 }
00296
00297 typename TransformType::InputPointType inputPoint =
00298 m_PhysicalToIndexTransform->TransformPoint(point) ;
00299
00300
00301 for (unsigned int i = 0 ; i < VImageDimension ; i++)
00302 {
00303 index[i] = static_cast<TCoordRep>(inputPoint[i]);
00304 }
00305
00306
00307 const bool isInside =
00308 this->GetLargestPossibleRegion().IsInside( index );
00309
00310 return isInside;
00311 }
00312
00323 template<class TCoordRep>
00324 bool TransformPhysicalPointToIndex(
00325 const Point<TCoordRep, VImageDimension>& point,
00326 IndexType & index ) const
00327 {
00328
00329 if ( !m_PhysicalToIndexTransform )
00330 {
00331 itkExceptionMacro("The Image lacks a PhysicalToIndexTransform");
00332 }
00333
00334 typename TransformType::InputPointType inputPoint =
00335 m_PhysicalToIndexTransform->TransformPoint(point) ;
00336
00337 typedef typename IndexType::IndexValueType IndexValueType;
00338
00339
00340 for (unsigned int i = 0 ; i < VImageDimension ; i++)
00341 {
00342 index[i] = static_cast<IndexValueType>(inputPoint[i]);
00343 }
00344
00345
00346 const bool isInside =
00347 this->GetLargestPossibleRegion().IsInside( index );
00348
00349 return isInside;
00350 }
00351
00363 template<class TCoordRep>
00364 void TransformContinuousIndexToPhysicalPoint(
00365 const ContinuousIndex<TCoordRep, VImageDimension>& index,
00366 Point<TCoordRep, VImageDimension>& point ) const
00367 {
00368
00369 if ( !m_IndexToPhysicalTransform )
00370 {
00371 itkExceptionMacro("The Image lacks a IndexToPhysicalTransform");
00372 }
00373
00374 typename TransformType::InputPointType inputPoint;
00375
00376
00377 for (unsigned int i = 0 ; i < VImageDimension ; i++)
00378 { inputPoint[i] = index[i]; }
00379
00380
00381 typename TransformType::OutputPointType outputPoint =
00382 m_IndexToPhysicalTransform->TransformPoint(inputPoint) ;
00383
00384
00385 point = outputPoint;
00386 }
00387
00400 template<class TCoordRep>
00401 void TransformIndexToPhysicalPoint(
00402 const IndexType & index,
00403 Point<TCoordRep, VImageDimension>& point ) const
00404 {
00405 if ( !m_IndexToPhysicalTransform )
00406 {
00407 itkExceptionMacro("The Image lacks a IndexToPhysicalTransform");
00408 }
00409
00410 typename TransformType::InputPointType inputPoint;
00411
00412
00413 for (unsigned int i = 0 ; i < VImageDimension ; i++)
00414 {
00415 inputPoint[i] = index[i];
00416 }
00417
00418
00419 typename TransformType::OutputPointType outputPoint =
00420 m_IndexToPhysicalTransform->TransformPoint(inputPoint) ;
00421
00422
00423 for (unsigned int i = 0 ; i < VImageDimension ; ++i)
00424 {
00425 point[i] = outputPoint[i];
00426 }
00427 }
00428
00445 virtual void CopyInformation(const DataObject *data);
00446
00447 protected:
00448 Image();
00449 virtual ~Image();
00450 void PrintSelf(std::ostream& os, Indent indent) const;
00451
00452 private:
00453 Image(const Self&);
00454 void operator=(const Self&);
00455
00457 PixelContainerPointer m_Buffer;
00458
00460 TransformPointer m_IndexToPhysicalTransform;
00461 TransformPointer m_PhysicalToIndexTransform;
00462 };
00463
00464 }
00465
00466 #ifndef ITK_MANUAL_INSTANTIATION
00467 #include "itkImage.txx"
00468 #endif
00469
00470 #endif
00471