ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkLabelImageGenericInterpolateImageFunction.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 itkLabelImageGenericInterpolateImageFunction_h
19 #define itkLabelImageGenericInterpolateImageFunction_h
20 
23 #include <vector>
24 #include <set>
25 
26 namespace itk
27 {
28 
40 template <typename TInputImage,template<class, typename> class TInterpolator, typename TCoordRep=double >
42  public InterpolateImageFunction<TInputImage, TCoordRep>
43 {
44 public:
45  ITK_DISALLOW_COPY_AND_ASSIGN(LabelImageGenericInterpolateImageFunction);
46 
52  using InputPixelType = typename TInputImage::PixelType;
53 
56 
58  itkNewMacro( Self );
59 
61  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
62 
64  using OutputType = typename Superclass::OutputType;
65 
67  using InputImageType = typename Superclass::InputImageType;
68 
70  using RealType = typename Superclass::RealType;
71 
73  using IndexType = typename Superclass::IndexType;
74 
76  using ContinuousIndexType = typename Superclass::ContinuousIndexType;
77 
79 
80  // The interpolator used for individual binary masks corresponding to each label
81  using InternalInterpolatorType = TInterpolator<LabelSelectionAdaptorType,TCoordRep>;
82 
87  const ContinuousIndexType & cindex ) const override
88  {
89  return this->EvaluateAtContinuousIndex( cindex, nullptr );
90  }
91 
92  void SetInputImage( const TInputImage *image ) override;
93 
94 protected:
97 
98  std::vector<typename InternalInterpolatorType::Pointer> m_InternalInterpolators;
99  std::vector<typename LabelSelectionAdaptorType::Pointer> m_LabelSelectionAdaptors;
100  using LabelSetType = std::set<typename TInputImage::PixelType>;
102 
103 private:
107  virtual OutputType EvaluateAtContinuousIndex(
108  const ContinuousIndexType &, OutputType * ) const;
109 };
110 
111 } // end namespace itk
112 
113 #ifndef ITK_MANUAL_INSTANTIATION
114 #include "itkLabelImageGenericInterpolateImageFunction.hxx"
115 #endif
116 
117 #endif
TInterpolator< LabelSelectionAdaptorType, TCoordRep > InternalInterpolatorType
Light weight base class for most itk classes.
typename Superclass::ContinuousIndexType ContinuousIndexType
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const override
typename Superclass::OutputType OutputType
Presents a label image as a binary image of one label.
typename NumericTraits< typename TInputImage::PixelType >::RealType RealType
std::vector< typename LabelSelectionAdaptorType::Pointer > m_LabelSelectionAdaptors
typename Superclass::IndexType IndexType
Base class for all image interpolaters.
typename Superclass::InputImageType InputImageType
Interpolation function for multi-label images that implicitly interpolates each unique value in the i...