ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkImage.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 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 {
74 template< typename TPixel, unsigned int VImageDimension = 2 >
75 class ITK_TEMPLATE_EXPORT Image:public ImageBase< VImageDimension >
76 {
77 public:
78  ITK_DISALLOW_COPY_AND_ASSIGN(Image);
79 
81  using Self = Image;
86 
88  itkNewMacro(Self);
89 
91  itkTypeMacro(Image, ImageBase);
92 
95  using PixelType = TPixel;
96 
98  using ValueType = TPixel;
99 
104  using InternalPixelType = TPixel;
105 
107 
112 
116 
118  using ImageDimensionType = typename Superclass::ImageDimensionType;
119 
123 
125  using OffsetType = typename Superclass::OffsetType;
126 
128  using SizeType = typename Superclass::SizeType;
130 
133 
136 
140 
143  using SpacingType = typename Superclass::SpacingType;
144  using SpacingValueType = typename Superclass::SpacingValueType;
145 
149 
153 
156 
163  template <typename UPixelType, unsigned int NUImageDimension = VImageDimension>
164  struct Rebind
165  {
167  };
168 
169 
170  template <typename UPixelType, unsigned int NUImageDimension = VImageDimension>
172 
173 
176  void Allocate(bool initializePixels = false) override;
177 
180  void Initialize() override;
181 
184  void FillBuffer(const TPixel & value);
185 
191  void SetPixel(const IndexType & index, const TPixel & value)
192  {
193  OffsetValueType offset = this->FastComputeOffset(index);
194  ( *m_Buffer )[offset] = value;
195  }
196 
201  const TPixel & GetPixel(const IndexType & index) const
202  {
203  OffsetValueType offset = this->FastComputeOffset(index);
204  return ( ( *m_Buffer )[offset] );
205  }
206 
211  TPixel & GetPixel(const IndexType & index)
212  {
213  OffsetValueType offset = this->FastComputeOffset(index);
214  return ( ( *m_Buffer )[offset] );
215  }
216 
221  TPixel & operator[](const IndexType & index)
222  { return this->GetPixel(index); }
223 
228  const TPixel & operator[](const IndexType & index) const
229  { return this->GetPixel(index); }
230 
233  virtual TPixel * GetBufferPointer()
234  { return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr; }
235  virtual const TPixel * GetBufferPointer() const
236  { return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr; }
238 
241  { return m_Buffer.GetPointer(); }
242 
244  { return m_Buffer.GetPointer(); }
245 
248  void SetPixelContainer(PixelContainer *container);
249 
260  virtual void Graft(const Self *data);
261 
264  { return AccessorType(); }
265 
268  { return AccessorType(); }
269 
272  { return NeighborhoodAccessorFunctorType(); }
273 
276  { return NeighborhoodAccessorFunctorType(); }
277 
278  unsigned int GetNumberOfComponentsPerPixel() const override;
279 
280 protected:
281  Image();
282  void PrintSelf(std::ostream & os, Indent indent) const override;
283  void Graft(const DataObject *data) override;
284 
285  ~Image() override = default;
286 
292  void ComputeIndexToPhysicalPointMatrices() override;
293  using Superclass::Graft;
294 private:
295 
298 };
299 } // end namespace itk
300 
301 #ifndef ITK_MANUAL_INSTANTIATION
302 #include "itkImage.hxx"
303 #endif
304 
305 #endif
const PixelContainer * GetPixelContainer() const
Definition: itkImage.h:243
unsigned int ImageDimensionType
Definition: itkImageBase.h:123
SpacePrecisionType SpacingValueType
Definition: itkImageBase.h:151
unsigned long SizeValueType
Definition: itkIntTypes.h:83
virtual TPixel * GetBufferPointer()
Definition: itkImage.h:233
An image region represents a structured region of data.
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
Definition: itkImage.h:211
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
Definition: itkImage.h:201
Implements a weak reference to an object.
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
Definition: itkImage.h:275
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Definition: itkImage.h:271
AccessorType GetPixelAccessor()
Definition: itkImage.h:263
Provides accessor interfaces to Get pixels and is meant to be used on pointers contained within Neigh...
signed long IndexValueType
Definition: itkIntTypes.h:90
const AccessorType GetPixelAccessor() const
Definition: itkImage.h:267
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
typename SizeType::SizeValueType SizeValueType
Definition: itkImageBase.h:142
Provides a common API for pixel accessors for Image and VectorImage.
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
PixelContainerPointer m_Buffer
Definition: itkImage.h:297
virtual const TPixel * GetBufferPointer() const
Definition: itkImage.h:235
Base class for templated image classes.
Definition: itkImageBase.h:105
Control indentation during Print() invocation.
Definition: itkIndent.h:49
PixelContainer * GetPixelContainer()
Definition: itkImage.h:240
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
Definition: itkImage.h:228
typename PixelContainer::Pointer PixelContainerPointer
Definition: itkImage.h:151
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition: itkImage.h:152
Give access to partial aspects a type.
Base class for most ITK classes.
Definition: itkObject.h:60
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
Definition: itkImage.h:221
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:133
signed long OffsetValueType
Definition: itkIntTypes.h:94
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
Definition: itkImage.h:191
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75
Defines an itk::Image front-end to a standard C-array.