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;
78
79 output = input[m_ComponentIdx];
80 return output;
81 }
82
83 inline ExternalType
84 Get(const InternalType & input, const SizeValueType offset) const
85 {
86 return Get(Superclass::Get(input, offset));
87 }
88
89 void
91 {
92 m_ComponentIdx = idx;
93 }
94
97 {
98 return m_ComponentIdx;
99 }
100
102 void
104 {
106 }
107
111 {
113 }
114
116
117protected:
119
120private:
122};
123} // end namespace Accessor
124
147template <typename TPixelType, unsigned int Dimension>
149 : public ImageAdaptor<VectorImage<TPixelType, Dimension>, Accessor::VectorImageToImagePixelAccessor<TPixelType>>
150{
151public:
152 ITK_DISALLOW_COPY_AND_MOVE(VectorImageToImageAdaptor);
153
158
161
163 itkNewMacro(Self);
164
166 itkOverrideGetNameOfClassMacro(VectorImageToImageAdaptor);
167
170 using typename Superclass::PixelContainer;
173 using typename Superclass::IOPixelType;
174
177
178 // Set/GetMethods to set the component to be extracted.
179 void
181 {
182 this->GetPixelAccessor().SetExtractComponentIdx(componentIdx);
183 }
184
185 // Set/GetMethods to set the component to be extracted.
188 {
190 }
191
192protected:
194 ~VectorImageToImageAdaptor() override = default;
195};
196} // end namespace itk
197
198#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