ITK  5.2.0
Insight Toolkit
itkImageBase.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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 #include "vnl/vnl_matrix_fixed.hxx" // Get the templates
43 
44 namespace itk
45 {
46 
104 template <unsigned int VImageDimension = 2>
105 class ITK_TEMPLATE_EXPORT ImageBase : public DataObject
106 {
107 public:
108  ITK_DISALLOW_COPY_AND_MOVE(ImageBase);
109 
111  using Self = ImageBase;
115 
117  itkNewMacro(Self);
118 
120  itkTypeMacro(ImageBase, DataObject);
121 
123  using ImageDimensionType = unsigned int;
124 
129  static constexpr ImageDimensionType ImageDimension = VImageDimension;
130 
134 
139 
143 
146 
153 
158 
163 
165  void
166  Initialize() override;
167 
169  static unsigned int
171  {
172  return VImageDimension;
173  }
174 
179  itkSetMacro(Origin, PointType);
180  virtual void
181  SetOrigin(const double origin[VImageDimension]);
182  virtual void
183  SetOrigin(const float origin[VImageDimension]);
185 
211  virtual void
212  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
244  Allocate(bool initialize = false);
245 
252  virtual void
253  SetLargestPossibleRegion(const RegionType & region);
254 
261  virtual const RegionType &
263  {
264  return m_LargestPossibleRegion;
265  }
266 
270  virtual void
271  SetBufferedRegion(const RegionType & region);
272 
276  virtual const RegionType &
278  {
279  return m_BufferedRegion;
280  }
281 
289  virtual void
290  SetRequestedRegion(const RegionType & region);
291 
299  void
300  SetRequestedRegion(const DataObject * data) override;
301 
306  virtual const RegionType &
308  {
309  return m_RequestedRegion;
310  }
311 
315  virtual void
316  SetRegions(const RegionType & region)
317  {
318  this->SetLargestPossibleRegion(region);
319  this->SetBufferedRegion(region);
320  this->SetRequestedRegion(region);
321  }
323 
324  virtual void
325  SetRegions(const SizeType & size)
326  {
327  RegionType region;
328  region.SetSize(size);
329 
330  this->SetLargestPossibleRegion(region);
331  this->SetBufferedRegion(region);
332  this->SetRequestedRegion(region);
333  }
334 
345  const OffsetValueType *
347  {
348  return m_OffsetTable;
349  }
351 
358  inline OffsetValueType
359  ComputeOffset(const IndexType & ind) const
360  {
361  OffsetValueType offset = 0;
362 
364  this->GetBufferedRegion().GetIndex(), ind, m_OffsetTable, offset);
365  return offset;
366  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
367  * Leaving here for documentation purposes
368  * OffsetValueType ComputeOffset(const IndexType & ind) const
369  * {
370  * // need to add bounds checking for the region/buffer?
371  * OffsetValueType offset = 0;
372  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
373  * // data is arranged as [][][][slice][row][col]
374  * // with Index[0] = col, Index[1] = row, Index[2] = slice
375  * for ( int i = VImageDimension - 1; i > 0; i-- )
376  * {
377  * offset += ( ind[i] - bufferedRegionIndex[i] ) * m_OffsetTable[i];
378  * }
379  * offset += ( ind[0] - bufferedRegionIndex[0] );
380  * return offset;
381  * }
382  */
383  }
384 
392  inline IndexType
394  {
395  IndexType index;
396  const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
398 
399  ImageHelper<VImageDimension, VImageDimension>::ComputeIndex(bufferedRegionIndex, offset, m_OffsetTable, index);
400  return index;
401  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
402  * Leaving here for documentation purposes
403  * IndexType ComputeIndex(OffsetValueType offset) const
404  * {
405  * IndexType index;
406  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
407  * for ( int i = VImageDimension - 1; i > 0; i-- )
408  * {
409  * index[i] = static_cast< IndexValueType >( offset / m_OffsetTable[i] );
410  * offset -= ( index[i] * m_OffsetTable[i] );
411  * index[i] += bufferedRegionIndex[i];
412  * }
413  * index[0] = bufferedRegionIndex[0] + static_cast< IndexValueType >( offset );
414  * return index;
415  * }
416  */
417  }
418 
425  virtual void
426  SetSpacing(const SpacingType & spacing);
427  virtual void
428  SetSpacing(const double spacing[VImageDimension]);
429  virtual void
430  SetSpacing(const float spacing[VImageDimension]);
432 
448  template <typename TCoordRep>
449  IndexType
451  {
452  IndexType index;
454 
455  for (unsigned int i = 0; i < VImageDimension; i++)
456  {
457  TCoordRep sum = NumericTraits<TCoordRep>::ZeroValue();
458  for (unsigned int j = 0; j < VImageDimension; j++)
459  {
460  sum += this->m_PhysicalPointToIndex[i][j] * (point[j] - this->m_Origin[j]);
461  }
462  index[i] = Math::RoundHalfIntegerUp<IndexValueType>(sum);
463  }
464  return index;
465  }
466 
471  template <typename TCoordRep>
472  bool
474  {
475  index = TransformPhysicalPointToIndex(point);
476 
477  // Now, check to see if the index is within allowed bounds
478  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
479  return isInside;
480  }
481 
482 
497  template <typename TIndexRep, typename TCoordRep>
500  {
503 
504  for (unsigned int k = 0; k < VImageDimension; ++k)
505  {
506  cvector[k] = point[k] - this->m_Origin[k];
507  }
508  cvector = m_PhysicalPointToIndex * cvector;
509  for (unsigned int i = 0; i < VImageDimension; ++i)
510  {
511  index[i] = static_cast<TIndexRep>(cvector[i]);
512  }
513  return index;
514  }
515 
520  template <typename TCoordRep, typename TIndexRep>
521  bool
524  {
525  index = TransformPhysicalPointToContinuousIndex<TIndexRep>(point);
526 
527  // Now, check to see if the index is within allowed bounds
528  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
529  return isInside;
530  }
531 
536  template <typename TCoordRep, typename TIndexRep>
537  void
539  Point<TCoordRep, VImageDimension> & point) const
540  {
541  for (unsigned int r = 0; r < VImageDimension; ++r)
542  {
543  TCoordRep sum = NumericTraits<TCoordRep>::ZeroValue();
544  for (unsigned int c = 0; c < VImageDimension; ++c)
545  {
546  sum += this->m_IndexToPhysicalPoint(r, c) * index[c];
547  }
548  point[r] = sum + this->m_Origin[r];
549  }
550  }
552 
557  template <typename TCoordRep, typename TIndexRep>
560  {
562  TransformContinuousIndexToPhysicalPoint(index, point);
563  return point;
564  }
566 
572  template <typename TCoordRep>
573  void
575  {
576  for (unsigned int i = 0; i < VImageDimension; ++i)
577  {
578  point[i] = this->m_Origin[i];
579  for (unsigned int j = 0; j < VImageDimension; ++j)
580  {
581  point[i] += m_IndexToPhysicalPoint[i][j] * index[j];
582  }
583  }
584  }
586 
592  template <typename TCoordRep>
595  {
597  TransformIndexToPhysicalPoint(index, point);
598  return point;
599  }
601 
616  template <typename TCoordRep>
617  void
619  FixedArray<TCoordRep, VImageDimension> & outputGradient) const
620  {
621  const DirectionType & direction = this->GetDirection();
622 
623  itkAssertInDebugAndIgnoreInReleaseMacro(inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
624 
625  for (unsigned int i = 0; i < VImageDimension; ++i)
626  {
627  using CoordSumType = typename NumericTraits<TCoordRep>::AccumulateType;
628  CoordSumType sum = NumericTraits<CoordSumType>::ZeroValue();
629  for (unsigned int j = 0; j < VImageDimension; ++j)
630  {
631  sum += direction[i][j] * inputGradient[j];
632  }
633  outputGradient[i] = static_cast<TCoordRep>(sum);
634  }
635  }
636 
643  template <typename TVector>
644  TVector
645  TransformLocalVectorToPhysicalVector(const TVector & inputGradient) const
646  {
647  TVector outputGradient;
648  TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
649  return outputGradient;
650  }
652 
653 
666  template <typename TCoordRep>
667  void
669  FixedArray<TCoordRep, VImageDimension> & outputGradient) const
670  {
671  const DirectionType & inverseDirection = this->GetInverseDirection();
672 
673  itkAssertInDebugAndIgnoreInReleaseMacro(inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
674 
675  for (unsigned int i = 0; i < VImageDimension; ++i)
676  {
677  using CoordSumType = typename NumericTraits<TCoordRep>::AccumulateType;
678  CoordSumType sum = NumericTraits<CoordSumType>::ZeroValue();
679  for (unsigned int j = 0; j < VImageDimension; ++j)
680  {
681  sum += inverseDirection[i][j] * inputGradient[j];
682  }
683  outputGradient[i] = static_cast<TCoordRep>(sum);
684  }
685  }
686 
693  template <typename TVector>
694  TVector
695  TransformPhysicalVectorToLocalVector(const TVector & inputGradient) const
696  {
697  TVector outputGradient;
698  TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
699  return outputGradient;
700  }
702 
712  void
713  CopyInformation(const DataObject * data) override;
714 
725  virtual void
726  Graft(const Self * image);
727 
735  void
736  UpdateOutputInformation() override;
737 
745  void
746  UpdateOutputData() override;
747 
751  void
752  SetRequestedRegionToLargestPossibleRegion() override;
753 
763  bool
764  RequestedRegionIsOutsideOfTheBufferedRegion() override;
765 
774  bool
775  VerifyRequestedRegion() override;
776 
795  virtual unsigned int
796  GetNumberOfComponentsPerPixel() const;
797  virtual void
798  SetNumberOfComponentsPerPixel(unsigned int);
800 
801 protected:
802  ImageBase();
803  ~ImageBase() override = default;
804  void
805  PrintSelf(std::ostream & os, Indent indent) const override;
806 
811  void
812  ComputeOffsetTable();
813 
819  virtual void
820  ComputeIndexToPhysicalPointMatrices();
821 
822 protected:
830 
835 
840  virtual void
841  InitializeBufferedRegion();
842 
854  FastComputeOffset(const IndexType & ind) const
855  {
856  OffsetValueType offset = 0;
858  Self::GetBufferedRegion().GetIndex(), ind, m_OffsetTable, offset);
859  return offset;
860  }
862 
874  IndexType
876  {
877  IndexType index;
878  const IndexType & bufferedRegionIndex = Self::GetBufferedRegion().GetIndex();
879  ImageHelper<VImageDimension, VImageDimension>::ComputeIndex(bufferedRegionIndex, offset, m_OffsetTable, index);
880  return index;
881  }
883 
884  void
885  Graft(const DataObject * data) override;
886 
887 private:
888  OffsetValueType m_OffsetTable[VImageDimension + 1]{};
889 
893 };
894 } // end namespace itk
895 
896 #ifndef ITK_MANUAL_INSTANTIATION
897 # include "itkImageBase.hxx"
898 #endif
899 
900 #endif
itk::ImageBase::TransformPhysicalPointToContinuousIndex
ContinuousIndex< TIndexRep, VImageDimension > TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, VImageDimension > &point) const
Returns the continuous index from a physical point.
Definition: itkImageBase.h:499
itk::ImageBase< TImage::ImageDimension >::OffsetValueType
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
itk::ImageBase::GetBufferedRegion
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:277
itk::ImageBase::GetImageDimension
static unsigned int GetImageDimension()
Definition: itkImageBase.h:170
itkObjectFactory.h
itk::ImageBase::TransformLocalVectorToPhysicalVector
TVector TransformLocalVectorToPhysicalVector(const TVector &inputGradient) const
Definition: itkImageBase.h:645
itk::Index< VImageDimension >
itk::ImageBase::SetRegions
virtual void SetRegions(const SizeType &size)
Definition: itkImageBase.h:325
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itkMatrix.h
itk::ImageBase::GetOffsetTable
const OffsetValueType * GetOffsetTable() const
Definition: itkImageBase.h:346
itkOffset.h
itkFloatTypes.h
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:105
itkImageHelper.h
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:69
itk::ImageBase< TImage::ImageDimension >::SizeValueType
typename SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:142
itk::ImageBase< TImage::ImageDimension >::PointValueType
SpacePrecisionType PointValueType
Definition: itkImageBase.h:156
itk::ImageBase::m_PhysicalPointToIndex
DirectionType m_PhysicalPointToIndex
Definition: itkImageBase.h:834
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::Vector< SpacingValueType, VImageDimension >
itk::NumericTraits::AccumulateType
double AccumulateType
Definition: itkNumericTraits.h:74
itk::ImageBase::TransformIndexToPhysicalPoint
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:574
itk::SmartPointer< Self >
itk::ImageBase
class ITK_TEMPLATE_EXPORT ImageBase
Definition: itkImageHelper.h:51
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itkDataObject.h
itk::ImageBase::m_LargestPossibleRegion
RegionType m_LargestPossibleRegion
Definition: itkImageBase.h:890
itk::ImageBase< TImage::ImageDimension >::IndexValueType
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:133
itk::FixedArray::GetDataPointer
ValueType * GetDataPointer()
Definition: itkFixedArray.h:313
itk::ImageBase::TransformContinuousIndexToPhysicalPoint
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, VImageDimension > &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:538
itkImageRegion.h
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageBase::FastComputeIndex
IndexType FastComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:875
itk::Index::GetIndex
const IndexValueType * GetIndex() const
Definition: itkIndex.h:228
itk::ImageBase::ComputeOffset
OffsetValueType ComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:359
itk::ImageBase::TransformContinuousIndexToPhysicalPoint
Point< TCoordRep, VImageDimension > TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, VImageDimension > &index) const
Definition: itkImageBase.h:559
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::SpacePrecisionType
double SpacePrecisionType
Definition: itkFloatTypes.h:30
itk::ImageBase::m_BufferedRegion
RegionType m_BufferedRegion
Definition: itkImageBase.h:892
itkFixedArray.h
itk::ImageBase::m_RequestedRegion
RegionType m_RequestedRegion
Definition: itkImageBase.h:891
itk::ImageBase::TransformPhysicalVectorToLocalVector
TVector TransformPhysicalVectorToLocalVector(const TVector &inputGradient) const
Definition: itkImageBase.h:695
itk::ImageBase::TransformPhysicalPointToContinuousIndex
bool TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const
Get the continuous index from a physical point.
Definition: itkImageBase.h:522
itk::ImageBase::SetRegions
virtual void SetRegions(const RegionType &region)
Definition: itkImageBase.h:316
itk::ImageBase::TransformPhysicalVectorToLocalVector
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:668
itk::DataObject
class ITK_FORWARD_EXPORT DataObject
Definition: itkDataObject.h:42
itk::FixedArray
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:52
itk::Matrix< SpacePrecisionType, VImageDimension, VImageDimension >
itk::ImageBase::GetLargestPossibleRegion
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:262
itk::Offset
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image.
Definition: itkOffset.h:67
itk::ImageBase::GetRequestedRegion
virtual const RegionType & GetRequestedRegion() const
Definition: itkImageBase.h:307
itk::ImageBase< TImage::ImageDimension >::SpacingValueType
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:151
itk::ImageBase::TransformPhysicalPointToIndex
IndexType TransformPhysicalPointToIndex(const Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:450
itk::NumericTraits::ZeroValue
static T ZeroValue()
Definition: itkNumericTraits.h:148
itk::ImageBase::m_Spacing
SpacingType m_Spacing
Definition: itkImageBase.h:826
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::OffsetValueType
signed long OffsetValueType
Definition: itkIntTypes.h:94
itk::ImageBase::TransformLocalVectorToPhysicalVector
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:618
itk::IndexValueType
signed long IndexValueType
Definition: itkIntTypes.h:90
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
itk::ImageBase::TransformIndexToPhysicalPoint
Point< TCoordRep, VImageDimension > TransformIndexToPhysicalPoint(const IndexType &index) const
Definition: itkImageBase.h:594
itk::ImageBase::FastComputeOffset
OffsetValueType FastComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:854
itk::Point< PointValueType, VImageDimension >
itk::ImageBase::ComputeIndex
IndexType ComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:393
itk::ImageBase::TransformPhysicalPointToIndex
bool TransformPhysicalPointToIndex(const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
Definition: itkImageBase.h:473
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:174
itk::ImageHelper::ComputeIndex
static void ComputeIndex(const IndexType &bufferedRegionIndex, OffsetValueType offset, const OffsetValueType offsetTable[], IndexType &index)
Definition: itkImageHelper.h:65
itk::ImageBase::m_IndexToPhysicalPoint
DirectionType m_IndexToPhysicalPoint
Definition: itkImageBase.h:833
itk::ImageBase::m_Origin
PointType m_Origin
Definition: itkImageBase.h:827
itk::ImageBase::m_InverseDirection
DirectionType m_InverseDirection
Definition: itkImageBase.h:829
itk::ImageBase::m_Direction
DirectionType m_Direction
Definition: itkImageBase.h:828
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::ImageHelper::ComputeOffset
static void ComputeOffset(const IndexType &bufferedRegionIndex, const IndexType &index, const OffsetValueType offsetTable[], OffsetValueType &offset)
Definition: itkImageHelper.h:103