ITK  5.0.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 #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_ASSIGN(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 Initialize() override;
166 
168  static unsigned int GetImageDimension()
169  { return VImageDimension; }
170 
175  itkSetMacro(Origin, PointType);
176  virtual void SetOrigin(const double origin[VImageDimension]);
177  virtual void SetOrigin(const float origin[VImageDimension]);
179 
205  virtual void SetDirection(const DirectionType & direction);
206 
210  itkGetConstReferenceMacro(Direction, DirectionType);
211 
215  itkGetConstReferenceMacro(InverseDirection, DirectionType);
216 
221  itkGetConstReferenceMacro(Spacing, SpacingType);
222 
227  itkGetConstReferenceMacro(Origin, PointType);
228 
236  virtual void Allocate(bool initialize=false);
237 
244  virtual void SetLargestPossibleRegion(const RegionType & region);
245 
252  virtual const RegionType & GetLargestPossibleRegion() const
253  { return m_LargestPossibleRegion; }
254 
258  virtual void SetBufferedRegion(const RegionType & region);
259 
263  virtual const RegionType & GetBufferedRegion() const
264  { return m_BufferedRegion; }
265 
273  virtual void SetRequestedRegion(const RegionType & region);
274 
282  void SetRequestedRegion( const DataObject *data ) override;
283 
288  virtual const RegionType & GetRequestedRegion() const
289  { return m_RequestedRegion; }
290 
294  virtual void SetRegions(const RegionType& region)
295  {
296  this->SetLargestPossibleRegion(region);
297  this->SetBufferedRegion(region);
298  this->SetRequestedRegion(region);
299  }
301 
302  virtual void SetRegions(const SizeType& size)
303  {
304  RegionType region; region.SetSize(size);
305 
306  this->SetLargestPossibleRegion(region);
307  this->SetBufferedRegion(region);
308  this->SetRequestedRegion(region);
309  }
310 
321  const OffsetValueType * GetOffsetTable() const { return m_OffsetTable; }
323 
330  inline OffsetValueType ComputeOffset(const IndexType & ind) const
331  {
332  OffsetValueType offset = 0;
333 
334  ImageHelper< VImageDimension, VImageDimension >::ComputeOffset(this->GetBufferedRegion().GetIndex(),
335  ind,
336  m_OffsetTable,
337  offset);
338  return offset;
339  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
340  * Leaving here for documentation purposes
341  * OffsetValueType ComputeOffset(const IndexType & ind) const
342  * {
343  * // need to add bounds checking for the region/buffer?
344  * OffsetValueType offset = 0;
345  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
346  * // data is arranged as [][][][slice][row][col]
347  * // with Index[0] = col, Index[1] = row, Index[2] = slice
348  * for ( int i = VImageDimension - 1; i > 0; i-- )
349  * {
350  * offset += ( ind[i] - bufferedRegionIndex[i] ) * m_OffsetTable[i];
351  * }
352  * offset += ( ind[0] - bufferedRegionIndex[0] );
353  * return offset;
354  * }
355  */
356  }
357 
366  {
367  IndexType index;
368  const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
370 
372  offset,
373  m_OffsetTable,
374  index);
375  return index;
376  /* NON TEMPLATE_META_PROGRAMMING_LOOP_UNROLLING data version
377  * Leaving here for documentation purposes
378  * IndexType ComputeIndex(OffsetValueType offset) const
379  * {
380  * IndexType index;
381  * const IndexType & bufferedRegionIndex = this->GetBufferedRegion().GetIndex();
382  * for ( int i = VImageDimension - 1; i > 0; i-- )
383  * {
384  * index[i] = static_cast< IndexValueType >( offset / m_OffsetTable[i] );
385  * offset -= ( index[i] * m_OffsetTable[i] );
386  * index[i] += bufferedRegionIndex[i];
387  * }
388  * index[0] = bufferedRegionIndex[0] + static_cast< IndexValueType >( offset );
389  * return index;
390  * }
391  */
392 
393  }
394 
401  virtual void SetSpacing(const SpacingType & spacing);
402  virtual void SetSpacing(const double spacing[VImageDimension]);
403  virtual void SetSpacing(const float spacing[VImageDimension]);
405 
410  template< typename TCoordRep >
413  IndexType & index) const
414  {
415  for ( unsigned int i = 0; i < VImageDimension; i++ )
416  {
417  TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
418  for ( unsigned int j = 0; j < VImageDimension; j++ )
419  {
420  sum += this->m_PhysicalPointToIndex[i][j] * ( point[j] - this->m_Origin[j] );
421  }
422  index[i] = Math::RoundHalfIntegerUp< IndexValueType >(sum);
423  }
424  // Now, check to see if the index is within allowed bounds
425  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
426  return isInside;
427  }
429 
434  template< typename TCoordRep, typename TIndexRep >
438  {
440 
441  for ( unsigned int k = 0; k < VImageDimension; ++k )
442  {
443  cvector[k] = point[k] - this->m_Origin[k];
444  }
445  cvector = m_PhysicalPointToIndex * cvector;
446  for ( unsigned int i = 0; i < VImageDimension; ++i )
447  {
448  index[i] = static_cast< TIndexRep >( cvector[i] );
449  }
450 
451  // Now, check to see if the index is within allowed bounds
452  const bool isInside = this->GetLargestPossibleRegion().IsInside(index);
453 
454  return isInside;
455  }
456 
461  template< typename TCoordRep, typename TIndexRep >
465  {
466  for ( unsigned int r = 0; r < VImageDimension; ++r )
467  {
468  TCoordRep sum = NumericTraits< TCoordRep >::ZeroValue();
469  for ( unsigned int c = 0; c < VImageDimension; ++c )
470  {
471  sum += this->m_IndexToPhysicalPoint(r, c) * index[c];
472  }
473  point[r] = sum + this->m_Origin[r];
474  }
475  }
477 
483  template< typename TCoordRep >
485  const IndexType & index,
487  {
488  for ( unsigned int i = 0; i < VImageDimension; ++i )
489  {
490  point[i] = this->m_Origin[i];
491  for ( unsigned int j = 0; j < VImageDimension; ++j )
492  {
493  point[i] += m_IndexToPhysicalPoint[i][j] * index[j];
494  }
495  }
496  }
498 
513  template< typename TCoordRep >
515  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
516  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
517  {
518  const DirectionType & direction = this->GetDirection();
519 
520  itkAssertInDebugAndIgnoreInReleaseMacro( inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
521 
522  for ( unsigned int i = 0; i < VImageDimension; ++i )
523  {
524  using CoordSumType = typename NumericTraits< TCoordRep >::AccumulateType;
525  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
526  for ( unsigned int j = 0; j < VImageDimension; ++j )
527  {
528  sum += direction[i][j] * inputGradient[j];
529  }
530  outputGradient[i] = static_cast< TCoordRep >( sum );
531  }
532  }
533 
546  template< typename TCoordRep >
548  const FixedArray< TCoordRep, VImageDimension > & inputGradient,
549  FixedArray< TCoordRep, VImageDimension > & outputGradient) const
550  {
551  const DirectionType & inverseDirection = this->GetInverseDirection();
552 
553  itkAssertInDebugAndIgnoreInReleaseMacro( inputGradient.GetDataPointer() != outputGradient.GetDataPointer());
554 
555  for ( unsigned int i = 0; i < VImageDimension; ++i )
556  {
557  using CoordSumType = typename NumericTraits< TCoordRep >::AccumulateType;
558  CoordSumType sum = NumericTraits< CoordSumType >::ZeroValue();
559  for ( unsigned int j = 0; j < VImageDimension; ++j )
560  {
561  sum += inverseDirection[i][j] * inputGradient[j];
562  }
563  outputGradient[i] = static_cast< TCoordRep >( sum );
564  }
565  }
566 
576  void CopyInformation(const DataObject *data) override;
577 
588  virtual void Graft(const Self *data);
589 
597  void UpdateOutputInformation() override;
598 
606  void UpdateOutputData() override;
607 
611  void SetRequestedRegionToLargestPossibleRegion() override;
612 
622  bool RequestedRegionIsOutsideOfTheBufferedRegion() override;
623 
632  bool VerifyRequestedRegion() override;
633 
652  virtual unsigned int GetNumberOfComponentsPerPixel() const;
653  virtual void SetNumberOfComponentsPerPixel(unsigned int);
655 
656 protected:
657  ImageBase();
658  ~ImageBase() override = default;
659  void PrintSelf(std::ostream & os, Indent indent) const override;
660 
665  void ComputeOffsetTable();
666 
672  virtual void ComputeIndexToPhysicalPointMatrices();
673 
674 protected:
682 
687 
692  virtual void InitializeBufferedRegion();
693 
705  {
706  OffsetValueType offset = 0;
707  ImageHelper<VImageDimension,VImageDimension>::ComputeOffset(Self::GetBufferedRegion().GetIndex(),
708  ind,
709  m_OffsetTable,
710  offset);
711  return offset;
712  }
714 
727  {
728  IndexType index;
729  const IndexType &bufferedRegionIndex = Self::GetBufferedRegion().GetIndex();
731  offset,
732  m_OffsetTable,
733  index);
734  return index;
735  }
737 
738  void Graft(const DataObject *data) override;
739 
740 private:
741  void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
742  {
743  SpacingType s(spacing);
744  this->SetSpacing(s);
745  }
746 
747  template <typename TSpacingValue>
748  void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
749  {
751  SpacingType s;
752  s.CastFrom(sf);
753  this->SetSpacing(s);
754  }
755 
756  OffsetValueType m_OffsetTable[VImageDimension + 1];
757 
761 };
762 } // end namespace itk
763 
764 #ifndef ITK_MANUAL_INSTANTIATION
765 #include "itkImageBase.hxx"
766 #endif
767 
768 #endif
void SetSize(const SizeType &size)
virtual void SetRegions(const RegionType &region)
Definition: itkImageBase.h:294
double SpacePrecisionType
Definition: itkFloatTypes.h:30
const OffsetValueType * GetOffsetTable() const
Definition: itkImageBase.h:321
static void ComputeOffset(const IndexType &bufferedRegionIndex, const IndexType &index, const OffsetValueType offsetTable[], OffsetValueType &offset)
Define numeric traits for std::vector.
RegionType m_BufferedRegion
Definition: itkImageBase.h:760
unsigned long SizeValueType
Definition: itkIntTypes.h:83
virtual const RegionType & GetLargestPossibleRegion() const
Definition: itkImageBase.h:252
IndexType FastComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:726
SpacingType m_Spacing
Definition: itkImageBase.h:678
An image region represents a structured region of data.
RegionType m_RequestedRegion
Definition: itkImageBase.h:759
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, VImageDimension > &point) const
Definition: itkImageBase.h:484
class ITK_FORWARD_EXPORT DataObject
Definition: itkDataObject.h:41
PointType m_Origin
Definition: itkImageBase.h:679
DirectionType m_IndexToPhysicalPoint
Definition: itkImageBase.h:685
static void ComputeIndex(const IndexType &bufferedRegionIndex, OffsetValueType offset, const OffsetValueType offsetTable[], IndexType &index)
DirectionType m_InverseDirection
Definition: itkImageBase.h:681
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:51
virtual void SetRegions(const SizeType &size)
Definition: itkImageBase.h:302
IndexType ComputeIndex(OffsetValueType offset) const
Definition: itkImageBase.h:365
DirectionType m_Direction
Definition: itkImageBase.h:680
bool TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, VImageDimension > &point, ContinuousIndex< TIndexRep, VImageDimension > &index) const
Get the continuous index from a physical point.
Definition: itkImageBase.h:435
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:263
signed long IndexValueType
Definition: itkIntTypes.h:90
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:514
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:68
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image...
Definition: itkOffset.h:67
OffsetValueType FastComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:704
void CastFrom(const Vector< TCoordRepB, NVectorDimension > &pa)
Definition: itkVector.h:237
typename SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:142
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, VImageDimension > &inputGradient, FixedArray< TCoordRep, VImageDimension > &outputGradient) const
Definition: itkImageBase.h:547
static unsigned int GetImageDimension()
Definition: itkImageBase.h:168
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
class ITK_TEMPLATE_EXPORT ImageBase
void InternalSetSpacing(const TSpacingValue spacing[VImageDimension])
Definition: itkImageBase.h:748
Base class for templated image classes.
Definition: itkImageBase.h:105
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:462
DirectionType m_PhysicalPointToIndex
Definition: itkImageBase.h:686
virtual const RegionType & GetRequestedRegion() const
Definition: itkImageBase.h:288
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void InternalSetSpacing(const SpacingValueType spacing[VImageDimension])
Definition: itkImageBase.h:741
OffsetValueType ComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:330
const IndexValueType * GetIndex() const
Definition: itkIndex.h:217
ValueType * GetDataPointer()
Base class for most ITK classes.
Definition: itkObject.h:60
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:133
signed long OffsetValueType
Definition: itkIntTypes.h:94
Base class for all data objects in ITK.
RegionType m_LargestPossibleRegion
Definition: itkImageBase.h:758
bool TransformPhysicalPointToIndex(const Point< TCoordRep, VImageDimension > &point, IndexType &index) const
Definition: itkImageBase.h:411