ITK  5.4.0
Insight Toolkit
itkVectorImage.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 itkVectorImage_h
19#define itkVectorImage_h
20
21#include "itkImageRegion.h"
26#include "itkWeakPointer.h"
27
28namespace itk
29{
80template <typename TPixel, unsigned int VImageDimension = 3>
81class ITK_TEMPLATE_EXPORT VectorImage : public ImageBase<VImageDimension>
82{
83public:
84 ITK_DISALLOW_COPY_AND_MOVE(VectorImage);
85
92
94 itkNewMacro(Self);
95
97 itkOverrideGetNameOfClassMacro(VectorImage);
98
104
108 using InternalPixelType = TPixel;
109
112
114
118
122
126
131 static constexpr unsigned int ImageDimension = VImageDimension;
132
134 using typename Superclass::IndexType;
135 using typename Superclass::IndexValueType;
136
138 using typename Superclass::OffsetType;
139
141 using typename Superclass::SizeType;
142
145
147 using typename Superclass::DirectionType;
148
151 using typename Superclass::RegionType;
152
155 using typename Superclass::SpacingType;
156
159 using typename Superclass::PointType;
160
164
166 using typename Superclass::OffsetValueType;
167
168 using VectorLengthType = unsigned int;
169
186 template <typename UPixelType, unsigned int VUImageDimension = VImageDimension>
187 struct Rebind
188 {
190 };
191
193 template <typename UElementType, unsigned int VUImageDimension>
194 struct Rebind<VariableLengthVector<UElementType>, VUImageDimension>
195 {
197 };
199
200 template <typename UPixelType, unsigned int VUImageDimension = VImageDimension>
202
205 void
206 Allocate(bool UseValueInitialization = false) override;
207
210 void
211 Initialize() override;
212
215 void
216 FillBuffer(const PixelType & value);
217
223 void
224 SetPixel(const IndexType & index, const PixelType & value)
225 {
226 OffsetValueType offset = m_VectorLength * this->FastComputeOffset(index);
227
228 for (VectorLengthType i = 0; i < m_VectorLength; ++i)
229 {
230 (*m_Buffer)[offset + i] = value[i];
231 }
232 }
233
239 const PixelType
240 GetPixel(const IndexType & index) const
241 {
242 OffsetValueType offset = m_VectorLength * this->FastComputeOffset(index);
243
244 // Do not create a local for this method, to use return value
245 // optimization.
246 return PixelType(&((*m_Buffer)[offset]), m_VectorLength);
247 }
248
258 PixelType
259 GetPixel(const IndexType & index)
260 {
261 OffsetValueType offset = m_VectorLength * this->FastComputeOffset(index);
262
263 // Correctness of this method relies of return value optimization, do
264 // not create a local for the value.
265 return PixelType(&((*m_Buffer)[offset]), m_VectorLength);
266 }
267
277 PixelType operator[](const IndexType & index) { return this->GetPixel(index); }
278
283 const PixelType operator[](const IndexType & index) const { return this->GetPixel(index); }
284
287 InternalPixelType *
289 {
290 return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
291 }
292 const InternalPixelType *
294 {
295 return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
296 }
300 PixelContainer *
302 {
303 return m_Buffer.GetPointer();
304 }
305
307 const PixelContainer *
309 {
310 return m_Buffer.GetPointer();
311 }
312
315 void
317
328 virtual void
329 Graft(const Self * image);
330
334 {
335 return AccessorType(m_VectorLength);
336 }
337
339 const AccessorType
341 {
342 return AccessorType(m_VectorLength);
343 }
344
346 NeighborhoodAccessorFunctorType
348 {
349 return NeighborhoodAccessorFunctorType(m_VectorLength);
350 }
351
353 const NeighborhoodAccessorFunctorType
355 {
356 return NeighborhoodAccessorFunctorType(m_VectorLength);
357 }
358
360 itkSetMacro(VectorLength, VectorLengthType);
361 itkGetConstReferenceMacro(VectorLength, VectorLengthType);
365 unsigned int
367
368 void
369 SetNumberOfComponentsPerPixel(unsigned int n) override;
370
371protected:
372 VectorImage() = default;
373 void
374 PrintSelf(std::ostream & os, Indent indent) const override;
375
376 ~VectorImage() override = default;
377 void
378 Graft(const DataObject * data) override;
379 using Superclass::Graft;
380
381private:
383 VectorLengthType m_VectorLength{ 0 };
384
387};
388} // end namespace itk
389
390#ifndef ITK_MANUAL_INSTANTIATION
391# include "itkVectorImage.hxx"
392#endif
393
394#endif
Base class for all data objects in ITK.
This class provides a common API for pixel accessors for Image and VectorImage. (between the DefaultV...
Give access to partial aspects of a type.
Base class for templated image classes.
Definition: itkImageBase.h:115
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:147
Defines an itk::Image front-end to a standard C-array.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for most ITK classes.
Definition: itkObject.h:62
Provides accessor interfaces to Access pixels and is meant to be used on pointers to pixels held by t...
Templated n-dimensional vector image class.
void Initialize() override
typename PixelContainer::ConstPointer PixelContainerConstPointer
const AccessorType GetPixelAccessor() const
typename Rebind< UPixelType, VUImageDimension >::Type RebindImageType
const PixelType operator[](const IndexType &index) const
Access a pixel.
unsigned int VectorLengthType
InternalPixelType IOPixelType
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
void SetPixelContainer(PixelContainer *container)
virtual void Graft(const Self *image)
void Graft(const DataObject *data) override
void FillBuffer(const PixelType &value)
const InternalPixelType * GetBufferPointer() const
PixelContainer * GetPixelContainer()
VectorImage()=default
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
void SetNumberOfComponentsPerPixel(unsigned int n) override
unsigned int GetNumberOfComponentsPerPixel() const override
PixelType GetPixel(const IndexType &index)
Get a "reference" to a pixel. This result cannot be used as an lvalue because the pixel is converted ...
~VectorImage() override=default
void Allocate(bool UseValueInitialization=false) override
typename PixelContainer::Pointer PixelContainerPointer
const PixelType GetPixel(const IndexType &index) const
Get a pixel (read only version).
void SetPixel(const IndexType &index, const PixelType &value)
Set a pixel value.
AccessorType GetPixelAccessor()
PixelType operator[](const IndexType &index)
Access a pixel. This result cannot be used as an lvalue because the pixel is converted on the fly to ...
const PixelContainer * GetPixelContainer() const
void PrintSelf(std::ostream &os, Indent indent) const override
InternalPixelType * GetBufferPointer()
Implements a weak reference to an object.
SmartPointer< const Self > ConstPointer
static Pointer New()
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long IndexValueType
Definition: itkIntTypes.h:90
long OffsetValueType
Definition: itkIntTypes.h:94
A structure which enable changing any image class' pixel type to another.