ITK  4.8.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  virtual void Initialize() ITK_OVERRIDE;
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 
238  virtual void Allocate(bool initialize=false);
239 
246  virtual void SetLargestPossibleRegion(const RegionType & region);
247 
254  virtual const RegionType & GetLargestPossibleRegion() const
255  { return m_LargestPossibleRegion; }
256 
260  virtual void SetBufferedRegion(const RegionType & region);
261 
265  virtual const RegionType & GetBufferedRegion() const
266  { return m_BufferedRegion; }
267 
275  virtual void SetRequestedRegion(const RegionType & region);
276 
284  virtual void SetRequestedRegion( const DataObject *data ) ITK_OVERRIDE;
285 
290  virtual const RegionType & GetRequestedRegion() const
291  { return m_RequestedRegion; }
292 
296  virtual void SetRegions(const RegionType& region)
297  {
298  this->SetLargestPossibleRegion(region);
299  this->SetBufferedRegion(region);
300  this->SetRequestedRegion(region);
301  }
303 
304  virtual void SetRegions(const SizeType& size)
305  {
306  RegionType region; region.SetSize(size);
307 
308  this->SetLargestPossibleRegion(region);
309  this->SetBufferedRegion(region);
310  this->SetRequestedRegion(region);
311  }
312 
323  const OffsetValueType * GetOffsetTable() const { return m_OffsetTable; }
325 
332  inline OffsetValueType ComputeOffset(const IndexType & ind) const
333  {
334  OffsetValueType offset = 0;
335 
337  ind,
339  offset);
340  return offset;
341  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
342  * Leaving here for documentation purposes
343  * OffsetValueType ComputeOffset(const IndexType & ind) const
344  * {
345  * // need to add bounds checking for the region/buffer?
346  * OffsetValueType offset = 0;
347  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
348  * // data is arranged as [][][][slice][row][col]
349  * // with Index[0] = col, Index[1] = row, Index[2] = slice
350  * for ( int i = VImageDimension - 1; i > 0; i-- )
351  * {
352  * offset += ( ind[i] - bufferedRegionIndex[i] ) * m_OffsetTable[i];
353  * }
354  * offset += ( ind[0] - bufferedRegionIndex[0] );
355  * return offset;
356  * }
357  */
358  }
359 
368  {
369  IndexType index;
370  const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
372 
374  offset,
376  index);
377  return index;
378  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
379  * Leaving here for documentation purposes
380  * IndexType ComputeIndex(OffsetValueType offset) const
381  * {
382  * IndexType index;
383  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
384  * for ( int i = VImageDimension - 1; i > 0; i-- )
385  * {
386  * index[i] = static_cast< IndexValueType >( offset / m_OffsetTable[i] );
387  * offset -= ( index[i] * m_OffsetTable[i] );
388  * index[i] += bufferedRegionIndex[i];
389  * }
390  * index[0] = bufferedRegionIndex[0] + static_cast< IndexValueType >( offset );
391  * return index;
392  * }
393  */
394 
395  }
396 
403  virtual void SetSpacing(const SpacingType & spacing);
404  virtual void SetSpacing(const double spacing[VImageDimension]);
405  virtual void SetSpacing(const float spacing[VImageDimension]);
407 
412  template< typename TCoordRep >
415  IndexType & index) const
416  {
419 
420  // Now, check to see if the index is within allowed bounds
421  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
422  return isInside;
423  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
424  * Leaving here for documentation purposes
425  * template< typename TCoordRep >
426  * bool TransformPhysicalPointToIndex(
427  * const Point< TCoordRep, VImageDimension > & point,
428  * IndexType & index) const
429  * {
430  * for ( unsigned int i = 0; i < VImageDimension; i++ )
431  * {
432  * TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
433  * for ( unsigned int j = 0; j < VImageDimension; j++ )
434  * {
435  * sum += this->m_PhysicalPointToIndex[i][j] * ( point[j] - this->m_Origin[j] );
436  * }
437  * index[i] = Math::RoundHalfIntegerUp< IndexValueType >(sum);
438  * }
439  * // Now, check to see if the index is within allowed bounds
440  * const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
441  * return isInside;
442  * }
443  */
444  }
445 
450  template< typename TCoordRep, typename TIndexRep >
454  {
456 
457  for ( unsigned int k = 0; k < VImageDimension; k++ )
458  {
459  cvector[k] = point[k] - this->m_Origin[k];
460  }
461  cvector = m_PhysicalPointToIndex * cvector;
462  for ( unsigned int i = 0; i < VImageDimension; i++ )
463  {
464  index[i] = static_cast< TIndexRep >( cvector[i] );
465  }
466 
467  // Now, check to see if the index is within allowed bounds
468  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
469 
470  return isInside;
471  }
472 
477  template< typename TCoordRep, typename TIndexRep >
481  {
482  for ( unsigned int r = 0; r < VImageDimension; r++ )
483  {
484  TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
485  for ( unsigned int c = 0; c < VImageDimension; c++ )
486  {
487  sum += this->m_IndexToPhysicalPoint(r, c) * index[c];
488  }
489  point[r] = sum + this->m_Origin[r];
490  }
491  }
493 
499  template< typename TCoordRep >
501  const IndexType & index,
503  {
506  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
507  * Leaving here for documentation purposes
508  * template< typename TCoordRep >
509  * void TransformIndexToPhysicalPoint(
510  * const IndexType & index,
511  * Point< TCoordRep, VImageDimension > & point) const
512  * {
513  * for ( unsigned int i = 0; i < VImageDimension; i++ )
514  * {
515  * point[i] = this->m_Origin[i];
516  * for ( unsigned int j = 0; j < VImageDimension; j++ )
517  * {
518  * point[i] += m_IndexToPhysicalPoint[i][j] * index[j];
519  * }
520  * }
521  * }
522  */
523  }
525 
537  template< typename TCoordRep >
539  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
540  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
541  {
542  //
543  //TODO: This temporary implementation should be replaced with Template
544  // MetaProgramming.
545  //
546  const DirectionType & direction = this->GetDirection();
547 
548  for ( unsigned int i = 0; i < VImageDimension; i++ )
549  {
550  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
551  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
552  for ( unsigned int j = 0; j < VImageDimension; j++ )
553  {
554  sum += direction[i][j] * inputGradient[j];
555  }
556  outputGradient[i] = static_cast< TCoordRep >( sum );
557  }
558  }
559 
568  template< typename TCoordRep >
570  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
571  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
572  {
573  //
574  //TODO: This temporary implementation should be replaced with Template
575  // MetaProgramming.
576  //
577  const DirectionType & inverseDirection = this->GetInverseDirection();
578 
579  for ( unsigned int i = 0; i < VImageDimension; i++ )
580  {
581  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
582  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
583  for ( unsigned int j = 0; j < VImageDimension; j++ )
584  {
585  sum += inverseDirection[i][j] * inputGradient[j];
586  }
587  outputGradient[i] = static_cast< TCoordRep >( sum );
588  }
589  }
590 
600  virtual void CopyInformation(const DataObject *data) ITK_OVERRIDE;
601 
612  virtual void Graft(const DataObject *data) ITK_OVERRIDE;
613 
621  virtual void UpdateOutputInformation() ITK_OVERRIDE;
622 
630  virtual void UpdateOutputData() ITK_OVERRIDE;
631 
635  virtual void SetRequestedRegionToLargestPossibleRegion() ITK_OVERRIDE;
636 
646  virtual bool RequestedRegionIsOutsideOfTheBufferedRegion() ITK_OVERRIDE;
647 
656  virtual bool VerifyRequestedRegion() ITK_OVERRIDE;
657 
676  virtual unsigned int GetNumberOfComponentsPerPixel() const;
677  virtual void SetNumberOfComponentsPerPixel(unsigned int);
679 
680 protected:
681  ImageBase();
682  ~ImageBase();
683  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
684 
689  void ComputeOffsetTable();
690 
697 
698 protected:
706 
711 
716  virtual void InitializeBufferedRegion();
717 
718 private:
719  ImageBase(const Self &); //purposely not implemented
720  void operator=(const Self &); //purposely not implemented
721 
722  void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
723  {
724  SpacingType s(spacing);
725  this->SetSpacing(s);
726  }
727 
728  template <typename TSpacingValue>
729  void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
730  {
732  SpacingType s;
733  s.CastFrom(sf);
734  this->SetSpacing(s);
735  }
736 
737  OffsetValueType m_OffsetTable[VImageDimension + 1];
738 
742 };
743 } // end namespace itk
744 
745 #ifndef ITK_MANUAL_INSTANTIATION
746 #include "itkImageBase.hxx"
747 #endif
748 
749 #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:296
virtual bool RequestedRegionIsOutsideOfTheBufferedRegion() override
OffsetValueType m_OffsetTable[VImageDimension+1]
Definition: itkImageBase.h:737
const OffsetValueType * GetOffsetTable() const
Definition: itkImageBase.h:323
Index< VImageDimension > IndexType
Definition: itkImageBase.h:134
const IndexType & GetIndex() const
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 CopyInformation(const DataObject *data) override
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
RegionType m_BufferedRegion
Definition: itkImageBase.h:741
virtual void SetDirection(const DirectionType &direction)
virtual void Allocate(bool initialize=false)
bool IsInside(const IndexType &index) const
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:254
SpacingType m_Spacing
Definition: itkImageBase.h:702
An image region represents a structured region of data.
RegionType m_RequestedRegion
Definition: itkImageBase.h:740
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:153
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:500
virtual void UpdateOutputInformation() override
PointType m_Origin
Definition: itkImageBase.h:703
DirectionType m_IndexToPhysicalPoint
Definition: itkImageBase.h:709
static void ComputeIndex(const IndexType &bufferedRegionIndex, OffsetValueType offset, const OffsetValueType offsetTable[], IndexType &index)
virtual void Initialize() override
Size< VImageDimension > SizeType
Definition: itkImageBase.h:143
DirectionType m_InverseDirection
Definition: itkImageBase.h:705
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:304
IndexType ComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:367
DirectionType m_Direction
Definition: itkImageBase.h:704
::itk::IndexValueType IndexValueType
Definition: itkIndex.h:79
virtual void InitializeBufferedRegion()
bool TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const
Get the continuous index from a physical point.
Definition: itkImageBase.h:451
virtual void SetNumberOfComponentsPerPixel(unsigned int)
void ComputeOffsetTable()
virtual void Graft(const DataObject *data) override
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:265
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:139
virtual const DirectionType & GetInverseDirection() const
SmartPointer< const Self > ConstPointer
Definition: itkImageBase.h:119
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:538
void CastFrom(const Vector< TCoordRepB, NVectorDimension > &pa)
Definition: itkVector.h:243
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:569
static unsigned int GetImageDimension()
Definition: itkImageBase.h:170
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:729
virtual void PrintSelf(std::ostream &os, Indent indent) const override
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:478
DirectionType m_PhysicalPointToIndex
Definition: itkImageBase.h:710
virtual const DirectionType & GetDirection() const
virtual const RegionType & GetRequestedRegion() const
Definition: itkImageBase.h:290
virtual void ComputeIndexToPhysicalPointMatrices()
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
Definition: itkImageBase.h:722
OffsetValueType ComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:332
double SpacePrecisionType
Definition: itkFloatTypes.h:30
virtual void UpdateOutputData() override
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:147
ImageBase Self
Definition: itkImageBase.h:116
itk::OffsetValueType OffsetValueType
Definition: itkOffset.h:69
SpacePrecisionType PointValueType
Definition: itkImageBase.h:158
static const unsigned int ImageDimension
Definition: itkImageBase.h:131
Base class for all data objects in ITK.
virtual void SetRequestedRegionToLargestPossibleRegion() override
IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:135
RegionType m_LargestPossibleRegion
Definition: itkImageBase.h:739
virtual bool VerifyRequestedRegion() override
bool TransformPhysicalPointToIndex(const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
Definition: itkImageBase.h:413
virtual void SetRequestedRegion(const RegionType &region)
virtual void SetBufferedRegion(const RegionType &region)
virtual void SetOrigin(PointType _arg)