ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
32namespace itk
33{
87template <typename TPixel, unsigned int VImageDimension = 2>
88class ITK_TEMPLATE_EXPORT Image : public ImageBase<VImageDimension>
89{
90public:
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>
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 }
268
269
271 PixelContainer *
273 {
274 return m_Buffer.GetPointer();
275 }
276
277 const PixelContainer *
279 {
280 return m_Buffer.GetPointer();
281 }
282
285 void
287
298 virtual void
299 Graft(const Self * image);
300
304 {
305 return AccessorType();
306 }
307
309 const AccessorType
311 {
312 return AccessorType();
313 }
314
316 NeighborhoodAccessorFunctorType
321
323 const NeighborhoodAccessorFunctorType
328
329 unsigned int
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
385protected:
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
401 using Superclass::Graft;
402
403private:
406};
407} // end namespace itk
408
409#ifndef ITK_MANUAL_INSTANTIATION
410# include "itkImage.hxx"
411#endif
412
413#endif
Base class for all data objects in ITK.
Provides a common API for pixel accessors for Image and VectorImage.
Give access to partial aspects a type.
Vector< SpacingValueType, VImageDimension > SpacingType
ImageRegion< VImageDimension > RegionType
SpacePrecisionType SpacingValueType
PointType m_Origin
virtual const RegionType & GetBufferedRegion() const
Index< VImageDimension > IndexType
Offset< VImageDimension > OffsetType
DirectionType m_Direction
OffsetValueType FastComputeOffset(const IndexType &ind) const
typename OffsetType::OffsetValueType OffsetValueType
typename IndexType::IndexValueType IndexValueType
virtual void Graft(const Self *image)
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Size< VImageDimension > SizeType
Point< PointValueType, VImageDimension > PointType
SpacingType m_Spacing
DirectionType m_InverseDirection
typename SizeType::SizeValueType SizeValueType
Templated n-dimensional image class.
Definition itkImage.h:89
void PrintSelf(std::ostream &os, Indent indent) const override
const AccessorType GetPixelAccessor() const
Definition itkImage.h:310
unsigned int GetNumberOfComponentsPerPixel() const override
AccessorType GetPixelAccessor()
Definition itkImage.h:303
SmartPointer< const Self > ConstPointer
Definition itkImage.h:97
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
virtual void Graft(const Self *image)
const PixelContainer * GetPixelContainer() const
Definition itkImage.h:278
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
WeakPointer< const Self > ConstWeakPointer
Definition itkImage.h:98
ImageBase< VImageDimension > Superclass
Definition itkImage.h:95
Image()=default
const TPixel & operator[](const IndexType &index) const
Access a pixel. This version can only be an rvalue.
Definition itkImage.h:251
void SetPixel(const IndexType &index, const TPixel &value)
Set a pixel value.
Definition itkImage.h:208
TPixel & GetPixel(const IndexType &index)
Get a reference to a pixel (e.g. for editing).
Definition itkImage.h:230
virtual const TPixel * GetBufferPointer() const
Definition itkImage.h:264
virtual TPixel * GetBufferPointer()
Definition itkImage.h:259
void Graft(const DataObject *data) override
void SetPixelContainer(PixelContainer *container)
NeighborhoodAccessorFunctor< Self > NeighborhoodAccessorFunctorType
Definition itkImage.h:128
DefaultPixelAccessor< PixelType > AccessorType
Definition itkImage.h:123
typename PixelContainer::ConstPointer PixelContainerConstPointer
Definition itkImage.h:165
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
Definition itkImage.h:324
void Allocate(bool initializePixels=false) override
itk::Image< UPixelType, VUImageDimension > RebindImageType
Definition itkImage.h:184
void FillBuffer(const TPixel &value)
typename PixelContainer::Pointer PixelContainerPointer
Definition itkImage.h:164
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Definition itkImage.h:317
void Initialize() override
~Image() override=default
ImportImageContainer< SizeValueType, PixelType > PixelContainer
Definition itkImage.h:145
DefaultPixelAccessorFunctor< Self > AccessorFunctorType
Definition itkImage.h:124
PixelContainer * GetPixelContainer()
Definition itkImage.h:272
const TPixel & GetPixel(const IndexType &index) const
Get a pixel (read only version).
Definition itkImage.h:219
void ComputeIndexToPhysicalPointMatrices() override
typename OffsetType::OffsetValueType OffsetValueType
TPixel & operator[](const IndexType &index)
Access a pixel. This version can be an lvalue.
Definition itkImage.h:241
Defines an itk::Image front-end to a standard C-array.
Control indentation during Print() invocation.
Definition itkIndent.h:50
Provides accessor interfaces to Get pixels and is meant to be used on pointers contained within Neigh...
Implements transparent reference counting.
Implements a weak reference to an object.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
itk::Image< UPixelType, VUImageDimension > Type
Definition itkImage.h:179