ITK  4.13.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:
50  typedef typename TInputImage::PixelType InputPixelType;
51 
54 
56  itkNewMacro( Self );
57 
59  itkStaticConstMacro( ImageDimension, unsigned int, TInputImage::ImageDimension );
60 
62  typedef typename Superclass::OutputType OutputType;
63 
65  typedef typename Superclass::InputImageType InputImageType;
66 
68  typedef typename Superclass::RealType RealType;
69 
71  typedef typename Superclass::IndexType IndexType;
72 
74  typedef typename Superclass::ContinuousIndexType ContinuousIndexType;
75 
77 
78  // The interpolator used for individual binary masks corresponding to each label
79  typedef TInterpolator<LabelSelectionAdaptorType,TCoordRep> InternalInterpolatorType;
80 
85  const ContinuousIndexType & cindex ) const
86  {
87  return this->EvaluateAtContinuousIndex( cindex, NULL );
88  }
89 
90  virtual void SetInputImage( const TInputImage *image );
91 
92 protected:
95 
96  std::vector<typename InternalInterpolatorType::Pointer> m_InternalInterpolators;
97  std::vector<typename LabelSelectionAdaptorType::Pointer> m_LabelSelectionAdaptors;
98  typedef std::set<typename TInputImage::PixelType> LabelSetType;
100 
101 private:
102  LabelImageGenericInterpolateImageFunction( const Self& ); //purposely not implemented
103  void operator=( const Self& ); //purposely not implemented
104 
108  virtual OutputType EvaluateAtContinuousIndex(
109  const ContinuousIndexType &, OutputType * ) const;
110 };
111 
112 } // end namespace itk
113 
114 #ifndef ITK_MANUAL_INSTANTIATION
115 #include "itkLabelImageGenericInterpolateImageFunction.hxx"
116 #endif
117 
118 #endif
LabelSelectionImageAdaptor< TInputImage, double > LabelSelectionAdaptorType
Light weight base class for most itk classes.
Presents a label image as a binary image of one label.
std::vector< typename LabelSelectionAdaptorType::Pointer > m_LabelSelectionAdaptors
virtual OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &cindex) const
Superclass::ContinuousIndexType ContinuousIndexType
Base class for all image interpolaters.
InterpolateImageFunction< TInputImage, TCoordRep > Superclass
Interpolation function for multi-label images that implicitly interpolates each unique value in the i...
TInterpolator< LabelSelectionAdaptorType, TCoordRep > InternalInterpolatorType