ITK  6.0.0
Insight Toolkit
itkWindowedSincInterpolateImageFunction.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 itkWindowedSincInterpolateImageFunction_h
19#define itkWindowedSincInterpolateImageFunction_h
20
24#include "itkMath.h"
25
26namespace itk
27{
28// clang-format off
29namespace Function
30{
38template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
39class ITK_TEMPLATE_EXPORT CosineWindowFunction
40{
41public:
42 inline TOutput
43 operator()(const TInput & A) const
44 {
45
47 static constexpr double factor = Math::pi / (2 * VRadius);
48 return static_cast<TOutput>(std::cos(A * factor));
49 }
50};
60template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
61class ITK_TEMPLATE_EXPORT HammingWindowFunction
62{
63public:
64 inline TOutput
65 operator()(const TInput & A) const
66 {
67
69 static constexpr double factor = Math::pi / VRadius;
70 return static_cast<TOutput>(0.54 + 0.46 * std::cos(A * factor));
71 }
72};
73
81template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
82class ITK_TEMPLATE_EXPORT WelchWindowFunction
83{
84public:
85 inline TOutput
86 operator()(const TInput & A) const
87 {
88
90 static constexpr double factor = 1.0 / (VRadius * VRadius);
91 return static_cast<TOutput>(1.0 - A * factor * A);
92 }
93};
105template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
106class ITK_TEMPLATE_EXPORT LanczosWindowFunction
107{
108public:
109 inline TOutput
110 operator()(const TInput & A) const
111 {
112 if (A == 0.0)
113 {
114 return static_cast<TOutput>(1.0);
115 } // namespace Function
116
118 static constexpr double factor = Math::pi / VRadius;
119 const double z = factor * A;
120 return static_cast<TOutput>(std::sin(z) / z);
121 } // namespace itk
122};
123
131template <unsigned int VRadius, typename TInput = double, typename TOutput = double>
132class ITK_TEMPLATE_EXPORT BlackmanWindowFunction
133{
134public:
135 inline TOutput
136 operator()(const TInput & A) const
137 {
138
140 static constexpr double factor1 = Math::pi / VRadius;
141
143 static constexpr double factor2 = 2.0 * Math::pi / VRadius;
144 return static_cast<TOutput>(0.42 + 0.5 * std::cos(A * factor1) + 0.08 * std::cos(A * factor2));
145 }
146};
147} // namespace Function
148// clang-format on
149
260template <typename TInputImage,
261 unsigned int VRadius,
262 typename TWindowFunction = Function::HammingWindowFunction<VRadius>,
263 class TBoundaryCondition = ZeroFluxNeumannBoundaryCondition<TInputImage, TInputImage>,
264 class TCoordinate = double>
266 : public InterpolateImageFunction<TInputImage, TCoordinate>
267{
268public:
269 ITK_DISALLOW_COPY_AND_MOVE(WindowedSincInterpolateImageFunction);
275
278
280 itkOverrideGetNameOfClassMacro(WindowedSincInterpolateImageFunction);
281
283 itkNewMacro(Self);
284
286 using typename Superclass::OutputType;
287
289 using typename Superclass::InputImageType;
290
292 using typename Superclass::RealType;
293
295 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
296
298 using typename Superclass::IndexType;
299 using typename Superclass::IndexValueType;
300
302 using typename Superclass::SizeType;
303
305 using ImageType = TInputImage;
306
308 using typename Superclass::ContinuousIndexType;
309
310 void
311 SetInputImage(const ImageType * image) override;
312
320 EvaluateAtContinuousIndex(const ContinuousIndexType & index) const override;
321
323 GetRadius() const override
324 {
325 constexpr auto radius = SizeType::Filled(VRadius);
326 return radius;
327 }
328
329protected:
332 void
333 PrintSelf(std::ostream & os, Indent indent) const override;
334
335private:
336 // Internal type alias
338
339 // Constant to store twice the radius
340 static constexpr unsigned int m_WindowSize{ 2 * VRadius };
341
343 TWindowFunction m_WindowFunction{};
344
346 static constexpr unsigned int m_OffsetTableSize = Math::UnsignedPower(m_WindowSize, ImageDimension);
347
350 unsigned int m_OffsetTable[m_OffsetTableSize]{};
351
353 unsigned int m_WeightOffsetTable[m_OffsetTableSize][ImageDimension]{};
354
356 static double
357 Sinc(const double x)
358 {
359 const double px = Math::pi * x;
360 return (x == 0.0) ? 1.0 : std::sin(px) / px;
361 }
362};
363} // namespace itk
366#ifndef ITK_MANUAL_INSTANTIATION
367# include "itkWindowedSincInterpolateImageFunction.hxx"
368#endif
369
370#endif // _itkWindowedSincInterpolateImageFunction_h
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Base class for all image interpolators.
Light weight base class for most itk classes.
Use the windowed sinc function to interpolate.
~WindowedSincInterpolateImageFunction() override=default
void PrintSelf(std::ostream &os, Indent indent) const override
OutputType EvaluateAtContinuousIndex(const ContinuousIndexType &index) const override
void SetInputImage(const ImageType *image) override
constexpr TReturnType UnsignedPower(const uintmax_t base, const uintmax_t exponent) noexcept
Definition: itkMath.h:793
static constexpr double pi
Definition: itkMath.h:66
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
long IndexValueType
Definition: itkIntTypes.h:93