ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkImageAdaptor.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 #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 > class VectorImage;
27 
54 template< typename TImage, typename TAccessor >
55 class ITK_TEMPLATE_EXPORT ImageAdaptor:public ImageBase< TImage::ImageDimension >
56 {
57 public:
58  ITK_DISALLOW_COPY_AND_ASSIGN(ImageAdaptor);
59 
64  static constexpr unsigned int ImageDimension = TImage::ImageDimension;
65 
67  using Self = ImageAdaptor;
72 
74  itkTypeMacro(ImageAdaptor, ImageBase);
75 
77  using InternalImageType = TImage;
78 
80  itkNewMacro(Self);
81 
84  using PixelType = typename TAccessor::ExternalType;
85 
88  using InternalPixelType = typename TAccessor::InternalType;
89 
91 
94  using AccessorType = TAccessor;
95 
98  using AccessorFunctorType = typename InternalImageType::AccessorFunctorType::template Rebind< Self >::Type;
99 
103 
105  using SizeType = typename Superclass::SizeType;
107 
109  using OffsetType = typename Superclass::OffsetType;
111 
115 
118  using SpacingType = typename Superclass::SpacingType;
119 
123 
128 
129 
136  template <typename UPixelType, unsigned int UImageDimension = TImage::ImageDimension>
137  struct Rebind
138  {
140  };
141 
142  template <typename UPixelType, unsigned int NUImageDimension = TImage::ImageDimension>
144 
145 
152  void SetLargestPossibleRegion(const RegionType & region) override;
153 
157  void SetBufferedRegion(const RegionType & region) override;
158 
162  void SetRequestedRegion(const RegionType & region) override;
163 
168  void SetRequestedRegion(const DataObject *data) override;
169 
176  const RegionType & GetRequestedRegion() const override;
177 
186  const RegionType & GetLargestPossibleRegion() const override;
187 
193  const RegionType & GetBufferedRegion() const override;
194 
196  void Allocate(bool initialize = false) override;
197 
200  void Initialize() override;
201 
203  void SetPixel(const IndexType & index, const PixelType & value)
204  { m_PixelAccessor.Set(m_Image->GetPixel(index), value); }
205 
207  PixelType GetPixel(const IndexType & index) const
208  { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
209 
211  PixelType operator[](const IndexType & index) const
212  { return m_PixelAccessor.Get( m_Image->GetPixel(index) ); }
213 
215  const OffsetValueType * GetOffsetTable() const;
216 
218  IndexType ComputeIndex(OffsetValueType offset) const;
219 
222  using PixelContainer = typename TImage::PixelContainer;
223  using PixelContainerPointer = typename TImage::PixelContainerPointer;
224  using PixelContainerConstPointer = typename TImage::PixelContainerConstPointer;
225 
228  { return m_Image->GetPixelContainer(); }
229 
231  { return m_Image->GetPixelContainer(); }
232 
235  void SetPixelContainer(PixelContainer *container);
236 
247  virtual void Graft(const Self *imgData);
248 
251 
254  InternalPixelType * GetBufferPointer();
255 
256  const InternalPixelType * GetBufferPointer() const;
257 
259  void SetSpacing(const SpacingType & values) override;
260 
261  void SetSpacing(const double *values /*[ImageDimension]*/) override;
262 
263  void SetSpacing(const float *values /*[ImageDimension]*/) override;
264 
268  const SpacingType & GetSpacing() const override;
269 
273  const PointType & GetOrigin() const override;
274 
276  void SetOrigin(const PointType values) override;
277 
278  void SetOrigin(const double *values /*[ImageDimension]*/) override;
279 
280  void SetOrigin(const float *values /*[ImageDimension]*/) override;
281 
283  void SetDirection(const DirectionType & direction) override;
284 
288  const DirectionType & GetDirection() const override;
289 
291  virtual void SetImage(TImage *);
292 
294  void Modified() const override;
295 
297  ModifiedTimeType GetMTime() const override;
298 
301  { return m_PixelAccessor; }
302 
305  { return m_PixelAccessor; }
306 
308  void SetPixelAccessor(const AccessorType & accessor)
309  { m_PixelAccessor = accessor; }
310 
312  void Update() override;
313 
314  void CopyInformation(const DataObject *data) override;
315 
318  void UpdateOutputInformation() override;
319 
320  void SetRequestedRegionToLargestPossibleRegion() override;
321 
322  void PropagateRequestedRegion() override;
323 
324  void UpdateOutputData() override;
325 
326  bool VerifyRequestedRegion() override;
327 
332  template< typename TCoordRep >
334  const Point< TCoordRep,
335  Self::ImageDimension > & point,
336  ContinuousIndex< TCoordRep,
337  Self::ImageDimension > & index) const
338  {
339  return m_Image->TransformPhysicalPointToContinuousIndex(point, index);
340  }
341 
346  template< typename TCoordRep >
348  const Point< TCoordRep,
349  Self::ImageDimension > & point,
350  IndexType & index) const
351  {
352  return m_Image->TransformPhysicalPointToIndex(point, index);
353  }
354 
359  template< typename TCoordRep >
361  const ContinuousIndex< TCoordRep,
362  Self::ImageDimension > & index,
363  Point< TCoordRep,
364  Self::ImageDimension > & point) const
365  {
366  m_Image->TransformContinuousIndexToPhysicalPoint(index, point);
367  }
368 
374  template< typename TCoordRep >
376  const IndexType & index,
377  Point< TCoordRep,
378  Self::ImageDimension > & point) const
379  {
380  m_Image->TransformIndexToPhysicalPoint(index, point);
381  }
382 
383  template< typename TCoordRep >
385  const FixedArray< TCoordRep, Self::ImageDimension > & inputGradient,
386  FixedArray< TCoordRep, Self::ImageDimension > & outputGradient) const
387  {
388  m_Image->TransformLocalVectorToPhysicalVector(inputGradient, outputGradient);
389  }
390 
391  template< typename TCoordRep >
393  const FixedArray< TCoordRep, Self::ImageDimension > & inputGradient,
394  FixedArray< TCoordRep, Self::ImageDimension > & outputGradient) const
395  {
396  m_Image->TransformPhysicalVectorToLocalVector(inputGradient, outputGradient);
397  }
398 
399 protected:
400 
401  ImageAdaptor();
402  ~ImageAdaptor() override = default;
403  void PrintSelf(std::ostream & os, Indent indent) const override;
404  void Graft(const DataObject *data) override;
405  using Superclass::Graft;
406 
407 private:
408 
409  // a specialized method to update PixelAccessors for VectorImages,
410  // to have the correct vector length of the image.
411  template< typename TPixelType >
412  void UpdateAccessor( typename ::itk::VectorImage< TPixelType, ImageDimension > * itkNotUsed( dummy ) )
413  {
414  this->m_PixelAccessor.SetVectorLength( this->m_Image->GetNumberOfComponentsPerPixel() );
415  }
416 
417  // The other image types don't expect an accessor which needs any updates
418  template< typename T > void UpdateAccessor( T *itkNotUsed( dummy ) ) { }
419 
420  // Adapted image, most of the calls to ImageAdaptor
421  // will be delegated to this image
422  typename TImage::Pointer m_Image;
423 
424  // Data accessor object,
425  // it converts the presentation of a pixel
427 };
428 } // end namespace itk
429 
430 #ifndef ITK_MANUAL_INSTANTIATION
431 #include "itkImageAdaptor.hxx"
432 #endif
433 
434 #endif
void UpdateAccessor(T *)
void TransformPhysicalVectorToLocalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
void TransformIndexToPhysicalPoint(const IndexType &index, Point< TCoordRep, Self::ImageDimension > &point) const
const PixelContainer * GetPixelContainer() const
unsigned long SizeValueType
Definition: itkIntTypes.h:83
An image region represents a structured region of data.
AccessorType m_PixelAccessor
Implements a weak reference to an object.
const AccessorType & GetPixelAccessor() const
bool TransformPhysicalPointToContinuousIndex(const Point< TCoordRep, Self::ImageDimension > &point, ContinuousIndex< TCoordRep, Self::ImageDimension > &index) const
Get the continuous index from a physical point.
void UpdateAccessor(typename::itk::VectorImage< TPixelType, ImageDimension > *)
PixelType GetPixel(const IndexType &index) const
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:51
PixelContainerPointer GetPixelContainer()
void SetPixelAccessor(const AccessorType &accessor)
typename InternalImageType::AccessorFunctorType::template Rebind< Self >::Type AccessorFunctorType
typename Accessor::AddPixelAccessor< TImage::PixelType >::ExternalType PixelType
signed long IndexValueType
Definition: itkIntTypes.h:90
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
void TransformContinuousIndexToPhysicalPoint(const ContinuousIndex< TCoordRep, Self::ImageDimension > &index, Point< TCoordRep, Self::ImageDimension > &point) const
typename SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:142
unsigned long ModifiedTimeType
Definition: itkIntTypes.h:104
void SetPixel(const IndexType &index, const PixelType &value)
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
PixelType operator[](const IndexType &index) const
bool TransformPhysicalPointToIndex(const Point< TCoordRep, Self::ImageDimension > &point, IndexType &index) const
Base class for templated image classes.
Definition: itkImageBase.h:105
A templated class holding a point in n-Dimensional image space.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
AccessorType & GetPixelAccessor()
Give access to partial aspects of voxels from an Image.
Base class for most ITK classes.
Definition: itkObject.h:60
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:133
typename Accessor::AddPixelAccessor< TImage::PixelType >::InternalType InternalPixelType
signed long OffsetValueType
Definition: itkIntTypes.h:94
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
void TransformLocalVectorToPhysicalVector(const FixedArray< TCoordRep, Self::ImageDimension > &inputGradient, FixedArray< TCoordRep, Self::ImageDimension > &outputGradient) const
TImage::Pointer m_Image