ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkVectorImage.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 itkVectorImage_h
19 #define itkVectorImage_h
20 
21 #include "itkImageRegion.h"
26 #include "itkWeakPointer.h"
27 
28 namespace itk
29 {
80 template< typename TPixel, unsigned int VImageDimension = 3 >
81 class ITK_TEMPLATE_EXPORT VectorImage:
82  public ImageBase< VImageDimension >
83 {
84 public:
85  ITK_DISALLOW_COPY_AND_ASSIGN(VectorImage);
86 
88  using Self = VectorImage;
93 
95  itkNewMacro(Self);
96 
98  itkTypeMacro(VectorImage, ImageBase);
99 
105 
109  using InternalPixelType = TPixel;
110 
113 
115 
119 
123 
127 
132  static constexpr unsigned int ImageDimension = VImageDimension;
133 
137 
139  using OffsetType = typename Superclass::OffsetType;
140 
142  using SizeType = typename Superclass::SizeType;
143 
146 
149 
153 
156  using SpacingType = typename Superclass::SpacingType;
157 
161 
165 
168 
169  using VectorLengthType = unsigned int;
170 
187  template <typename UPixelType, unsigned int NUImageDimension = VImageDimension>
188  struct Rebind
189  {
191  };
192 
194  template <typename UElementType, unsigned int NUImageDimension>
195  struct Rebind< VariableLengthVector< UElementType >, NUImageDimension>
196  {
198  };
200 
201  template <typename UPixelType, unsigned int NUImageDimension = VImageDimension>
202  using RebindImageType = typename Rebind<UPixelType, NUImageDimension>::Type;
203 
206  void Allocate(bool UseDefaultConstructor = false) override;
207 
210  void Initialize() override;
211 
214  void FillBuffer(const PixelType & value);
215 
221  void SetPixel(const IndexType & index, const PixelType & value)
222  {
223  OffsetValueType offset = m_VectorLength * this->FastComputeOffset(index);
224 
225  for ( VectorLengthType i = 0; i < m_VectorLength; i++ )
226  {
227  ( *m_Buffer )[offset + i] = value[i];
228  }
229  }
230 
236  const PixelType GetPixel(const IndexType & index) const
237  {
238  OffsetValueType offset = m_VectorLength * this->FastComputeOffset(index);
239 
240  // Do not create a local for this method, to use return value
241  // optimization.
242  return PixelType(&( ( *m_Buffer )[offset] ), m_VectorLength);
243  }
244 
254  PixelType GetPixel(const IndexType & index)
255  {
256  OffsetValueType offset = m_VectorLength * this->FastComputeOffset(index);
257 
258  // Correctness of this method relies of return value optimization, do
259  // not create a local for the value.
260  return PixelType(&( ( *m_Buffer )[offset] ), m_VectorLength);
261  }
262 
272  PixelType operator[](const IndexType & index) { return this->GetPixel(index); }
273 
278  const PixelType operator[](const IndexType & index) const { return this->GetPixel(index); }
279 
283  {
284  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
285  }
287  {
288  return m_Buffer ? m_Buffer->GetBufferPointer() : nullptr;
289  }
291 
293  PixelContainer * GetPixelContainer() { return m_Buffer.GetPointer(); }
294 
296  const PixelContainer * GetPixelContainer() const { return m_Buffer.GetPointer(); }
297 
300  void SetPixelContainer(PixelContainer *container);
301 
312  virtual void Graft(const Self *data);
313 
315  AccessorType GetPixelAccessor() { return AccessorType(m_VectorLength); }
316 
318  const AccessorType GetPixelAccessor() const { return AccessorType(m_VectorLength); }
319 
322  {
323  return NeighborhoodAccessorFunctorType(m_VectorLength);
324  }
325 
328  {
329  return NeighborhoodAccessorFunctorType(m_VectorLength);
330  }
331 
333  itkSetMacro(VectorLength, VectorLengthType);
334  itkGetConstReferenceMacro(VectorLength, VectorLengthType);
336 
338  unsigned int GetNumberOfComponentsPerPixel() const override;
339 
340  void SetNumberOfComponentsPerPixel(unsigned int n) override;
341 
342 protected:
343  VectorImage();
344  void PrintSelf(std::ostream & os, Indent indent) const override;
345 
346  ~VectorImage() override = default;
347  void Graft(const DataObject *data) override;
348  using Superclass::Graft;
349 private:
351  VectorLengthType m_VectorLength{0};
352 
355 };
356 } // end namespace itk
357 
358 #ifndef ITK_MANUAL_INSTANTIATION
359 #include "itkVectorImage.hxx"
360 #endif
361 
362 #endif
const InternalPixelType * GetBufferPointer() const
typename Rebind< UPixelType, NUImageDimension >::Type RebindImageType
typename PixelContainer::ConstPointer PixelContainerConstPointer
PixelContainer * GetPixelContainer()
PixelType GetPixel(const IndexType &index)
Get a &quot;reference&quot; to a pixel. This result cannot be used as an lvalue because the pixel is converted ...
An image region represents a structured region of data.
Templated n-dimensional vector image class.
Implements a weak reference to an object.
A structure which enable changing any image class&#39; pixel type to another.
const AccessorType GetPixelAccessor() const
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Give access to partial aspects of a type.
This class provides a common API for pixel accessors for Image and VectorImage. (between the DefaultV...
typename PixelContainer::Pointer PixelContainerPointer
Represents an array whose length can be defined at run-time.
InternalPixelType * GetBufferPointer()
signed long IndexValueType
Definition: itkIntTypes.h:90
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
AccessorType GetPixelAccessor()
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:138
class ITK_TEMPLATE_EXPORT ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:105
void SetPixel(const IndexType &index, const PixelType &value)
Set a pixel value.
const PixelType GetPixel(const IndexType &index) const
Get a pixel (read only version).
Control indentation during Print() invocation.
Definition: itkIndent.h:49
PixelContainerPointer m_Buffer
const PixelType operator[](const IndexType &index) const
Access a pixel.
Base class for most ITK classes.
Definition: itkObject.h:60
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 ...
typename IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:133
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
signed long OffsetValueType
Definition: itkIntTypes.h:94
const PixelContainer * GetPixelContainer() const
Base class for all data objects in ITK.
Provides accessor interfaces to Access pixels and is meant to be used on pointers to pixels held by t...
Defines an itk::Image front-end to a standard C-array.