ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkVectorImageToImageAdaptor.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 itkVectorImageToImageAdaptor_h
19 #define itkVectorImageToImageAdaptor_h
20 
21 #include "itkImageAdaptor.h"
22 #include "itkVectorImage.h"
23 
24 namespace itk
25 {
26 namespace Accessor
27 {
46 template< typename TType >
48  : private DefaultVectorPixelAccessor< TType >
49 {
50 public:
51 
52  typedef unsigned int VectorLengthType;
53 
56  typedef TType ExternalType;
57 
59  typedef TType InternalType;
60 
62 
63  inline void Set(ActualPixelType output, const ExternalType & input) const
64  {
65  output[m_ComponentIdx] = input;
66  }
67 
68  inline void Set(InternalType &output, const ExternalType & input,
69  const unsigned long offset) const
70  {
71  return Set( Superclass::Get( output, offset ), input );
72  }
73 
74  inline ExternalType Get(const ActualPixelType & input) const
75  {
76  ExternalType output;
77 
78  output = input[m_ComponentIdx];
79  return output;
80  }
81 
82  inline ExternalType Get(const InternalType &input, const SizeValueType offset) const
83  {
84  return Get( Superclass::Get(input, offset) );
85  }
86 
88  {
89  m_ComponentIdx = idx;
90  }
91 
93  {
94  return m_ComponentIdx;
95  }
96 
99  {
101  }
102 
105 
106  VectorImageToImagePixelAccessor( unsigned int length = 1)
107  :m_ComponentIdx(0)
108  {
109  Superclass::SetVectorLength( length );
110  }
111 
112 protected:
114 
115 private:
117 };
118 } // end namespace Accessor
119 
141 template< typename TPixelType, unsigned int Dimension >
143  ImageAdaptor< VectorImage< TPixelType, Dimension >,
144  Accessor::VectorImageToImagePixelAccessor< TPixelType > >
145 {
146 public:
150  typedef ImageAdaptor< VectorImageType,
152 
155 
157  itkNewMacro(Self);
158 
161 
168 
171 
172  // Set/GetMethods to set the component to be extracted.
174  {
175  this->GetPixelAccessor().SetExtractComponentIdx(componentIdx);
176  }
177 
178  // Set/GetMethods to set the component to be extracted.
180  {
181  return this->GetPixelAccessor().GetExtractComponentIdx();
182  }
183 
184 protected:
187 
188 private:
189  VectorImageToImageAdaptor(const Self &) ITK_DELETE_FUNCTION;
190  void operator=(const Self &) ITK_DELETE_FUNCTION;
191 };
192 } // end namespace itk
193 
194 #endif
Presents a VectorImage and extracts a component from it into an image.
Superclass::PixelContainerPointer PixelContainerPointer
Superclass::PixelContainerConstPointer PixelContainerConstPointer
TImage::PixelContainer PixelContainer
ExternalType Get(const InternalType &input, const SizeValueType offset) const
ExternalType Get(const ActualPixelType &input) const
VectorLengthType GetExtractComponentIndex() const
unsigned long SizeValueType
Definition: itkIntTypes.h:143
Give access to partial aspects of a type.
TImage::PixelContainerConstPointer PixelContainerConstPointer
Represents an array whose length can be defined at run-time.
VectorImage< TPixelType, Dimension > VectorImageType
void SetExtractComponentIndex(VectorLengthType componentIdx)
ExternalType Get(const InternalType &input, const SizeValueType offset) const
TImage::PixelContainerPointer PixelContainerPointer
VectorImageType::VectorLengthType VectorLengthType
Give access to partial aspects of voxels from an Image.
void Set(ActualPixelType output, const ExternalType &input) const
void Set(InternalType &output, const ExternalType &input, const unsigned long offset) const
Base class for all data objects in ITK.
ImageAdaptor< VectorImageType, Accessor::VectorImageToImagePixelAccessor< TPixelType > > Superclass