RTK  2.5.0
Reconstruction Toolkit
rtkLookupTableImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright RTK 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  * 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 
19 #ifndef rtkLookupTableImageFilter_h
20 #define rtkLookupTableImageFilter_h
21 
24 
25 namespace rtk
26 {
27 
43 namespace Functor
44 {
45 
46 template <class TInput, class TOutput>
47 class ITK_TEMPLATE_EXPORT LUT
48 {
49 public:
54  using InterpolatorPointer = typename InterpolatorType::Pointer;
55 
56  LUT()
57  : m_LookupTableDataPointer(nullptr)
58  , m_Interpolator(InterpolatorType::New()){};
59  ~LUT() = default;
60 
62  LookupTablePointer
64  {
65  return m_LookupTablePointer;
66  }
67  void
69  {
70  m_LookupTablePointer = lut;
71  m_LookupTableDataPointer = lut->GetBufferPointer();
72  m_InverseLUTSpacing = 1. / m_LookupTablePointer->GetSpacing()[0];
73  m_Interpolator->SetInputImage(lut);
74  }
76 
77  LookupTableDataPointer
79  {
80  return m_LookupTableDataPointer;
81  }
82 
83  bool
84  operator!=(const LUT & lut) const
85  {
86  return m_LookupTableDataPointer != lut.GetLookupTable();
87  }
88  bool
89  operator==(const LUT & lut) const
90  {
91  return m_LookupTableDataPointer == lut.GetLookupTable();
92  }
93 
94  inline TOutput
95  operator()(const TInput & val) const;
96 
97 private:
102 };
103 
104 template <class TInput, class TOutput>
105 TOutput
106 LUT<TInput, TOutput>::operator()(const TInput & val) const
107 {
108  return m_LookupTableDataPointer[val];
109 }
110 
111 template <>
112 inline float
113 LUT<float, float>::operator()(const float & val) const
114 {
115  InterpolatorType::ContinuousIndexType index;
116  index[0] = m_InverseLUTSpacing * (val - m_LookupTablePointer->GetOrigin()[0]);
117  return float(m_Interpolator->EvaluateAtContinuousIndex(index));
118 }
119 
120 template <>
121 inline double
122 LUT<double, double>::operator()(const double & val) const
123 {
124  InterpolatorType::ContinuousIndexType index;
125  index[0] = m_InverseLUTSpacing * (val - m_LookupTablePointer->GetOrigin()[0]);
126  return double(m_Interpolator->EvaluateAtContinuousIndex(index));
127 }
128 
129 } // end namespace Functor
130 
144 template <class TInputImage, class TOutputImage>
145 class ITK_TEMPLATE_EXPORT LookupTableImageFilter
146  : public itk::UnaryFunctorImageFilter<TInputImage,
147  TOutputImage,
148  Functor::LUT<typename TInputImage::PixelType, typename TOutputImage::PixelType>>
149 {
150 public:
151  ITK_DISALLOW_COPY_AND_MOVE(LookupTableImageFilter);
152 
156 
162 
164  itkNewMacro(Self);
165 
167 #ifdef itkOverrideGetNameOfClassMacro
168  itkOverrideGetNameOfClassMacro(LookupTableImageFilter);
169 #else
171 #endif
172 
173 
175  virtual void
177  {
178  // Idem as itkSetObjectMacro + call to functor SetLookupTableDataPointer
179  itkDebugMacro("setting "
180  << "LookupTable"
181  " to "
182  << _arg);
183  if (this->m_LookupTable != _arg || (_arg && _arg->GetTimeStamp() > this->GetTimeStamp()))
184  {
185  this->m_LookupTable = _arg;
186  this->Modified();
187  this->GetFunctor().SetLookupTable(_arg);
188  }
189  }
191 
193  itkGetModifiableObjectMacro(LookupTable, LookupTableType);
194 
197  void
198  BeforeThreadedGenerateData() override;
199 
200 protected:
201  LookupTableImageFilter() = default;
202  ~LookupTableImageFilter() override = default;
203  typename LookupTableType::Pointer m_LookupTable;
204 };
205 
206 } // end namespace rtk
207 
208 #ifndef ITK_MANUAL_INSTANTIATION
209 # include "rtkLookupTableImageFilter.hxx"
210 #endif
211 
212 #endif
TPixel PixelType
typename itk::LinearInterpolateImageFunction< LookupTableType, double > InterpolatorType
LookupTablePointer GetLookupTable()
void SetLookupTable(LookupTablePointer lut)
virtual void SetLookupTable(LookupTableType *_arg)
Function to do the lookup operation.
LookupTableDataPointer m_LookupTableDataPointer
TOutput operator()(const TInput &val) const
typename LookupTableType::PixelType * LookupTableDataPointer
LookupTableType::Pointer m_LookupTable
LookupTablePointer m_LookupTablePointer
InterpolatorPointer m_Interpolator
typename InterpolatorType::Pointer InterpolatorPointer
bool operator==(const LUT &lut) const
typename FunctorType::LookupTableType LookupTableType
LookupTableDataPointer GetLookupTable() const
typename LookupTableType::Pointer LookupTablePointer
bool operator!=(const LUT &lut) const