ITK  6.0.0
Insight Toolkit
itkImageAdaptor.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  * https://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 #ifndef itkImageAdaptor_h
19 #define itkImageAdaptor_h
20 
21 #include "itkImage.h"
22 
23 namespace itk
24 {
25 
26 template <typename TPixelType, unsigned int VImageDimension>
27 class VectorImage;
28 
55 template <typename TImage, typename TAccessor>
56 class ITK_TEMPLATE_EXPORT ImageAdaptor : public ImageBase<TImage::ImageDimension>
57 {
58 public:
59  ITK_DISALLOW_COPY_AND_MOVE(ImageAdaptor);
60 
65  static constexpr unsigned int ImageDimension = TImage::ImageDimension;
66 
68  using Self = ImageAdaptor;
73 
75  itkOverrideGetNameOfClassMacro(ImageAdaptor);
76 
78  using InternalImageType = TImage;
79 
81  itkNewMacro(Self);
82 
85  using PixelType = typename TAccessor::ExternalType;
86 
89  using InternalPixelType = typename TAccessor::InternalType;
90 
92 
95  using AccessorType = TAccessor;
96 
99  using AccessorFunctorType = typename InternalImageType::AccessorFunctorType::template Rebind<Self>::Type;
100 
102  using typename Superclass::IndexType;
104 
106  using typename Superclass::SizeType;
108 
110  using typename Superclass::OffsetType;
112 
115  using typename Superclass::RegionType;
116 
119  using typename Superclass::SpacingType;
120 
123  using typename Superclass::PointType;
124 
128  using typename Superclass::DirectionType;
129 
130 
137  template <typename UPixelType, unsigned int UImageDimension = TImage::ImageDimension>
138  struct Rebind
139  {
141  };
142 
143  template <typename UPixelType, unsigned int VUImageDimension = TImage::ImageDimension>
145 
146 
153  void
154  SetLargestPossibleRegion(const RegionType & region) override;
155 
159  void
160  SetBufferedRegion(const RegionType & region) override;
161 
165  void
166  SetRequestedRegion(const RegionType & region) override;
167 
172  void
173  SetRequestedRegion(const DataObject * data) override;
174 
181  const RegionType &
182  GetRequestedRegion() const override;
183 
192  const RegionType &
193  GetLargestPossibleRegion() const override;
194 
200  const RegionType &
201  GetBufferedRegion() const override;
202 
204  void
205  Allocate(bool initialize = false) override;
206 
209  void
210  Initialize() override;
211 
213  void
214  SetPixel(const IndexType & index, const PixelType & value)
215  {
216  m_PixelAccessor.Set(m_Image->GetPixel(index), value);
217  }
218 
220  PixelType
221  GetPixel(const IndexType & index) const
222  {
223  return m_PixelAccessor.Get(m_Image->GetPixel(index));
224  }
225 
227  PixelType operator[](const IndexType & index) const { return m_PixelAccessor.Get(m_Image->GetPixel(index)); }
228 
230  const OffsetValueType *
231  GetOffsetTable() const;
232 
234  IndexType
235  ComputeIndex(OffsetValueType offset) const;
236 
239  using PixelContainer = typename TImage::PixelContainer;
240  using PixelContainerPointer = typename TImage::PixelContainerPointer;
241  using PixelContainerConstPointer = typename TImage::PixelContainerConstPointer;
242 
246  {
247  return m_Image->GetPixelContainer();
248  }
249 
250  const PixelContainer *
252  {
253  return m_Image->GetPixelContainer();
254  }
255 
258  void
259  SetPixelContainer(PixelContainer * container);
260 
271  virtual void
272  Graft(const Self * imgData);
273 
274 #ifndef ITK_FUTURE_LEGACY_REMOVE
275 
276  using InternalPixelPointerType [[deprecated("Please just use `InternalPixelType *` instead!")]] = InternalPixelType *;
277 #endif
278 
281  InternalPixelType *
282  GetBufferPointer();
283 
284  const InternalPixelType *
285  GetBufferPointer() const;
286 
288  void
289  SetSpacing(const SpacingType & spacing) override;
290 
291  void
292  SetSpacing(const double * spacing /*[ImageDimension]*/) override;
293 
294  void
295  SetSpacing(const float * spacing /*[ImageDimension]*/) override;
296 
300  const SpacingType &
301  GetSpacing() const override;
302 
306  const PointType &
307  GetOrigin() const override;
308 
310  void
311  SetOrigin(const PointType origin) override;
312 
313  void
314  SetOrigin(const double * origin /*[ImageDimension]*/) override;
315 
316  void
317  SetOrigin(const float * origin /*[ImageDimension]*/) override;
318 
320  void
321  SetDirection(const DirectionType & direction) override;
322 
326  const DirectionType &
327  GetDirection() const override;
328 
330  virtual void
331  SetImage(TImage *);
332 
334  void
335  Modified() const override;
336 
339  GetMTime() const override;
340 
342  AccessorType &
344  {
345  return m_PixelAccessor;
346  }
347 
349  const AccessorType &
351  {
352  return m_PixelAccessor;
353  }
354 
356  void
357  SetPixelAccessor(const AccessorType & accessor)
358  {
359  m_PixelAccessor = accessor;
360  }
361 
363  void
364  Update() override;
365 
366  void
367  CopyInformation(const DataObject * data) override;
368 
371  void
372  UpdateOutputInformation() override;
373 
374  void
375  SetRequestedRegionToLargestPossibleRegion() override;
376 
377  void
378  PropagateRequestedRegion() override;
379 
380  void
381  UpdateOutputData() override;
382 
383  bool
384  VerifyRequestedRegion() override;
385 
387  template <typename TIndexRep, typename TCoordRep>
390  {
391  return m_Image->template TransformPhysicalPointToContinuousIndex<TIndexRep>(point);
392  }
393 
402  template <typename TCoordRep>
403  ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
404  bool TransformPhysicalPointToContinuousIndex(const Point<TCoordRep, Self::ImageDimension> & point,
405  ContinuousIndex<TCoordRep, Self::ImageDimension> & index) const
406  {
407  return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
408  }
409 
411  template <typename TCoordRep>
412  [[nodiscard]] IndexType
414  {
415  return m_Image->TransformPhysicalPointToIndex(point);
416  }
417 
426  template <typename TCoordRep>
427  ITK_NODISCARD("Call the overload which has the point as the only parameter and returns the index")
428  bool TransformPhysicalPointToIndex(const Point<TCoordRep, Self::ImageDimension> & point, IndexType & index) const
429  {
430  return m_Image->TransformPhysicalPointToIndex(point, index);
431  }
432 
437  template <typename TCoordRep>
438  void
441  {
442  m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
443  }
444 
446  template <typename TCoordRep, typename TIndexRep>
449  {
450  return m_Image->template TransformContinuousIndexToPhysicalPoint<TIndexRep>(index);
451  }
452 
458  template <typename TCoordRep>
459  void
461  {
462  m_Image->TransformIndexToPhysicalPoint(index, point);
463  }
464 
466  template <typename TCoordRep>
469  {
470  return m_Image->template TransformIndexToPhysicalPoint<TCoordRep>(index);
471  }
472 
473  template <typename TCoordRep>
474  void
476  FixedArray<TCoordRep, Self::ImageDimension> & outputGradient) const
477  {
478  m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
479  }
480 
481  template <typename TVector>
482  [[nodiscard]] TVector
483  TransformLocalVectorToPhysicalVector(const TVector & inputGradient) const
484  {
485  TVector outputGradient;
486  TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
487  return outputGradient;
488  }
489 
490  template <typename TCoordRep>
491  void
493  FixedArray<TCoordRep, Self::ImageDimension> & outputGradient) const
494  {
495  m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
496  }
497 
498  template <typename TVector>
499  [[nodiscard]] TVector
500  TransformPhysicalVectorToLocalVector(const TVector & inputGradient) const
501  {
502  TVector outputGradient;
503  TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
504  return outputGradient;
505  }
506 
507 protected:
508  ImageAdaptor();
509  ~ImageAdaptor() override = default;
510  void
511  PrintSelf(std::ostream & os, Indent indent) const override;
512  void
513  Graft(const DataObject * data) override;
514  using Superclass::Graft;
515 
516 private:
517  // a specialized method to update PixelAccessors for VectorImages,
518  // to have the correct vector length of the image.
519  template <typename TPixelType>
520  void
522  {
523  this->m_PixelAccessor.SetVectorLength(this->m_Image->GetNumberOfComponentsPerPixel());
524  }
525 
526  // The other image types don't expect an accessor which needs any updates
527  template <typename T>
528  void
529  UpdateAccessor(T * itkNotUsed(dummy))
530  {}
531 
532  // Adapted image, most of the calls to ImageAdaptor
533  // will be delegated to this image
534  typename TImage::Pointer m_Image{};
535 
536  // Data accessor object,
537  // it converts the presentation of a pixel
538  AccessorType m_PixelAccessor{};
539 };
540 } // end namespace itk
541 
542 #ifndef ITK_MANUAL_INSTANTIATION
543 # include "itkImageAdaptor.hxx"
544 #endif
545 
546 #endif
itk::ImageAdaptor::GetPixel
PixelType GetPixel(const IndexType &index) const
Definition: itkImageAdaptor.h:221
itk::ImageBase< TImage::ImageDimension >::OffsetValueType
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:147
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::PixelContainerConstPointer
typename TImage::PixelContainerConstPointer PixelContainerConstPointer
Definition: itkImageAdaptor.h:241
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Index< VImageDimension >
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::InternalImageType
TImage InternalImageType
Definition: itkImageAdaptor.h:78
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::ModifiedTimeType
SizeValueType ModifiedTimeType
Definition: itkIntTypes.h:105
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::PixelContainer
typename TImage::PixelContainer PixelContainer
Definition: itkImageAdaptor.h:239
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:114
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:80
itk::ImageAdaptor::TransformLocalVectorToPhysicalVector
TVector TransformLocalVectorToPhysicalVector(const TVector &inputGradient) const
Definition: itkImageAdaptor.h:483
itk::ImageBase< TImage::ImageDimension >::SizeValueType
typename SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:151
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::InternalPixelType
typename Accessor::AddPixelAccessor< TImage::PixelType > ::InternalType InternalPixelType
Definition: itkImageAdaptor.h:89
itk::VectorImage
Templated n-dimensional vector image class.
Definition: itkImageAlgorithm.h:29
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::ImageAdaptor::UpdateAccessor
void UpdateAccessor(T *)
Definition: itkImageAdaptor.h:529
itk::ImageAdaptor::TransformLocalVectorToPhysicalVector
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
Definition: itkImageAdaptor.h:475
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageAdaptor::GetPixelContainer
const PixelContainer * GetPixelContainer() const
Definition: itkImageAdaptor.h:251
itk::ImageBase< TImage::ImageDimension >::IndexValueType
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:142
itk::ImageAdaptor::operator[]
PixelType operator[](const IndexType &index) const
Definition: itkImageAdaptor.h:227
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:93
itk::ImageAdaptor
Give access to partial aspects of voxels from an Image.
Definition: itkImageAdaptor.h:56
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageAdaptor::TransformPhysicalPointToContinuousIndex
ContinuousIndex< TIndexRep, TImage::ImageDimension > TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, TImage::ImageDimension > &point) const
Definition: itkImageAdaptor.h:389
itk::ImageAdaptor::TransformPhysicalPointToIndex
IndexType TransformPhysicalPointToIndex(const Point< TCoordRep, Self::ImageDimension > &point) const
Definition: itkImageAdaptor.h:413
itk::ImageAdaptor::UpdateAccessor
void UpdateAccessor(typename itk::VectorImage< TPixelType, ImageDimension > *)
Definition: itkImageAdaptor.h:521
itk::ImageAdaptor::Rebind
Definition: itkImageAdaptor.h:138
itk::ImageAdaptor::TransformPhysicalVectorToLocalVector
TVector TransformPhysicalVectorToLocalVector(const TVector &inputGradient) const
Definition: itkImageAdaptor.h:500
itk::point
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::AccessorFunctorType
typename InternalImageType::AccessorFunctorType::template Rebind< Self >::Type AccessorFunctorType
Definition: itkImageAdaptor.h:99
itk::ImageAdaptor::SetPixelAccessor
void SetPixelAccessor(const AccessorType &accessor)
Definition: itkImageAdaptor.h:357
itk::ImageAdaptor::TransformIndexToPhysicalPoint
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, Self::ImageDimension > &point) const
Definition: itkImageAdaptor.h:460
itk::OffsetValueType
long OffsetValueType
Definition: itkIntTypes.h:97
itk::FixedArray
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:53
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::PixelContainerPointer
typename TImage::PixelContainerPointer PixelContainerPointer
Definition: itkImageAdaptor.h:240
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::PixelType
typename Accessor::AddPixelAccessor< TImage::PixelType > ::ExternalType PixelType
Definition: itkImageAdaptor.h:85
itk::ImageAdaptor::TransformPhysicalVectorToLocalVector
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
Definition: itkImageAdaptor.h:492
itk::WeakPointer
Implements a weak reference to an object.
Definition: itkWeakPointer.h:44
itk::ImageAdaptor::TransformIndexToPhysicalPoint
Point< TCoordRep, Self::ImageDimension > TransformIndexToPhysicalPoint(const IndexType &index) const
Definition: itkImageAdaptor.h:468
itk::ImageAdaptor::TransformContinuousIndexToPhysicalPoint
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TCoordRep, Self::ImageDimension > &index, Point< TCoordRep, Self::ImageDimension > &point) const
Definition: itkImageAdaptor.h:439
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::ImageAdaptor::SetPixel
void SetPixel(const IndexType &index, const PixelType &value)
Definition: itkImageAdaptor.h:214
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::Accessor::AddPixelAccessor< TImage::PixelType >
AddImageFilter
Definition: itkAddImageFilter.h:81
itk::ImageAdaptor::GetPixelAccessor
AccessorType & GetPixelAccessor()
Definition: itkImageAdaptor.h:343
itk::ImageAdaptor::GetPixelAccessor
const AccessorType & GetPixelAccessor() const
Definition: itkImageAdaptor.h:350
itk::ImageAdaptor::GetPixelContainer
PixelContainerPointer GetPixelContainer()
Definition: itkImageAdaptor.h:245
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itk::ImageAdaptor< TImage, Accessor::AddPixelAccessor< TImage::PixelType > >::IOPixelType
PixelType IOPixelType
Definition: itkImageAdaptor.h:91
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::ImageAdaptor::TransformContinuousIndexToPhysicalPoint
Point< TCoordRep, TImage::ImageDimension > TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TIndexRep, Self::ImageDimension > &index) const
Definition: itkImageAdaptor.h:448