ITK  4.8.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 VectorImage:
82  public ImageBase< VImageDimension >
83 {
84 public:
86  typedef VectorImage Self;
91 
93  itkNewMacro(Self);
94 
96  itkTypeMacro(VectorImage, ImageBase);
97 
103 
107  typedef TPixel InternalPixelType;
108 
111 
113 
117 
121 
126 
131  itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
132 
136 
139 
141  typedef typename Superclass::SizeType SizeType;
142 
145 
148 
152 
156 
160 
164 
167 
168  typedef unsigned int VectorLengthType;
169 
185  template <typename UPixelType, unsigned int UImageDimension = VImageDimension>
186  struct Rebind
187  {
189  };
190 
192  template <typename UElementType, unsigned int UImageDimension>
193  struct Rebind< VariableLengthVector< UElementType >, UImageDimension>
194  {
196  };
197 
202  virtual void Allocate(bool UseDefaultConstructor = false) ITK_OVERRIDE;
203 
206  virtual void Initialize() ITK_OVERRIDE;
207 
210  void FillBuffer(const PixelType & value);
211 
217  void SetPixel(const IndexType & index, const PixelType & value)
218  {
219  OffsetValueType offset = m_VectorLength * this->ComputeOffset(index);
220 
221  for ( VectorLengthType i = 0; i < m_VectorLength; i++ )
222  {
223  ( *m_Buffer )[offset + i] = value[i];
224  }
225  }
226 
232  const PixelType GetPixel(const IndexType & index) const
233  {
234  OffsetValueType offset = m_VectorLength * this->ComputeOffset(index);
235 
236  // Do not create a local for this method, to use return value
237  // optimization.
238  return PixelType(&( ( *m_Buffer )[offset] ), m_VectorLength);
239  }
240 
250  PixelType GetPixel(const IndexType & index)
251  {
252  OffsetValueType offset = m_VectorLength * this->ComputeOffset(index);
253 
254  // Correctness of this method relies of return value optimization, do
255  // not create a local for the value.
256  return PixelType(&( ( *m_Buffer )[offset] ), m_VectorLength);
257  }
258 
268  PixelType operator[](const IndexType & index) { return this->GetPixel(index); }
269 
274  const PixelType operator[](const IndexType & index) const { return this->GetPixel(index); }
275 
279  {
280  return m_Buffer ? m_Buffer->GetBufferPointer() : ITK_NULLPTR;
281  }
283  {
284  return m_Buffer ? m_Buffer->GetBufferPointer() : ITK_NULLPTR;
285  }
287 
290 
292  const PixelContainer * GetPixelContainer() const { return m_Buffer.GetPointer(); }
293 
296  void SetPixelContainer(PixelContainer *container);
297 
308  virtual void Graft(const DataObject *data) ITK_OVERRIDE;
309 
312 
315 
318  {
320  }
321 
324  {
326  }
327 
329  itkSetMacro(VectorLength, VectorLengthType);
330  itkGetConstReferenceMacro(VectorLength, VectorLengthType);
332 
334  virtual unsigned int GetNumberOfComponentsPerPixel() const ITK_OVERRIDE;
335 
336  virtual void SetNumberOfComponentsPerPixel(unsigned int n) ITK_OVERRIDE;
337 
338 protected:
339  VectorImage();
340  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
341 
342  virtual ~VectorImage() {}
343 
344 private:
345  VectorImage(const Self &); // purposely not implementated
346  void operator=(const Self &); //purposely not implemented
347 
350 
353 };
354 } // end namespace itk
355 
356 #ifndef ITK_MANUAL_INSTANTIATION
357 #include "itkVectorImage.hxx"
358 #endif
359 
360 #endif
VectorImage Self
const InternalPixelType * GetBufferPointer() const
Superclass::OffsetValueType OffsetValueType
static const unsigned int ImageDimension
virtual void Initialize() override
virtual void SetNumberOfComponentsPerPixel(unsigned int n) override
ImportImageContainer< SizeValueType, InternalPixelType > PixelContainer
Index< VImageDimension > IndexType
Definition: itkImageBase.h:134
PixelContainer::Pointer PixelContainerPointer
virtual void Allocate(bool UseDefaultConstructor=false) override
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 ...
Superclass::IndexType IndexType
Superclass::SizeType SizeType
Templated n-dimensional vector image class.
ObjectType * GetPointer() const
SmartPointer< const Self > ConstPointer
AccessorType GetPixelAccessor(void)
Implements a weak reference to an object.
A structure which enable changing any image class&#39; pixel type to another.
Size< VImageDimension > SizeType
Definition: itkImageBase.h:143
Superclass::OffsetType OffsetType
Point< PointValueType, VImageDimension > PointType
Definition: itkImageBase.h:159
NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
Matrix< SpacePrecisionType, VImageDimension, VImageDimension > DirectionType
Definition: itkImageBase.h:164
Superclass::DirectionType DirectionType
Give access to partial aspects of a type.
This class provides a common API for pixel accessors for Image and VectorImage. (between the DefaultV...
void operator=(const Self &)
Superclass::PointType PointType
Represents an array whose length can be defined at run-time.
Offset< VImageDimension > OffsetType
Definition: itkImageBase.h:139
InternalPixelType * GetBufferPointer()
DefaultVectorPixelAccessor< InternalPixelType > AccessorType
Vector< SpacingValueType, VImageDimension > SpacingType
Definition: itkImageBase.h:154
Superclass::SpacingType SpacingType
VariableLengthVector< TPixel > PixelType
itk::VectorImage< UPixelType, UImageDimension > Type
InternalPixelType IOPixelType
OffsetType::OffsetValueType OffsetValueType
Definition: itkImageBase.h:140
VectorImageNeighborhoodAccessorFunctor< Self > NeighborhoodAccessorFunctorType
unsigned int VectorLengthType
WeakPointer< const Self > ConstWeakPointer
Base class for templated image classes.
Definition: itkImageBase.h:112
const AccessorType GetPixelAccessor(void) const
void PrintSelf(std::ostream &os, Indent indent) const override
void SetPixel(const IndexType &index, const PixelType &value)
Set a pixel value.
Superclass::RegionType RegionType
VectorLengthType m_VectorLength
const PixelType GetPixel(const IndexType &index) const
Get a pixel (read only version).
Control indentation during Print() invocation.
Definition: itkIndent.h:49
virtual unsigned int GetNumberOfComponentsPerPixel() const override
OffsetValueType ComputeOffset(const IndexType &ind) const
Definition: itkImageBase.h:332
PixelContainerPointer m_Buffer
void FillBuffer(const PixelType &value)
ImageRegion< VImageDimension > RegionType
Definition: itkImageBase.h:147
const PixelType operator[](const IndexType &index) const
Access a pixel.
Superclass::IndexValueType IndexValueType
DefaultVectorPixelAccessorFunctor< Self > AccessorFunctorType
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 ...
PixelContainer::ConstPointer PixelContainerConstPointer
SmartPointer< Self > Pointer
const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
void SetPixelContainer(PixelContainer *container)
virtual void Graft(const DataObject *data) override
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...
IndexType::IndexValueType IndexValueType
Definition: itkImageBase.h:135
ImageBase< VImageDimension > Superclass
Defines an itk::Image front-end to a standard C-array.