ITK  5.2.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  * 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 {
47 template <typename TType>
49 {
50 public:
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 
115  VectorImageToImagePixelAccessor(unsigned int length = 1) { Superclass::SetVectorLength(length); }
116 
117 protected:
119 
120 private:
122 };
123 } // end namespace Accessor
124 
147 template <typename TPixelType, unsigned int Dimension>
149  : public ImageAdaptor<VectorImage<TPixelType, Dimension>, Accessor::VectorImageToImagePixelAccessor<TPixelType>>
150 {
151 public:
152  ITK_DISALLOW_COPY_AND_MOVE(VectorImageToImageAdaptor);
153 
158 
161 
163  itkNewMacro(Self);
164 
167 
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  {
189  return this->GetPixelAccessor().GetExtractComponentIdx();
190  }
191 
192 protected:
193  VectorImageToImageAdaptor() = default;
194  ~VectorImageToImageAdaptor() override = default;
195 };
196 } // end namespace itk
197 
198 #endif
itk::Accessor::VectorImageToImagePixelAccessor::Set
void Set(ActualPixelType output, const ExternalType &input) const
Definition: itkVectorImageToImageAdaptor.h:63
itk::ImageAdaptor::PixelContainerConstPointer
typename TImage::PixelContainerConstPointer PixelContainerConstPointer
Definition: itkImageAdaptor.h:241
itk::Accessor::VectorImageToImagePixelAccessor::m_ComponentIdx
VectorLengthType m_ComponentIdx
Definition: itkVectorImageToImageAdaptor.h:121
itk::VectorImageToImageAdaptor::~VectorImageToImageAdaptor
~VectorImageToImageAdaptor() override=default
itk::Accessor::VectorImageToImagePixelAccessor::GetVectorLength
VectorLengthType GetVectorLength() const
Definition: itkVectorImageToImageAdaptor.h:110
itk::VectorImageToImageAdaptor
Presents a VectorImage and extracts a component from it into an image.
Definition: itkVectorImageToImageAdaptor.h:148
itk::ImageAdaptor::PixelContainer
typename TImage::PixelContainer PixelContainer
Definition: itkImageAdaptor.h:239
itk::Accessor::VectorImageToImagePixelAccessor
Extract components from a VectorImage.
Definition: itkVectorImageToImageAdaptor.h:48
itk::VectorImage< TPixelType, Dimension >
itk::VectorImageToImageAdaptor::IOPixelType
typename Superclass::IOPixelType IOPixelType
Definition: itkVectorImageToImageAdaptor.h:173
itk::VectorImageToImageAdaptor::SetExtractComponentIndex
void SetExtractComponentIndex(VectorLengthType componentIdx)
Definition: itkVectorImageToImageAdaptor.h:180
itk::Accessor::VectorImageToImagePixelAccessor::VectorImageToImagePixelAccessor
VectorImageToImagePixelAccessor(unsigned int length=1)
Definition: itkVectorImageToImageAdaptor.h:115
itk::DefaultVectorPixelAccessor< TPixelType >::InternalType
TPixelType InternalType
Definition: itkDefaultVectorPixelAccessor.h:62
itk::VectorImageToImageAdaptor::VectorImageToImageAdaptor
VectorImageToImageAdaptor()=default
itk::SmartPointer< Self >
itkImageAdaptor.h
itk::VectorImageToImageAdaptor::PixelContainer
typename Superclass::PixelContainer PixelContainer
Definition: itkVectorImageToImageAdaptor.h:170
itk::ImageAdaptor
Give access to partial aspects of voxels from an Image.
Definition: itkImageAdaptor.h:56
itkVectorImage.h
itk::DefaultVectorPixelAccessor::Get
ExternalType Get(const InternalType &input, const SizeValueType offset) const
Definition: itkDefaultVectorPixelAccessor.h:78
itk::Accessor::VectorImageToImagePixelAccessor::VectorLengthType
unsigned int VectorLengthType
Definition: itkVectorImageToImageAdaptor.h:51
itk::Accessor::VectorImageToImagePixelAccessor::Get
ExternalType Get(const InternalType &input, const SizeValueType offset) const
Definition: itkVectorImageToImageAdaptor.h:84
itk::VectorImageToImageAdaptor::PixelContainerConstPointer
typename Superclass::PixelContainerConstPointer PixelContainerConstPointer
Definition: itkVectorImageToImageAdaptor.h:172
itk::Accessor::VectorImageToImagePixelAccessor::ExternalType
TType ExternalType
Definition: itkVectorImageToImageAdaptor.h:55
itk::Accessor::VectorImageToImagePixelAccessor::Set
void Set(InternalType &output, const ExternalType &input, const unsigned long offset) const
Definition: itkVectorImageToImageAdaptor.h:69
itk::DefaultVectorPixelAccessor
Give access to partial aspects of a type.
Definition: itkDefaultVectorPixelAccessor.h:50
itk::VectorImageToImageAdaptor::VectorLengthType
typename VectorImageType::VectorLengthType VectorLengthType
Definition: itkVectorImageToImageAdaptor.h:176
itk::Accessor::VectorImageToImagePixelAccessor::Get
ExternalType Get(const ActualPixelType &input) const
Definition: itkVectorImageToImageAdaptor.h:75
itk::VariableLengthVector
Represents an array whose length can be defined at run-time.
Definition: itkConstantBoundaryCondition.h:28
itk::Accessor::VectorImageToImagePixelAccessor::SetVectorLength
void SetVectorLength(VectorLengthType l)
Definition: itkVectorImageToImageAdaptor.h:103
itk::DefaultVectorPixelAccessor::SetVectorLength
void SetVectorLength(VectorLengthType l)
Definition: itkDefaultVectorPixelAccessor.h:87
itk::ImageAdaptor::PixelContainerPointer
typename TImage::PixelContainerPointer PixelContainerPointer
Definition: itkImageAdaptor.h:240
itk::VectorImageToImageAdaptor::GetExtractComponentIndex
VectorLengthType GetExtractComponentIndex() const
Definition: itkVectorImageToImageAdaptor.h:187
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::DefaultVectorPixelAccessor< TPixelType >::VectorLengthType
unsigned int VectorLengthType
Definition: itkDefaultVectorPixelAccessor.h:53
itk::Accessor::VectorImageToImagePixelAccessor::SetExtractComponentIdx
void SetExtractComponentIdx(VectorLengthType idx)
Definition: itkVectorImageToImageAdaptor.h:90
itk::DefaultVectorPixelAccessor::GetVectorLength
VectorLengthType GetVectorLength() const
Definition: itkDefaultVectorPixelAccessor.h:96
itk::Object
Base class for most ITK classes.
Definition: itkObject.h:62
itk::Accessor::VectorImageToImagePixelAccessor::GetExtractComponentIdx
VectorLengthType GetExtractComponentIdx() const
Definition: itkVectorImageToImageAdaptor.h:96
itk::ImageAdaptor< VectorImage< TPixelType, Dimension >, Accessor::VectorImageToImagePixelAccessor< TPixelType > >::GetPixelAccessor
AccessorType & GetPixelAccessor()
Definition: itkImageAdaptor.h:341
itk::VectorImageToImageAdaptor::PixelContainerPointer
typename Superclass::PixelContainerPointer PixelContainerPointer
Definition: itkVectorImageToImageAdaptor.h:171
itk::VectorImage< TPixelType, Dimension >::VectorLengthType
unsigned int VectorLengthType
Definition: itkVectorImage.h:168
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::ImageAdaptor::IOPixelType
PixelType IOPixelType
Definition: itkImageAdaptor.h:91
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293