ITK  6.0.0
Insight Toolkit
itkImage.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 itkImage_h
19 #define itkImage_h
20 
21 #include "itkImageRegion.h"
25 #include "itkPoint.h"
26 #include "itkFixedArray.h"
27 #include "itkWeakPointer.h"
29 
30 #include <type_traits> // For is_same
31 
32 namespace itk
33 {
87 template <typename TPixel, unsigned int VImageDimension = 2>
88 class ITK_TEMPLATE_EXPORT Image : public ImageBase<VImageDimension>
89 {
90 public:
91  ITK_DISALLOW_COPY_AND_MOVE(Image);
92 
94  using Self = Image;
99 
101  itkNewMacro(Self);
102 
104  itkOverrideGetNameOfClassMacro(Image);
105 
108  using PixelType = TPixel;
109 
111  using ValueType = TPixel;
112 
117  using InternalPixelType = TPixel;
118 
120 
125 
129 
131  using typename Superclass::ImageDimensionType;
132 
134  using typename Superclass::IndexType;
135  using typename Superclass::IndexValueType;
136 
138  using typename Superclass::OffsetType;
139 
141  using typename Superclass::SizeType;
142  using typename Superclass::SizeValueType;
143 
146 
148  using typename Superclass::DirectionType;
149 
152  using typename Superclass::RegionType;
153 
156  using typename Superclass::SpacingType;
157  using typename Superclass::SpacingValueType;
158 
161  using typename Superclass::PointType;
162 
166 
168  using typename Superclass::OffsetValueType;
169 
176  template <typename UPixelType, unsigned int VUImageDimension = VImageDimension>
177  struct Rebind
178  {
180  };
181 
182 
183  template <typename UPixelType, unsigned int VUImageDimension = VImageDimension>
185 
186 
189  void
190  Allocate(bool initializePixels = false) override;
191 
194  void
195  Initialize() override;
196 
199  void
200  FillBuffer(const TPixel & value);
201 
207  void
208  SetPixel(const IndexType & index, const TPixel & value)
209  {
210  const OffsetValueType offset = this->FastComputeOffset(index);
211  (*m_Buffer)[offset] = value;
212  }
213 
218  const TPixel &
219  GetPixel(const IndexType & index) const
220  {
221  const OffsetValueType offset = this->FastComputeOffset(index);
222  return ((*m_Buffer)[offset]);
223  }
224 
229  TPixel &
230  GetPixel(const IndexType & index)
231  {
232  const OffsetValueType offset = this->FastComputeOffset(index);
233  return ((*m_Buffer)[offset]);
234  }
235 
240  TPixel &
241  operator[](const IndexType & index)
242  {
243  return this->GetPixel(index);
244  }
245 
250  const TPixel &
251  operator[](const IndexType & index) const
252  {
253  return this->GetPixel(index);
254  }
255 
258  virtual TPixel *
260  {
261  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
262  }
263  virtual const TPixel *
265  {
266  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
267  }
271  PixelContainer *
273  {
274  return m_Buffer.GetPointer();
275  }
276 
277  const PixelContainer *
279  {
280  return m_Buffer.GetPointer();
281  }
282 
285  void
286  SetPixelContainer(PixelContainer * container);
287 
298  virtual void
299  Graft(const Self * image);
300 
302  AccessorType
304  {
305  return AccessorType();
306  }
307 
309  const AccessorType
311  {
312  return AccessorType();
313  }
314 
316  NeighborhoodAccessorFunctorType
318  {
320  }
321 
323  const NeighborhoodAccessorFunctorType
325  {
327  }
328 
329  unsigned int
330  GetNumberOfComponentsPerPixel() const override;
331 
337  template <typename TEqualityComparable>
338  friend std::enable_if_t<std::is_same_v<TEqualityComparable, TPixel>, bool>
341  {
342  if ((lhs.GetBufferedRegion() != rhs.GetBufferedRegion()) || (lhs.m_Spacing != rhs.m_Spacing) ||
343  (lhs.m_Origin != rhs.m_Origin) || (lhs.m_Direction != rhs.m_Direction) ||
345  {
346  return false;
347  }
348 
349  if (lhs.m_Buffer == rhs.m_Buffer)
350  {
351  return true;
352  }
353 
354  if ((lhs.m_Buffer == nullptr) || (rhs.m_Buffer == nullptr))
355  {
356  return false;
357  }
358 
359  auto & lhsBuffer = *(lhs.m_Buffer);
360  auto & rhsBuffer = *(rhs.m_Buffer);
361 
362  const auto bufferSize = lhsBuffer.Size();
363 
364  if (bufferSize != rhsBuffer.Size())
365  {
366  return false;
367  }
368 
369  const TEqualityComparable * const lhsBufferPointer = lhsBuffer.GetBufferPointer();
370  const TEqualityComparable * const rhsBufferPointer = rhsBuffer.GetBufferPointer();
371 
372  return ((lhsBufferPointer == rhsBufferPointer) ||
373  std::equal(lhsBufferPointer, lhsBufferPointer + bufferSize, rhsBufferPointer));
374  }
375 
377  template <typename TEqualityComparable>
378  friend std::enable_if_t<std::is_same_v<TEqualityComparable, TPixel>, bool>
381  {
382  return !(lhs == rhs);
383  }
384 
385 protected:
386  Image() = default;
387  void
388  PrintSelf(std::ostream & os, Indent indent) const override;
389  void
390  Graft(const DataObject * data) override;
391 
392  ~Image() override = default;
393 
399  void
400  ComputeIndexToPhysicalPointMatrices() override;
401  using Superclass::Graft;
402 
403 private:
406 };
407 } // end namespace itk
408 
409 #ifndef ITK_MANUAL_INSTANTIATION
410 # include "itkImage.hxx"
411 #endif
412 
413 #endif
itk::ImageBase::OffsetValueType
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:147
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Image::GetNeighborhoodAccessor
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Definition: itkImage.h:317
itk::ImageBase::GetBufferedRegion
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:294
itk::Image::GetPixelAccessor
const AccessorType GetPixelAccessor() const
Definition: itkImage.h:310
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::Index< VImageDimension >
itk::Image< TNode *, VImageDimension >::IOPixelType
PixelType IOPixelType
Definition: itkImage.h:119
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itkWeakPointer.h
itk::Image::m_Buffer
PixelContainerPointer m_Buffer
Definition: itkImage.h:405
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:114
itkPoint.h
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::Image::GetPixel
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
Definition: itkImage.h:219
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::Image::operator!=
friend std::enable_if_t< std::is_same_v< TEqualityComparable, TPixel >, bool > operator!=(const Image< TEqualityComparable, VImageDimension > &lhs, const Image< TEqualityComparable, VImageDimension > &rhs)
Definition: itkImage.h:379
itkNeighborhoodAccessorFunctor.h
itk::SmartPointer< Self >
itk::DefaultPixelAccessor
Give access to partial aspects a type.
Definition: itkDefaultPixelAccessor.h:54
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::Image::GetBufferPointer
virtual const TPixel * GetBufferPointer() const
Definition: itkImage.h:264
itk::NeighborhoodAccessorFunctor
Provides accessor interfaces to Get pixels and is meant to be used on pointers contained within Neigh...
Definition: itkNeighborhoodAccessorFunctor.h:41
itkImageRegion.h
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:93
itk::Image::operator==
friend std::enable_if_t< std::is_same_v< TEqualityComparable, TPixel >, bool > operator==(const Image< TEqualityComparable, VImageDimension > &lhs, const Image< TEqualityComparable, VImageDimension > &rhs)
Definition: itkImage.h:339
itk::Image::GetPixelAccessor
AccessorType GetPixelAccessor()
Definition: itkImage.h:303
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itkDefaultPixelAccessor.h
itk::Image< TNode *, VImageDimension >::ValueType
TNode * ValueType
Definition: itkImage.h:111
itk::DefaultPixelAccessorFunctor
Provides a common API for pixel accessors for Image and VectorImage.
Definition: itkDefaultPixelAccessorFunctor.h:47
itk::Image::Rebind
Definition: itkImage.h:177
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::Image::SetPixel
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
Definition: itkImage.h:208
itkFixedArray.h
itk::Image< TNode *, VImageDimension >::InternalPixelType
TNode * InternalPixelType
Definition: itkImage.h:117
itk::Image::GetPixel
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
Definition: itkImage.h:230
itkImportImageContainer.h
itk::ImportImageContainer
Defines an itk::Image front-end to a standard C-array.
Definition: itkImportImageContainer.h:45
itk::OffsetValueType
long OffsetValueType
Definition: itkIntTypes.h:97
itk::Image< TNode *, VImageDimension >::PixelType
TNode * PixelType
Definition: itkImage.h:108
itk::WeakPointer
Implements a weak reference to an object.
Definition: itkWeakPointer.h:44
itk::ImageBase::m_Spacing
SpacingType m_Spacing
Definition: itkImageBase.h:867
itk::Image::GetPixelContainer
const PixelContainer * GetPixelContainer() const
Definition: itkImage.h:278
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::Image< TNode *, VImageDimension >::PixelContainerConstPointer
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition: itkImage.h:165
itk::Image::operator[]
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
Definition: itkImage.h:251
itk::Image::GetNeighborhoodAccessor
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
Definition: itkImage.h:324
itkDefaultPixelAccessorFunctor.h
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:61
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::Image::GetBufferPointer
virtual TPixel * GetBufferPointer()
Definition: itkImage.h:259
itk::Image::operator[]
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
Definition: itkImage.h:241
New
static Pointer New()
itk::Image< TNode *, VImageDimension >::PixelContainerPointer
typename PixelContainer::Pointer PixelContainerPointer
Definition: itkImage.h:164
AddImageFilter
Definition: itkAddImageFilter.h:81
itk::Image::GetPixelContainer
PixelContainer * GetPixelContainer()
Definition: itkImage.h:272
itk::ImageBase::m_Origin
PointType m_Origin
Definition: itkImageBase.h:868
itk::ImageBase::m_InverseDirection
DirectionType m_InverseDirection
Definition: itkImageBase.h:870
itk::ImageBase::m_Direction
DirectionType m_Direction
Definition: itkImageBase.h:869
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293