ITK  6.0.0
Insight Toolkit
itkVectorImageToImageAdaptor.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 itkVectorImageToImageAdaptor_h
19#define itkVectorImageToImageAdaptor_h
20
21#include "itkImageAdaptor.h"
22#include "itkVectorImage.h"
23
24namespace itk
25{
26namespace Accessor
27{
47template <typename TType>
49{
50public:
51 using VectorLengthType = unsigned int;
52
55 using ExternalType = TType;
56
58 using InternalType = TType;
59
61
62 inline void
63 Set(ActualPixelType output, const ExternalType & input) const
64 {
65 output[m_ComponentIdx] = input;
66 }
67
68 inline void
69 Set(InternalType & output, const ExternalType & input, const unsigned long offset) const
70 {
71 return Set(Superclass::Get(output, offset), input);
72 }
73
74 inline ExternalType
75 Get(const ActualPixelType & input) const
76 {
77 ExternalType output = input[m_ComponentIdx];
78 return output;
79 }
80
81 inline ExternalType
82 Get(const InternalType & input, const SizeValueType offset) const
83 {
84 return Get(Superclass::Get(input, offset));
85 }
86
87 void
89 {
90 m_ComponentIdx = idx;
91 }
92
95 {
96 return m_ComponentIdx;
97 }
98
100 void
102 {
104 }
105
109 {
111 }
112
114
115protected:
117
118private:
120};
121} // end namespace Accessor
122
145template <typename TPixelType, unsigned int Dimension>
147 : public ImageAdaptor<VectorImage<TPixelType, Dimension>, Accessor::VectorImageToImagePixelAccessor<TPixelType>>
148{
149public:
150 ITK_DISALLOW_COPY_AND_MOVE(VectorImageToImageAdaptor);
151
156
159
161 itkNewMacro(Self);
162
164 itkOverrideGetNameOfClassMacro(VectorImageToImageAdaptor);
165
168 using typename Superclass::PixelContainer;
171 using typename Superclass::IOPixelType;
172
175
176 // Set/GetMethods to set the component to be extracted.
177 void
179 {
180 this->GetPixelAccessor().SetExtractComponentIdx(componentIdx);
181 }
182
183 // Set/GetMethods to set the component to be extracted.
186 {
188 }
189
190protected:
192 ~VectorImageToImageAdaptor() override = default;
193};
194} // end namespace itk
195
196#endif
ExternalType Get(const InternalType &input, const SizeValueType offset) const
ExternalType Get(const ActualPixelType &input) const
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.
Give access to partial aspects of a type.
ExternalType Get(const InternalType &input, const SizeValueType offset) const
Give access to partial aspects of voxels from an Image.
typename TImage::PixelContainerPointer PixelContainerPointer
typename TImage::PixelContainer PixelContainer
typename TImage::PixelContainerConstPointer PixelContainerConstPointer
Base class for most ITK classes.
Definition: itkObject.h:62
Presents a VectorImage and extracts a component from it into an image.
typename VectorImageType::VectorLengthType VectorLengthType
VectorLengthType GetExtractComponentIndex() const
~VectorImageToImageAdaptor() override=default
void SetExtractComponentIndex(VectorLengthType componentIdx)
Templated n-dimensional vector image class.
unsigned int VectorLengthType
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86