ITK  4.4.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 {
81 template< class TPixel, unsigned int VImageDimension = 3 >
82 class ITK_EXPORT VectorImage:
83  public ImageBase< VImageDimension >
84 {
85 public:
87  typedef VectorImage Self;
92 
94  itkNewMacro(Self);
95 
97  itkTypeMacro(VectorImage, ImageBase);
98 
104 
108  typedef TPixel InternalPixelType;
109 
112 
114 
118 
122 
127 
132  itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
133 
135  typedef typename Superclass::IndexType IndexType;
137 
139  typedef typename Superclass::OffsetType OffsetType;
140 
142  typedef typename Superclass::SizeType SizeType;
143 
146 
148  typedef typename Superclass::DirectionType DirectionType;
149 
152  typedef typename Superclass::RegionType RegionType;
153 
156  typedef typename Superclass::SpacingType SpacingType;
157 
160  typedef typename Superclass::PointType PointType;
161 
165 
168 
169  typedef unsigned int VectorLengthType;
170 
186  template <typename UPixelType, unsigned int UImageDimension = VImageDimension>
187  struct Rebind
188  {
190  };
191 
193  template <typename UElementType, unsigned int UImageDimension>
194  struct Rebind< VariableLengthVector< UElementType >, UImageDimension>
195  {
197  };
198 
203  void Allocate();
204 
207  virtual void Initialize();
208 
211  void FillBuffer(const PixelType & value);
212 
218  void SetPixel(const IndexType & index, const PixelType & value)
219  {
220  OffsetValueType offset = m_VectorLength * this->ComputeOffset(index);
221 
222  for ( VectorLengthType i = 0; i < m_VectorLength; i++ )
223  {
224  ( *m_Buffer )[offset + i] = value[i];
225  }
226  }
227 
233  const PixelType GetPixel(const IndexType & index) const
234  {
235  OffsetValueType offset = m_VectorLength * this->ComputeOffset(index);
236  PixelType p(&( ( *m_Buffer )[offset] ), m_VectorLength);
237 
238  return p;
239  }
240 
245  PixelType GetPixel(const IndexType & index)
246  {
247  OffsetValueType offset = m_VectorLength * this->ComputeOffset(index);
248  PixelType p(&( ( *m_Buffer )[offset] ), m_VectorLength);
249 
250  return p;
251  }
252 
258  PixelType operator[](const IndexType & index) { return this->GetPixel(index); }
259 
265  PixelType operator[](const IndexType & index) const { return this->GetPixel(index); }
266 
269  InternalPixelType * GetBufferPointer()
270  {
271  return m_Buffer ? m_Buffer->GetBufferPointer() : 0;
272  }
273  const InternalPixelType * GetBufferPointer() const
274  {
275  return m_Buffer ? m_Buffer->GetBufferPointer() : 0;
276  }
278 
280  PixelContainer * GetPixelContainer() { return m_Buffer.GetPointer(); }
281 
283  const PixelContainer * GetPixelContainer() const { return m_Buffer.GetPointer(); }
284 
287  void SetPixelContainer(PixelContainer *container);
288 
299  virtual void Graft(const DataObject *data);
300 
302  AccessorType GetPixelAccessor(void) { return AccessorType(m_VectorLength); }
303 
305  const AccessorType GetPixelAccessor(void) const { return AccessorType(m_VectorLength); }
306 
308  NeighborhoodAccessorFunctorType GetNeighborhoodAccessor()
309  {
310  return NeighborhoodAccessorFunctorType(m_VectorLength);
311  }
312 
314  const NeighborhoodAccessorFunctorType GetNeighborhoodAccessor() const
315  {
316  return NeighborhoodAccessorFunctorType(m_VectorLength);
317  }
318 
320  itkSetMacro(VectorLength, VectorLengthType);
321  itkGetConstReferenceMacro(VectorLength, VectorLengthType);
323 
325  virtual unsigned int GetNumberOfComponentsPerPixel() const;
326 
327  virtual void SetNumberOfComponentsPerPixel(unsigned int n);
328 
329 protected:
330  VectorImage();
331  void PrintSelf(std::ostream & os, Indent indent) const;
332 
333  virtual ~VectorImage() {}
334 
335 private:
336  VectorImage(const Self &); // purposely not implementated
337  void operator=(const Self &); //purposely not implemented
338 
341 
344 };
345 } // end namespace itk
346 
347 #ifndef ITK_MANUAL_INSTANTIATION
348 #include "itkVectorImage.hxx"
349 #endif
350 
351 #endif
352