ITK  5.2.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  * 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 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 namespace itk
31 {
85 template <typename TPixel, unsigned int VImageDimension = 2>
86 class ITK_TEMPLATE_EXPORT Image : public ImageBase<VImageDimension>
87 {
88 public:
89  ITK_DISALLOW_COPY_AND_MOVE(Image);
90 
92  using Self = Image;
97 
99  itkNewMacro(Self);
100 
102  itkTypeMacro(Image, ImageBase);
103 
106  using PixelType = TPixel;
107 
109  using ValueType = TPixel;
110 
115  using InternalPixelType = TPixel;
116 
118 
123 
127 
129  using ImageDimensionType = typename Superclass::ImageDimensionType;
130 
134 
136  using OffsetType = typename Superclass::OffsetType;
137 
139  using SizeType = typename Superclass::SizeType;
141 
144 
147 
151 
154  using SpacingType = typename Superclass::SpacingType;
155  using SpacingValueType = typename Superclass::SpacingValueType;
156 
160 
164 
167 
174  template <typename UPixelType, unsigned int NUImageDimension = VImageDimension>
175  struct Rebind
176  {
178  };
179 
180 
181  template <typename UPixelType, unsigned int NUImageDimension = VImageDimension>
183 
184 
187  void
188  Allocate(bool initializePixels = false) override;
189 
192  void
193  Initialize() override;
194 
197  void
198  FillBuffer(const TPixel & value);
199 
205  void
206  SetPixel(const IndexType & index, const TPixel & value)
207  {
208  OffsetValueType offset = this->FastComputeOffset(index);
209  (*m_Buffer)[offset] = value;
210  }
211 
216  const TPixel &
217  GetPixel(const IndexType & index) const
218  {
219  OffsetValueType offset = this->FastComputeOffset(index);
220  return ((*m_Buffer)[offset]);
221  }
222 
227  TPixel &
228  GetPixel(const IndexType & index)
229  {
230  OffsetValueType offset = this->FastComputeOffset(index);
231  return ((*m_Buffer)[offset]);
232  }
233 
238  TPixel & operator[](const IndexType & index) { return this->GetPixel(index); }
239 
244  const TPixel & operator[](const IndexType & index) const { return this->GetPixel(index); }
245 
248  virtual TPixel *
250  {
251  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
252  }
253  virtual const TPixel *
255  {
256  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
257  }
259 
261  PixelContainer *
263  {
264  return m_Buffer.GetPointer();
265  }
266 
267  const PixelContainer *
269  {
270  return m_Buffer.GetPointer();
271  }
272 
275  void
276  SetPixelContainer(PixelContainer * container);
277 
288  virtual void
289  Graft(const Self * image);
290 
292  AccessorType
294  {
295  return AccessorType();
296  }
297 
299  const AccessorType
301  {
302  return AccessorType();
303  }
304 
306  NeighborhoodAccessorFunctorType
308  {
310  }
311 
313  const NeighborhoodAccessorFunctorType
315  {
317  }
318 
319  unsigned int
320  GetNumberOfComponentsPerPixel() const override;
321 
323  friend bool
324  operator==(const Image & lhs, const Image & rhs)
325  {
326  if ((lhs.GetBufferedRegion() != rhs.GetBufferedRegion()) || (lhs.m_Spacing != rhs.m_Spacing) ||
327  (lhs.m_Origin != rhs.m_Origin) || (lhs.m_Direction != rhs.m_Direction) ||
329  {
330  return false;
331  }
332 
333  if (lhs.m_Buffer == rhs.m_Buffer)
334  {
335  return true;
336  }
337 
338  if ((lhs.m_Buffer == nullptr) || (rhs.m_Buffer == nullptr))
339  {
340  return false;
341  }
342 
343  auto & lhsBuffer = *(lhs.m_Buffer);
344  auto & rhsBuffer = *(rhs.m_Buffer);
345 
346  const auto bufferSize = lhsBuffer.Size();
347 
348  if (bufferSize != rhsBuffer.Size())
349  {
350  return false;
351  }
352 
353  const TPixel * const lhsBufferPointer = lhsBuffer.GetBufferPointer();
354  const TPixel * const rhsBufferPointer = rhsBuffer.GetBufferPointer();
355 
356  return ((lhsBufferPointer == rhsBufferPointer) ||
357  std::equal(lhsBufferPointer, lhsBufferPointer + bufferSize, rhsBufferPointer));
358  }
359 
361  friend bool
362  operator!=(const Image & lhs, const Image & rhs)
363  {
364  return !(lhs == rhs);
365  }
366 
367 protected:
368  Image();
369  void
370  PrintSelf(std::ostream & os, Indent indent) const override;
371  void
372  Graft(const DataObject * data) override;
373 
374  ~Image() override = default;
375 
381  void
382  ComputeIndexToPhysicalPointMatrices() override;
383  using Superclass::Graft;
384 
385 private:
388 };
389 } // end namespace itk
390 
391 #ifndef ITK_MANUAL_INSTANTIATION
392 # include "itkImage.hxx"
393 #endif
394 
395 #endif
itk::ImageBase::OffsetValueType
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
itk::Image::GetNeighborhoodAccessor
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Definition: itkImage.h:307
itk::Image::operator==
friend bool operator==(const Image &lhs, const Image &rhs)
Definition: itkImage.h:324
itk::ImageBase::GetBufferedRegion
virtual const RegionType & GetBufferedRegion() const
Definition: itkImageBase.h:277
itk::Image::GetPixelAccessor
const AccessorType GetPixelAccessor() const
Definition: itkImage.h:300
itk::Index< VImageDimension >
itk::Image< TNode *, VImageDimension >::IOPixelType
PixelType IOPixelType
Definition: itkImage.h:117
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itkWeakPointer.h
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itk::Image::m_Buffer
PixelContainerPointer m_Buffer
Definition: itkImage.h:387
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:105
itkPoint.h
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:69
itk::ImageBase::SizeValueType
typename SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:142
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::Vector< SpacingValueType, VImageDimension >
itk::Image::GetPixel
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
Definition: itkImage.h:217
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
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:254
itk::ImageBase::IndexValueType
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:133
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::Image::GetPixelAccessor
AccessorType GetPixelAccessor()
Definition: itkImage.h:293
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itkDefaultPixelAccessor.h
itk::Image< TNode *, VImageDimension >::ValueType
TNode * ValueType
Definition: itkImage.h:109
itk::DefaultPixelAccessorFunctor
Provides a common API for pixel accessors for Image and VectorImage.
Definition: itkDefaultPixelAccessorFunctor.h:47
itk::ImageBase::ImageDimensionType
unsigned int ImageDimensionType
Definition: itkImageBase.h:123
itk::Image::Rebind
Definition: itkImage.h:175
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:206
itkFixedArray.h
itk::Image< TNode *, VImageDimension >::InternalPixelType
TNode * InternalPixelType
Definition: itkImage.h:115
itk::Image::GetPixel
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
Definition: itkImage.h:228
itkImportImageContainer.h
itk::ImportImageContainer
Defines an itk::Image front-end to a standard C-array.
Definition: itkImportImageContainer.h:45
itk::Matrix< SpacePrecisionType, VImageDimension, VImageDimension >
itk::Image< TNode *, VImageDimension >::PixelType
TNode * PixelType
Definition: itkImage.h:106
itk::Offset
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image.
Definition: itkOffset.h:67
itk::ImageBase::SpacingValueType
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:151
itk::WeakPointer
Implements a weak reference to an object.
Definition: itkWeakPointer.h:44
itk::ImageBase::m_Spacing
SpacingType m_Spacing
Definition: itkImageBase.h:826
itk::Image::GetPixelContainer
const PixelContainer * GetPixelContainer() const
Definition: itkImage.h:268
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::Image< TNode *, VImageDimension >::PixelContainerConstPointer
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition: itkImage.h:163
itk::OffsetValueType
signed long OffsetValueType
Definition: itkIntTypes.h:94
itk::Image::operator[]
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
Definition: itkImage.h:244
itk::Image::GetNeighborhoodAccessor
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
Definition: itkImage.h:314
itk::Image::operator!=
friend bool operator!=(const Image &lhs, const Image &rhs)
Definition: itkImage.h:362
itk::IndexValueType
signed long IndexValueType
Definition: itkIntTypes.h:90
itkDefaultPixelAccessorFunctor.h
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
itk::Point< PointValueType, VImageDimension >
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::Image::GetBufferPointer
virtual TPixel * GetBufferPointer()
Definition: itkImage.h:249
itk::Image::operator[]
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
Definition: itkImage.h:238
itk::Image< TNode *, VImageDimension >::PixelContainerPointer
typename PixelContainer::Pointer PixelContainerPointer
Definition: itkImage.h:162
itk::Image::GetPixelContainer
PixelContainer * GetPixelContainer()
Definition: itkImage.h:262
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