ITK  4.13.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 #include <vxl_version.h>
42 #if VXL_VERSION_DATE_FULL < 20160229
43 #include "vnl/vnl_matrix_fixed.txx" // Get the templates
44 #else
45 #include "vnl/vnl_matrix_fixed.hxx" // Get the templates
46 #endif
47 
49 
50 namespace itk
51 {
52 /* Forward declaration (ImageTransformHelper include's ImageBase) */
53 template< unsigned int NImageDimension, unsigned int R, unsigned int C, typename TPointValue, typename TMatrixValue >
54 class ITK_FORWARD_EXPORT ImageTransformHelper;
55 
113 template< unsigned int VImageDimension = 2 >
114 class ITK_TEMPLATE_EXPORT ImageBase:public DataObject
115 {
116 public:
118  typedef ImageBase Self;
122 
124  itkNewMacro(Self);
125 
127  itkTypeMacro(ImageBase, DataObject);
128 
130  typedef unsigned int ImageDimensionType;
131 
136  itkStaticConstMacro(ImageDimension, ImageDimensionType, VImageDimension);
137 
141 
146 
150 
153 
160 
165 
170 
172  virtual void Initialize() ITK_OVERRIDE;
173 
175  static unsigned int GetImageDimension()
176  { return VImageDimension; }
177 
182  itkSetMacro(Origin, PointType);
183  virtual void SetOrigin(const double origin[VImageDimension]);
184  virtual void SetOrigin(const float origin[VImageDimension]);
186 
212  virtual void SetDirection(const DirectionType & direction);
213 
217  itkGetConstReferenceMacro(Direction, DirectionType);
218 
222  itkGetConstReferenceMacro(InverseDirection, DirectionType);
223 
228  itkGetConstReferenceMacro(Spacing, SpacingType);
229 
234  itkGetConstReferenceMacro(Origin, PointType);
235 
243  virtual void Allocate(bool initialize=false);
244 
251  virtual void SetLargestPossibleRegion(const RegionType & region);
252 
259  virtual const RegionType & GetLargestPossibleRegion() const
260  { return m_LargestPossibleRegion; }
261 
265  virtual void SetBufferedRegion(const RegionType & region);
266 
270  virtual const RegionType & GetBufferedRegion() const
271  { return m_BufferedRegion; }
272 
280  virtual void SetRequestedRegion(const RegionType & region);
281 
289  virtual void SetRequestedRegion( const DataObject *data ) ITK_OVERRIDE;
290 
295  virtual const RegionType & GetRequestedRegion() const
296  { return m_RequestedRegion; }
297 
301  virtual void SetRegions(const RegionType& region)
302  {
303  this->SetLargestPossibleRegion(region);
304  this->SetBufferedRegion(region);
305  this->SetRequestedRegion(region);
306  }
308 
309  virtual void SetRegions(const SizeType& size)
310  {
311  RegionType region; region.SetSize(size);
312 
313  this->SetLargestPossibleRegion(region);
314  this->SetBufferedRegion(region);
315  this->SetRequestedRegion(region);
316  }
317 
328  const OffsetValueType * GetOffsetTable() const { return m_OffsetTable; }
330 
337  inline OffsetValueType ComputeOffset(const IndexType & ind) const
338  {
339  OffsetValueType offset = 0;
340 
341  ImageHelper< VImageDimension, VImageDimension >::ComputeOffset(this->GetBufferedRegion().GetIndex(),
342  ind,
343  m_OffsetTable,
344  offset);
345  return offset;
346  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
347  * Leaving here for documentation purposes
348  * OffsetValueType ComputeOffset(const IndexType & ind) const
349  * {
350  * // need to add bounds checking for the region/buffer?
351  * OffsetValueType offset = 0;
352  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
353  * // data is arranged as [][][][slice][row][col]
354  * // with Index[0] = col, Index[1] = row, Index[2] = slice
355  * for ( int i = VImageDimension - 1; i > 0; i-- )
356  * {
357  * offset += ( ind[i] - bufferedRegionIndex[i] ) * m_OffsetTable[i];
358  * }
359  * offset += ( ind[0] - bufferedRegionIndex[0] );
360  * return offset;
361  * }
362  */
363  }
364 
373  {
374  IndexType index;
375  const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
377 
379  offset,
380  m_OffsetTable,
381  index);
382  return index;
383  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
384  * Leaving here for documentation purposes
385  * IndexType ComputeIndex(OffsetValueType offset) const
386  * {
387  * IndexType index;
388  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
389  * for ( int i = VImageDimension - 1; i > 0; i-- )
390  * {
391  * index[i] = static_cast< IndexValueType >( offset / m_OffsetTable[i] );
392  * offset -= ( index[i] * m_OffsetTable[i] );
393  * index[i] += bufferedRegionIndex[i];
394  * }
395  * index[0] = bufferedRegionIndex[0] + static_cast< IndexValueType >( offset );
396  * return index;
397  * }
398  */
399 
400  }
401 
408  virtual void SetSpacing(const SpacingType & spacing);
409  virtual void SetSpacing(const double spacing[VImageDimension]);
410  virtual void SetSpacing(const float spacing[VImageDimension]);
412 
417  template< typename TCoordRep >
420  IndexType & index) const
421  {
423  ::TransformPhysicalPointToIndex(this->m_PhysicalPointToIndex, this->m_Origin, point, index);
424 
425  // Now, check to see if the index is within allowed bounds
426  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
427  return isInside;
428  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
429  * Leaving here for documentation purposes
430  * template< typename TCoordRep >
431  * bool TransformPhysicalPointToIndex(
432  * const Point< TCoordRep, VImageDimension > & point,
433  * IndexType & index) const
434  * {
435  * for ( unsigned int i = 0; i < VImageDimension; i++ )
436  * {
437  * TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
438  * for ( unsigned int j = 0; j < VImageDimension; j++ )
439  * {
440  * sum += this->m_PhysicalPointToIndex[i][j] * ( point[j] - this->m_Origin[j] );
441  * }
442  * index[i] = Math::RoundHalfIntegerUp< IndexValueType >(sum);
443  * }
444  * // Now, check to see if the index is within allowed bounds
445  * const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
446  * return isInside;
447  * }
448  */
449  }
450 
455  template< typename TCoordRep, typename TIndexRep >
459  {
461 
462  for ( unsigned int k = 0; k < VImageDimension; ++k )
463  {
464  cvector[k] = point[k] - this->m_Origin[k];
465  }
466  cvector = m_PhysicalPointToIndex * cvector;
467  for ( unsigned int i = 0; i < VImageDimension; ++i )
468  {
469  index[i] = static_cast< TIndexRep >( cvector[i] );
470  }
471 
472  // Now, check to see if the index is within allowed bounds
473  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
474 
475  return isInside;
476  }
477 
482  template< typename TCoordRep, typename TIndexRep >
486  {
487  for ( unsigned int r = 0; r < VImageDimension; ++r )
488  {
489  TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
490  for ( unsigned int c = 0; c < VImageDimension; ++c )
491  {
492  sum += this->m_IndexToPhysicalPoint(r, c) * index[c];
493  }
494  point[r] = sum + this->m_Origin[r];
495  }
496  }
498 
504  template< typename TCoordRep >
506  const IndexType & index,
508  {
510  TransformIndexToPhysicalPoint(this->m_IndexToPhysicalPoint, this->m_Origin, index, point);
511  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
512  * Leaving here for documentation purposes
513  * template< typename TCoordRep >
514  * void TransformIndexToPhysicalPoint(
515  * const IndexType & index,
516  * Point< TCoordRep, VImageDimension > & point) const
517  * {
518  * for ( unsigned int i = 0; i < VImageDimension; ++i )
519  * {
520  * point[i] = this->m_Origin[i];
521  * for ( unsigned int j = 0; j < VImageDimension; ++j )
522  * {
523  * point[i] += m_IndexToPhysicalPoint[i][j] * index[j];
524  * }
525  * }
526  * }
527  */
528  }
530 
545  template< typename TCoordRep >
547  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
548  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
549  {
550  //
551  //TODO: This temporary implementation should be replaced with Template
552  // MetaProgramming.
553  //
554  const DirectionType & direction = this->GetDirection();
555 
556  itkAssertInDebugAndIgnoreInReleaseMacro( inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
557 
558  for ( unsigned int i = 0; i < VImageDimension; ++i )
559  {
560  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
561  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
562  for ( unsigned int j = 0; j < VImageDimension; ++j )
563  {
564  sum += direction[i][j] * inputGradient[j];
565  }
566  outputGradient[i] = static_cast< TCoordRep >( sum );
567  }
568  }
569 
582  template< typename TCoordRep >
584  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
585  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
586  {
587  //
588  //TODO: This temporary implementation should be replaced with Template
589  // MetaProgramming.
590  //
591  const DirectionType & inverseDirection = this->GetInverseDirection();
592 
593  itkAssertInDebugAndIgnoreInReleaseMacro( inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
594 
595  for ( unsigned int i = 0; i < VImageDimension; ++i )
596  {
597  typedef typename NumericTraits< TCoordRep >::AccumulateType CoordSumType;
598  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
599  for ( unsigned int j = 0; j < VImageDimension; ++j )
600  {
601  sum += inverseDirection[i][j] * inputGradient[j];
602  }
603  outputGradient[i] = static_cast< TCoordRep >( sum );
604  }
605  }
606 
616  virtual void CopyInformation(const DataObject *data) ITK_OVERRIDE;
617 
628  virtual void Graft(const Self *data);
629 
637  virtual void UpdateOutputInformation() ITK_OVERRIDE;
638 
646  virtual void UpdateOutputData() ITK_OVERRIDE;
647 
651  virtual void SetRequestedRegionToLargestPossibleRegion() ITK_OVERRIDE;
652 
662  virtual bool RequestedRegionIsOutsideOfTheBufferedRegion() ITK_OVERRIDE;
663 
672  virtual bool VerifyRequestedRegion() ITK_OVERRIDE;
673 
692  virtual unsigned int GetNumberOfComponentsPerPixel() const;
693  virtual void SetNumberOfComponentsPerPixel(unsigned int);
695 
696 protected:
697  ImageBase();
698  ~ImageBase() ITK_OVERRIDE;
699  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
700 
705  void ComputeOffsetTable();
706 
712  virtual void ComputeIndexToPhysicalPointMatrices();
713 
714 protected:
718  SpacingType m_Spacing;
719  PointType m_Origin;
720  DirectionType m_Direction;
721  DirectionType m_InverseDirection;
722 
725  DirectionType m_IndexToPhysicalPoint;
726  DirectionType m_PhysicalPointToIndex;
727 
732  virtual void InitializeBufferedRegion();
733 
744  OffsetValueType FastComputeOffset(const IndexType &ind) const
745  {
746  OffsetValueType offset = 0;
747  ImageHelper<VImageDimension,VImageDimension>::ComputeOffset(Self::GetBufferedRegion().GetIndex(),
748  ind,
749  m_OffsetTable,
750  offset);
751  return offset;
752  }
754 
767  {
768  IndexType index;
769  const IndexType &bufferedRegionIndex = Self::GetBufferedRegion().GetIndex();
771  offset,
772  m_OffsetTable,
773  index);
774  return index;
775  }
777 
778  virtual void Graft(const DataObject *data) ITK_OVERRIDE;
779 
780 private:
781  ITK_DISALLOW_COPY_AND_ASSIGN(ImageBase);
782 
783  void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
784  {
785  SpacingType s(spacing);
786  this->SetSpacing(s);
787  }
788 
789  template <typename TSpacingValue>
790  void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
791  {
793  SpacingType s;
794  s.CastFrom(sf);
795  this->SetSpacing(s);
796  }
797 
798  OffsetValueType m_OffsetTable[VImageDimension + 1];
799 
803 };
804 } // end namespace itk
805 
806 #ifndef ITK_MANUAL_INSTANTIATION
807 #include "itkImageBase.hxx"
808 #endif
809 
810 #endif
void SetSize(const SizeType &size)
const IndexValueType * GetIndex() const
Definition: itkIndex.h:255
static void TransformIndexToPhysicalPoint(const MatrixType &matrix, const OriginType &origin, const IndexType &index, DoublePoint &point)
virtual void SetRegions(const RegionType &region)
Definition: itkImageBase.h:301
const OffsetValueType * GetOffsetTable() const
Definition: itkImageBase.h:328
Index< VImageDimension > IndexType
Definition: itkImageBase.h:139
Represent the offset between two n-dimensional indexes in a n-dimensional image.
Definition: itkOffset.h:56
Fast index/physical index computation.
static void ComputeOffset(const IndexType &bufferedRegionIndex, const IndexType &index, const OffsetValueType offsetTable[], OffsetValueType &offset)
itk::SizeValueType SizeValueType
Definition: itkSize.h:60
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
RegionType m_BufferedRegion
Definition: itkImageBase.h:802
signed long OffsetValueType
Definition: itkIntTypes.h:154
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:259
IndexType FastComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:766
An image region represents a structured region of data.
RegionType m_RequestedRegion
Definition: itkImageBase.h:801
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:158
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:505
static void ComputeIndex(const IndexType &bufferedRegionIndex, OffsetValueType offset, const OffsetValueType offsetTable[], IndexType &index)
Size< VImageDimension > SizeType
Definition: itkImageBase.h:148
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:164
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Definition: itkImageBase.h:169
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:50
virtual void SetRegions(const SizeType &size)
Definition: itkImageBase.h:309
unsigned int ImageDimensionType
Definition: itkImageBase.h:127
IndexType ComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:372
::itk::IndexValueType IndexValueType
Definition: itkIndex.h:80
bool TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const
Get the continuous index from a physical point.
Definition: itkImageBase.h:456
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:270
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:144
SmartPointer< const Self > ConstPointer
Definition: itkImageBase.h:121
Vector< SpacingValueType, VImageDimension > SpacingType
Definition: itkImageBase.h:159
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:546
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:583
SmartPointer< Self > Pointer
Definition: itkImageBase.h:120
SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:149
OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:145
void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
Definition: itkImageBase.h:790
static void TransformPhysicalPointToIndex(const MatrixType &matrix, const OriginType &origin, const DoublePoint &point, IndexType &index)
DataObject Superclass
Definition: itkImageBase.h:119
Base class for templated image classes.
Definition: itkImageBase.h:114
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:483
virtual const RegionType & GetRequestedRegion() const
Definition: itkImageBase.h:295
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
Definition: itkImageBase.h:783
OffsetValueType ComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:337
double SpacePrecisionType
Definition: itkFloatTypes.h:30
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:152
ImageBase Self
Definition: itkImageBase.h:118
ValueType * GetDataPointer()
SpacePrecisionType PointValueType
Definition: itkImageBase.h:163
Base class for all data objects in ITK.
IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:140
RegionType m_LargestPossibleRegion
Definition: itkImageBase.h:800
bool TransformPhysicalPointToIndex(const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
Definition: itkImageBase.h:418