ITK  4.5.0
Insight Segmentation and Registration Toolkit
itkImageBase.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 /*=========================================================================
19  *
20  * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21  *
22  * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23  *
24  * For complete copyright, license and disclaimer of warranty information
25  * please refer to the NOTICE file at the top of the ITK source tree.
26  *
27  *=========================================================================*/
28 #ifndef __itkImageBase_h
29 #define __itkImageBase_h
30 
31 #include "itkDataObject.h"
32 
33 #include "itkImageRegion.h"
34 #include "itkMatrix.h"
35 #include "itkObjectFactory.h"
36 #include "itkOffset.h"
37 #include "itkFixedArray.h"
38 #include "itkImageHelper.h"
39 #include "itkFloatTypes.h"
40 
41 //HACK: vnl/vnl_matrix_fixed.txx is needed here?
42 // to avoid undefined symbol vnl_matrix_fixed<double, 8u, 8u>::set_identity()", referenced from
43 #include "vnl/vnl_matrix_fixed.txx"
44 
46 
47 namespace itk
48 {
49 
50 /* Forward declaration (ImageTransformHelper include's ImageBase) */
51 template< unsigned int NImageDimension, unsigned int R, unsigned int C, typename TPointValue, typename TMatrixValue >
53 
111 template< unsigned int VImageDimension = 2 >
112 class ImageBase:public DataObject
113 {
114 public:
116  typedef ImageBase Self;
120 
122  itkNewMacro(Self);
123 
125  itkTypeMacro(ImageBase, DataObject);
126 
131  itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
132 
136 
141 
145 
148 
155 
160 
165 
167  void Initialize();
168 
170  static unsigned int GetImageDimension()
171  { return VImageDimension; }
172 
177  itkSetMacro(Origin, PointType);
178  virtual void SetOrigin(const double origin[VImageDimension]);
179  virtual void SetOrigin(const float origin[VImageDimension]);
181 
207  virtual void SetDirection(const DirectionType & direction);
208 
212  itkGetConstReferenceMacro(Direction, DirectionType);
213 
217  itkGetConstReferenceMacro(InverseDirection, DirectionType);
218 
223  itkGetConstReferenceMacro(Spacing, SpacingType);
224 
229  itkGetConstReferenceMacro(Origin, PointType);
230 
237  virtual void Allocate() {}
238 
245  virtual void SetLargestPossibleRegion(const RegionType & region);
246 
253  virtual const RegionType & GetLargestPossibleRegion() const
254  { return m_LargestPossibleRegion; }
255 
259  virtual void SetBufferedRegion(const RegionType & region);
260 
264  virtual const RegionType & GetBufferedRegion() const
265  { return m_BufferedRegion; }
266 
274  virtual void SetRequestedRegion(const RegionType & region);
275 
283  virtual void SetRequestedRegion( const DataObject *data );
284 
289  virtual const RegionType & GetRequestedRegion() const
290  { return m_RequestedRegion; }
291 
295  virtual void SetRegions(const RegionType& region)
296  {
297  this->SetLargestPossibleRegion(region);
298  this->SetBufferedRegion(region);
299  this->SetRequestedRegion(region);
300  }
302 
303  virtual void SetRegions(const SizeType& size)
304  {
305  RegionType region; region.SetSize(size);
306 
307  this->SetLargestPossibleRegion(region);
308  this->SetBufferedRegion(region);
309  this->SetRequestedRegion(region);
310  }
311 
322  const OffsetValueType * GetOffsetTable() const { return m_OffsetTable; }
324 
331  inline OffsetValueType ComputeOffset(const IndexType & ind) const
332  {
333  OffsetValueType offset = 0;
334 
336  ind,
338  offset);
339  return offset;
340  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
341  * Leaving here for documentation purposes
342  * OffsetValueType ComputeOffset(const IndexType & ind) const
343  * {
344  * // need to add bounds checking for the region/buffer?
345  * OffsetValueType offset = 0;
346  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
347  * // data is arranged as [][][][slice][row][col]
348  * // with Index[0] = col, Index[1] = row, Index[2] = slice
349  * for ( int i = VImageDimension - 1; i > 0; i-- )
350  * {
351  * offset += ( ind[i] - bufferedRegionIndex[i] ) * m_OffsetTable[i];
352  * }
353  * offset += ( ind[0] - bufferedRegionIndex[0] );
354  * return offset;
355  * }
356  */
357  }
358 
367  {
368  IndexType index;
369  const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
371 
373  offset,
375  index);
376  return index;
377  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
378  * Leaving here for documentation purposes
379  * IndexType ComputeIndex(OffsetValueType offset) const
380  * {
381  * IndexType index;
382  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
383  * for ( int i = VImageDimension - 1; i > 0; i-- )
384  * {
385  * index[i] = static_cast< IndexValueType >( offset / m_OffsetTable[i] );
386  * offset -= ( index[i] * m_OffsetTable[i] );
387  * index[i] += bufferedRegionIndex[i];
388  * }
389  * index[0] = bufferedRegionIndex[0] + static_cast< IndexValueType >( offset );
390  * return index;
391  * }
392  */
393 
394  }
395 
402  virtual void SetSpacing(const SpacingType & spacing);
403  virtual void SetSpacing(const double spacing[VImageDimension]);
404  virtual void SetSpacing(const float spacing[VImageDimension]);
406 
411  template< typename TCoordRep >
414  IndexType & index) const
415  {
418 
419  // Now, check to see if the index is within allowed bounds
420  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
421  return isInside;
422  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
423  * Leaving here for documentation purposes
424  * template< typename TCoordRep >
425  * bool TransformPhysicalPointToIndex(
426  * const Point< TCoordRep, VImageDimension > & point,
427  * IndexType & index) const
428  * {
429  * for ( unsigned int i = 0; i < VImageDimension; i++ )
430  * {
431  * TCoordRep sum = NumericTraits< TCoordRep >::Zero;
432  * for ( unsigned int j = 0; j < VImageDimension; j++ )
433  * {
434  * sum += this->m_PhysicalPointToIndex[i][j] * ( point[j] - this->m_Origin[j] );
435  * }
436  * index[i] = Math::RoundHalfIntegerUp< IndexValueType >(sum);
437  * }
438  * // Now, check to see if the index is within allowed bounds
439  * const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
440  * return isInside;
441  * }
442  */
443  }
444 
449  template< typename TCoordRep, typename TIndexRep >
453  {
455 
456  for ( unsigned int k = 0; k < VImageDimension; k++ )
457  {
458  cvector[k] = point[k] - this->m_Origin[k];
459  }
460  cvector = m_PhysicalPointToIndex * cvector;
461  for ( unsigned int i = 0; i < VImageDimension; i++ )
462  {
463  index[i] = static_cast< TIndexRep >( cvector[i] );
464  }
465 
466  // Now, check to see if the index is within allowed bounds
467  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
468 
469  return isInside;
470  }
471 
476  template< typename TCoordRep, typename TIndexRep >
480  {
481  for ( unsigned int r = 0; r < VImageDimension; r++ )
482  {
483  TCoordRep sum = NumericTraits< TCoordRep >::Zero;
484  for ( unsigned int c = 0; c < VImageDimension; c++ )
485  {
486  sum += this->m_IndexToPhysicalPoint(r, c) * index[c];
487  }
488  point[r] = sum + this->m_Origin[r];
489  }
490  }
492 
498  template< typename TCoordRep >
500  const IndexType & index,
502  {
505  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
506  * Leaving here for documentation purposes
507  * template< typename TCoordRep >
508  * void TransformIndexToPhysicalPoint(
509  * const IndexType & index,
510  * Point< TCoordRep, VImageDimension > & point) const
511  * {
512  * for ( unsigned int i = 0; i < VImageDimension; i++ )
513  * {
514  * point[i] = this->m_Origin[i];
515  * for ( unsigned int j = 0; j < VImageDimension; j++ )
516  * {
517  * point[i] += m_IndexToPhysicalPoint[i][j] * index[j];
518  * }
519  * }
520  * }
521  */
522  }
524 
536  template< typename TCoordRep >
538  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
539  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
540  {
541  //
542  //TODO: This temporary implementation should be replaced with Template
543  // MetaProgramming.
544  //
545  const DirectionType & direction = this->GetDirection();
546 
547  for ( unsigned int i = 0; i < VImageDimension; i++ )
548  {
549  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
550  CoordSumType sum = NumericTraits< CoordSumType >::Zero;
551  for ( unsigned int j = 0; j < VImageDimension; j++ )
552  {
553  sum += direction[i][j] * inputGradient[j];
554  }
555  outputGradient[i] = static_cast< TCoordRep >( sum );
556  }
557  }
558 
567  template< typename TCoordRep >
569  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
570  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
571  {
572  //
573  //TODO: This temporary implementation should be replaced with Template
574  // MetaProgramming.
575  //
576  const DirectionType & inverseDirection = this->GetInverseDirection();
577 
578  for ( unsigned int i = 0; i < VImageDimension; i++ )
579  {
580  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
581  CoordSumType sum = NumericTraits< CoordSumType >::Zero;
582  for ( unsigned int j = 0; j < VImageDimension; j++ )
583  {
584  sum += inverseDirection[i][j] * inputGradient[j];
585  }
586  outputGradient[i] = static_cast< TCoordRep >( sum );
587  }
588  }
589 
599  virtual void CopyInformation(const DataObject *data);
600 
611  virtual void Graft(const DataObject *data);
612 
620  virtual void UpdateOutputInformation();
621 
629  virtual void UpdateOutputData();
630 
635 
646 
655  virtual bool VerifyRequestedRegion();
656 
675  virtual unsigned int GetNumberOfComponentsPerPixel() const;
676  virtual void SetNumberOfComponentsPerPixel(unsigned int);
678 
679 protected:
680  ImageBase();
681  ~ImageBase();
682  virtual void PrintSelf(std::ostream & os, Indent indent) const;
683 
688  void ComputeOffsetTable();
689 
696 
697 protected:
705 
710 
715  virtual void InitializeBufferedRegion(void);
716 
717 private:
718  ImageBase(const Self &); //purposely not implemented
719  void operator=(const Self &); //purposely not implemented
720 
721  void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
722  {
723  SpacingType s(spacing);
724  this->SetSpacing(s);
725  }
726 
727  template <typename TSpacingValue>
728  void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
729  {
731  SpacingType s;
732  s.CastFrom(sf);
733  this->SetSpacing(s);
734  }
735 
736  OffsetValueType m_OffsetTable[VImageDimension + 1];
737 
741 };
742 } // end namespace itk
743 
744 #ifndef ITK_MANUAL_INSTANTIATION
745 #include "itkImageBase.hxx"
746 #endif
747 
748 #endif
void SetSize(const SizeType &size)
static void TransformIndexToPhysicalPoint(const MatrixType &matrix, const OriginType &origin, const IndexType &index, DoublePoint &point)
virtual void SetRegions(const RegionType &region)
Definition: itkImageBase.h:295
virtual const DirectionType & GetInverseDirection()
OffsetValueType m_OffsetTable[VImageDimension+1]
Definition: itkImageBase.h:736
const OffsetValueType * GetOffsetTable() const
Definition: itkImageBase.h:322
Index< VImageDimension > IndexType
Definition: itkImageBase.h:134
const IndexType & GetIndex() const
virtual void Graft(const DataObject *data)
void operator=(const Self &)
Represent the offset between two n-dimensional indexes in a n-dimensional image.
Definition: itkOffset.h:55
Fast index/physical index computation.
Definition: itkImageBase.h:52
static void ComputeOffset(const IndexType &bufferedRegionIndex, const IndexType &index, const OffsetValueType offsetTable[], OffsetValueType &offset)
itk::SizeValueType SizeValueType
Definition: itkSize.h:60
virtual void PrintSelf(std::ostream &os, Indent indent) const
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
RegionType m_BufferedRegion
Definition: itkImageBase.h:740
virtual void SetDirection(const DirectionType &direction)
bool IsInside(const IndexType &index) const
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:253
virtual const DirectionType & GetDirection()
SpacingType m_Spacing
Definition: itkImageBase.h:701
An image region represents a structured region of data.
RegionType m_RequestedRegion
Definition: itkImageBase.h:739
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:153
virtual void SetRequestedRegionToLargestPossibleRegion()
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:499
PointType m_Origin
Definition: itkImageBase.h:702
virtual bool VerifyRequestedRegion()
DirectionType m_IndexToPhysicalPoint
Definition: itkImageBase.h:708
static void ComputeIndex(const IndexType &bufferedRegionIndex, OffsetValueType offset, const OffsetValueType offsetTable[], IndexType &index)
Size< VImageDimension > SizeType
Definition: itkImageBase.h:143
DirectionType m_InverseDirection
Definition: itkImageBase.h:704
virtual bool RequestedRegionIsOutsideOfTheBufferedRegion()
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:159
virtual void SetSpacing(const SpacingType &spacing)
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Definition: itkImageBase.h:164
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
virtual void SetRegions(const SizeType &size)
Definition: itkImageBase.h:303
IndexType ComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:366
DirectionType m_Direction
Definition: itkImageBase.h:703
::itk::IndexValueType IndexValueType
Definition: itkIndex.h:79
bool TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const
Get the continuous index from a physical point.
Definition: itkImageBase.h:450
virtual void Allocate()
Definition: itkImageBase.h:237
virtual void SetNumberOfComponentsPerPixel(unsigned int)
void ComputeOffsetTable()
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:264
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:139
SmartPointer< const Self > ConstPointer
Definition: itkImageBase.h:119
virtual void CopyInformation(const DataObject *data)
virtual unsigned int GetNumberOfComponentsPerPixel() const
Vector< SpacingValueType, VImageDimension > SpacingType
Definition: itkImageBase.h:154
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:537
void CastFrom(const Vector< TCoordRepB, NVectorDimension > &pa)
Definition: itkVector.h:232
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:568
static unsigned int GetImageDimension()
Definition: itkImageBase.h:170
virtual void UpdateOutputData()
virtual void InitializeBufferedRegion(void)
SmartPointer< Self > Pointer
Definition: itkImageBase.h:118
SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:144
OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:140
void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
Definition: itkImageBase.h:728
static void TransformPhysicalPointToIndex(const MatrixType &matrix, const OriginType &origin, const DoublePoint &point, IndexType &index)
DataObject Superclass
Definition: itkImageBase.h:117
virtual void SetLargestPossibleRegion(const RegionType &region)
Base class for templated image classes.
Definition: itkImageBase.h:112
A templated class holding a point in n-Dimensional image space.
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, VImageDimension > &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:477
DirectionType m_PhysicalPointToIndex
Definition: itkImageBase.h:709
virtual const RegionType & GetRequestedRegion() const
Definition: itkImageBase.h:289
virtual void ComputeIndexToPhysicalPointMatrices()
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
Definition: itkImageBase.h:721
OffsetValueType ComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:331
double SpacePrecisionType
Definition: itkFloatTypes.h:30
Define additional traits for native types such as int or float.
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:147
ImageBase Self
Definition: itkImageBase.h:116
itk::OffsetValueType OffsetValueType
Definition: itkOffset.h:69
virtual void UpdateOutputInformation()
SpacePrecisionType PointValueType
Definition: itkImageBase.h:158
static const unsigned int ImageDimension
Definition: itkImageBase.h:131
Base class for all data objects in ITK.
IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:135
RegionType m_LargestPossibleRegion
Definition: itkImageBase.h:738
bool TransformPhysicalPointToIndex(const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
Definition: itkImageBase.h:412
virtual void SetRequestedRegion(const RegionType &region)
virtual void SetBufferedRegion(const RegionType &region)
virtual void SetOrigin(PointType _arg)