ITK  4.9.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 
128  typedef unsigned int ImageDimensionType;
129 
134  itkStaticConstMacro(ImageDimension, ImageDimensionType, VImageDimension);
135 
139 
144 
148 
151 
158 
163 
168 
170  virtual void Initialize() ITK_OVERRIDE;
171 
173  static unsigned int GetImageDimension()
174  { return VImageDimension; }
175 
180  itkSetMacro(Origin, PointType);
181  virtual void SetOrigin(const double origin[VImageDimension]);
182  virtual void SetOrigin(const float origin[VImageDimension]);
184 
210  virtual void SetDirection(const DirectionType & direction);
211 
215  itkGetConstReferenceMacro(Direction, DirectionType);
216 
220  itkGetConstReferenceMacro(InverseDirection, DirectionType);
221 
226  itkGetConstReferenceMacro(Spacing, SpacingType);
227 
232  itkGetConstReferenceMacro(Origin, PointType);
233 
241  virtual void Allocate(bool initialize=false);
242 
249  virtual void SetLargestPossibleRegion(const RegionType & region);
250 
257  virtual const RegionType & GetLargestPossibleRegion() const
258  { return m_LargestPossibleRegion; }
259 
263  virtual void SetBufferedRegion(const RegionType & region);
264 
268  virtual const RegionType & GetBufferedRegion() const
269  { return m_BufferedRegion; }
270 
278  virtual void SetRequestedRegion(const RegionType & region);
279 
287  virtual void SetRequestedRegion( const DataObject *data ) ITK_OVERRIDE;
288 
293  virtual const RegionType & GetRequestedRegion() const
294  { return m_RequestedRegion; }
295 
299  virtual void SetRegions(const RegionType& region)
300  {
301  this->SetLargestPossibleRegion(region);
302  this->SetBufferedRegion(region);
303  this->SetRequestedRegion(region);
304  }
306 
307  virtual void SetRegions(const SizeType& size)
308  {
309  RegionType region; region.SetSize(size);
310 
311  this->SetLargestPossibleRegion(region);
312  this->SetBufferedRegion(region);
313  this->SetRequestedRegion(region);
314  }
315 
326  const OffsetValueType * GetOffsetTable() const { return m_OffsetTable; }
328 
335  inline OffsetValueType ComputeOffset(const IndexType & ind) const
336  {
337  OffsetValueType offset = 0;
338 
340  ind,
342  offset);
343  return offset;
344  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
345  * Leaving here for documentation purposes
346  * OffsetValueType ComputeOffset(const IndexType & ind) const
347  * {
348  * // need to add bounds checking for the region/buffer?
349  * OffsetValueType offset = 0;
350  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
351  * // data is arranged as [][][][slice][row][col]
352  * // with Index[0] = col, Index[1] = row, Index[2] = slice
353  * for ( int i = VImageDimension - 1; i > 0; i-- )
354  * {
355  * offset += ( ind[i] - bufferedRegionIndex[i] ) * m_OffsetTable[i];
356  * }
357  * offset += ( ind[0] - bufferedRegionIndex[0] );
358  * return offset;
359  * }
360  */
361  }
362 
371  {
372  IndexType index;
373  const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
375 
377  offset,
379  index);
380  return index;
381  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
382  * Leaving here for documentation purposes
383  * IndexType ComputeIndex(OffsetValueType offset) const
384  * {
385  * IndexType index;
386  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
387  * for ( int i = VImageDimension - 1; i > 0; i-- )
388  * {
389  * index[i] = static_cast< IndexValueType >( offset / m_OffsetTable[i] );
390  * offset -= ( index[i] * m_OffsetTable[i] );
391  * index[i] += bufferedRegionIndex[i];
392  * }
393  * index[0] = bufferedRegionIndex[0] + static_cast< IndexValueType >( offset );
394  * return index;
395  * }
396  */
397 
398  }
399 
406  virtual void SetSpacing(const SpacingType & spacing);
407  virtual void SetSpacing(const double spacing[VImageDimension]);
408  virtual void SetSpacing(const float spacing[VImageDimension]);
410 
415  template< typename TCoordRep >
418  IndexType & index) const
419  {
422 
423  // Now, check to see if the index is within allowed bounds
424  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
425  return isInside;
426  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
427  * Leaving here for documentation purposes
428  * template< typename TCoordRep >
429  * bool TransformPhysicalPointToIndex(
430  * const Point< TCoordRep, VImageDimension > & point,
431  * IndexType & index) const
432  * {
433  * for ( unsigned int i = 0; i < VImageDimension; i++ )
434  * {
435  * TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
436  * for ( unsigned int j = 0; j < VImageDimension; j++ )
437  * {
438  * sum += this->m_PhysicalPointToIndex[i][j] * ( point[j] - this->m_Origin[j] );
439  * }
440  * index[i] = Math::RoundHalfIntegerUp< IndexValueType >(sum);
441  * }
442  * // Now, check to see if the index is within allowed bounds
443  * const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
444  * return isInside;
445  * }
446  */
447  }
448 
453  template< typename TCoordRep, typename TIndexRep >
457  {
459 
460  for ( unsigned int k = 0; k < VImageDimension; ++k )
461  {
462  cvector[k] = point[k] - this->m_Origin[k];
463  }
464  cvector = m_PhysicalPointToIndex * cvector;
465  for ( unsigned int i = 0; i < VImageDimension; ++i )
466  {
467  index[i] = static_cast< TIndexRep >( cvector[i] );
468  }
469 
470  // Now, check to see if the index is within allowed bounds
471  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
472 
473  return isInside;
474  }
475 
480  template< typename TCoordRep, typename TIndexRep >
484  {
485  for ( unsigned int r = 0; r < VImageDimension; ++r )
486  {
487  TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
488  for ( unsigned int c = 0; c < VImageDimension; ++c )
489  {
490  sum += this->m_IndexToPhysicalPoint(r, c) * index[c];
491  }
492  point[r] = sum + this->m_Origin[r];
493  }
494  }
496 
502  template< typename TCoordRep >
504  const IndexType & index,
506  {
509  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
510  * Leaving here for documentation purposes
511  * template< typename TCoordRep >
512  * void TransformIndexToPhysicalPoint(
513  * const IndexType & index,
514  * Point< TCoordRep, VImageDimension > & point) const
515  * {
516  * for ( unsigned int i = 0; i < VImageDimension; ++i )
517  * {
518  * point[i] = this->m_Origin[i];
519  * for ( unsigned int j = 0; j < VImageDimension; ++j )
520  * {
521  * point[i] += m_IndexToPhysicalPoint[i][j] * index[j];
522  * }
523  * }
524  * }
525  */
526  }
528 
540  template< typename TCoordRep >
542  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
543  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
544  {
545  //
546  //TODO: This temporary implementation should be replaced with Template
547  // MetaProgramming.
548  //
549  const DirectionType & direction = this->GetDirection();
550 
551  for ( unsigned int i = 0; i < VImageDimension; ++i )
552  {
553  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
554  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
555  for ( unsigned int j = 0; j < VImageDimension; ++j )
556  {
557  sum += direction[i][j] * inputGradient[j];
558  }
559  outputGradient[i] = static_cast< TCoordRep >( sum );
560  }
561  }
562 
571  template< typename TCoordRep >
573  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
574  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
575  {
576  //
577  //TODO: This temporary implementation should be replaced with Template
578  // MetaProgramming.
579  //
580  const DirectionType & inverseDirection = this->GetInverseDirection();
581 
582  for ( unsigned int i = 0; i < VImageDimension; ++i )
583  {
584  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
585  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
586  for ( unsigned int j = 0; j < VImageDimension; ++j )
587  {
588  sum += inverseDirection[i][j] * inputGradient[j];
589  }
590  outputGradient[i] = static_cast< TCoordRep >( sum );
591  }
592  }
593 
603  virtual void CopyInformation(const DataObject *data) ITK_OVERRIDE;
604 
615  virtual void Graft(const DataObject *data) ITK_OVERRIDE;
616 
624  virtual void UpdateOutputInformation() ITK_OVERRIDE;
625 
633  virtual void UpdateOutputData() ITK_OVERRIDE;
634 
638  virtual void SetRequestedRegionToLargestPossibleRegion() ITK_OVERRIDE;
639 
649  virtual bool RequestedRegionIsOutsideOfTheBufferedRegion() ITK_OVERRIDE;
650 
659  virtual bool VerifyRequestedRegion() ITK_OVERRIDE;
660 
679  virtual unsigned int GetNumberOfComponentsPerPixel() const;
680  virtual void SetNumberOfComponentsPerPixel(unsigned int);
682 
683 protected:
684  ImageBase();
685  ~ImageBase();
686  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
687 
692  void ComputeOffsetTable();
693 
700 
701 protected:
709 
714 
719  virtual void InitializeBufferedRegion();
720 
732  {
733  OffsetValueType offset = 0;
735  ind,
737  offset);
738  return offset;
739  }
741 
754  {
755  IndexType index;
756  const IndexType &bufferedRegionIndex = Self::GetBufferedRegion().GetIndex();
758  offset,
760  index);
761  return index;
762  }
764 
765 private:
766  ImageBase(const Self &) ITK_DELETE_FUNCTION;
767  void operator=(const Self &) ITK_DELETE_FUNCTION;
768 
769  void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
770  {
771  SpacingType s(spacing);
772  this->SetSpacing(s);
773  }
774 
775  template <typename TSpacingValue>
776  void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
777  {
779  SpacingType s;
780  s.CastFrom(sf);
781  this->SetSpacing(s);
782  }
783 
784  OffsetValueType m_OffsetTable[VImageDimension + 1];
785 
789 };
790 } // end namespace itk
791 
792 #ifndef ITK_MANUAL_INSTANTIATION
793 #include "itkImageBase.hxx"
794 #endif
795 
796 #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:299
virtual bool RequestedRegionIsOutsideOfTheBufferedRegion() override
OffsetValueType m_OffsetTable[VImageDimension+1]
Definition: itkImageBase.h:784
const OffsetValueType * GetOffsetTable() const
Definition: itkImageBase.h:326
Index< VImageDimension > IndexType
Definition: itkImageBase.h:137
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:788
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:257
IndexType FastComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:753
SpacingType m_Spacing
Definition: itkImageBase.h:705
An image region represents a structured region of data.
RegionType m_RequestedRegion
Definition: itkImageBase.h:787
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:156
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:503
virtual void UpdateOutputInformation() override
PointType m_Origin
Definition: itkImageBase.h:706
DirectionType m_IndexToPhysicalPoint
Definition: itkImageBase.h:712
static void ComputeIndex(const IndexType &bufferedRegionIndex, OffsetValueType offset, const OffsetValueType offsetTable[], IndexType &index)
virtual void Initialize() override
Size< VImageDimension > SizeType
Definition: itkImageBase.h:146
DirectionType m_InverseDirection
Definition: itkImageBase.h:708
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:162
virtual void SetSpacing(const SpacingType &spacing)
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Definition: itkImageBase.h:167
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
virtual void SetRegions(const SizeType &size)
Definition: itkImageBase.h:307
unsigned int ImageDimensionType
Definition: itkImageBase.h:125
IndexType ComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:370
DirectionType m_Direction
Definition: itkImageBase.h:707
::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:454
virtual void SetNumberOfComponentsPerPixel(unsigned int)
void ComputeOffsetTable()
virtual void Graft(const DataObject *data) override
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:268
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:142
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:157
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:541
OffsetValueType FastComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:731
void CastFrom(const Vector< TCoordRepB, NVectorDimension > &pa)
Definition: itkVector.h:243
static const ImageDimensionType ImageDimension
Definition: itkImageBase.h:134
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:572
static unsigned int GetImageDimension()
Definition: itkImageBase.h:173
SmartPointer< Self > Pointer
Definition: itkImageBase.h:118
SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:147
OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:143
void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
Definition: itkImageBase.h:776
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:481
DirectionType m_PhysicalPointToIndex
Definition: itkImageBase.h:713
virtual const DirectionType & GetDirection() const
virtual const RegionType & GetRequestedRegion() const
Definition: itkImageBase.h:293
virtual void ComputeIndexToPhysicalPointMatrices()
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
Definition: itkImageBase.h:769
OffsetValueType ComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:335
double SpacePrecisionType
Definition: itkFloatTypes.h:30
virtual void UpdateOutputData() override
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:150
ImageBase Self
Definition: itkImageBase.h:116
itk::OffsetValueType OffsetValueType
Definition: itkOffset.h:69
SpacePrecisionType PointValueType
Definition: itkImageBase.h:161
Base class for all data objects in ITK.
virtual void SetRequestedRegionToLargestPossibleRegion() override
IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:138
RegionType m_LargestPossibleRegion
Definition: itkImageBase.h:786
virtual bool VerifyRequestedRegion() override
bool TransformPhysicalPointToIndex(const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
Definition: itkImageBase.h:416
virtual void SetRequestedRegion(const RegionType &region)
virtual void SetBufferedRegion(const RegionType &region)
virtual void SetOrigin(PointType _arg)