ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkLabelImageGenericInterpolateImageFunction.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 itkLabelImageGenericInterpolateImageFunction_h
19#define itkLabelImageGenericInterpolateImageFunction_h
20
23#include <vector>
24#include <set>
25
26namespace itk
27{
28
39
40template <typename TInputImage,
41 template <typename TInterpInputImage, typename TCoordRep> class TInterpolator,
42 typename TCoordRep = double>
43class ITK_EXPORT LabelImageGenericInterpolateImageFunction : public InterpolateImageFunction<TInputImage, TCoordRep>
44{
45public:
46 ITK_DISALLOW_COPY_AND_MOVE(LabelImageGenericInterpolateImageFunction);
47
53 using InputPixelType = typename TInputImage::PixelType;
54
56 itkOverrideGetNameOfClassMacro(LabelImageGenericInterpolateImageFunction);
57
59 itkNewMacro(Self);
60
62 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
63
66
69
72
75
78
81
83
84 // The interpolator used for individual binary masks corresponding to each label
85 using InternalInterpolatorType = TInterpolator<LabelSelectionAdaptorType, TCoordRep>;
86
91 EvaluateAtContinuousIndex(const ContinuousIndexType & cindex) const override
92 {
93 return this->EvaluateAtContinuousIndex(cindex, nullptr);
94 }
95
96 void
97 SetInputImage(const TInputImage * image) override;
98
105 GetRadius() const override
106 {
107 return SizeType::Filled(1);
108 }
109
110protected:
113
114 std::vector<typename InternalInterpolatorType::Pointer> m_InternalInterpolators;
115 std::vector<typename LabelSelectionAdaptorType::Pointer> m_LabelSelectionAdaptors;
116 using LabelSetType = std::set<typename TInputImage::PixelType>;
118
119private:
123 virtual OutputType
125};
126
127} // end namespace itk
128
129#ifndef ITK_MANUAL_INSTANTIATION
130# include "itkLabelImageGenericInterpolateImageFunction.hxx"
131#endif
132
133#endif
typename NumericTraits< typename TInputImage::PixelType >::RealType RealType
typename InputImageType::IndexType IndexType
ContinuousIndex< TCoordRep, Self::ImageDimension > ContinuousIndexType
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
std::vector< typename LabelSelectionAdaptorType::Pointer > m_LabelSelectionAdaptors
InterpolateImageFunction< TInputImage, TCoordRep > Superclass
LabelSelectionImageAdaptor< TInputImage, double > LabelSelectionAdaptorType
std::vector< typename InternalInterpolatorType::Pointer > m_InternalInterpolators
TInterpolator< LabelSelectionAdaptorType, TCoordRep > InternalInterpolatorType
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &, OutputType *) const
~LabelImageGenericInterpolateImageFunction() override=default
void SetInputImage(const TInputImage *image) override
Presents a label image as a binary image of one label.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....