File: itkMinMaxCurvatureFlowImageFilter.h

package info (click to toggle)
insighttoolkit5 5.4.4-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 704,404 kB
  • sloc: cpp: 783,697; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 461; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (155 lines) | stat: -rw-r--r-- 6,422 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkMinMaxCurvatureFlowImageFilter_h
#define itkMinMaxCurvatureFlowImageFilter_h

#include "itkCurvatureFlowImageFilter.h"
#include "itkMinMaxCurvatureFlowFunction.h"

namespace itk
{
/**
 * \class MinMaxCurvatureFlowImageFilter
 * \brief Denoise an image using min/max curvature flow.
 *
 * MinMaxCurvatureFlowImageFilter implements a curvature driven image denoising
 * algorithm. Iso-brightness contours in the grayscale input image are viewed
 * as a level set. The level set is then evolved using a curvature-based speed
 * function:
 *
 * \f[  I_t = F_{\mbox{minmax}} |\nabla I| \f]
 *
 * where \f$ F_{\mbox{minmax}} = \max(\kappa,0) \f$ if
 * \f$ \mbox{Avg}_{\mbox{stencil}}(x) \f$ is less than or equal to \f$ T_{threshold} \f$
 * and \f$ \min(\kappa,0) \f$, otherwise.
 * \f$ \kappa \f$ is the mean curvature of the iso-brightness contour
 * at point \f$ x \f$.
 *
 * In min/max curvature flow, movement is turned on or off depending
 * on the scale of the noise one wants to remove. Switching depends on
 * the average image value of a region of radius \f$ R \f$ around each
 * point. The choice of \f$ R \f$, the stencil radius, governs the scale of
 * the noise to be removed.
 *
 * The threshold value \f$ T_{threshold} \f$ is the average intensity obtained
 * in the direction perpendicular to the gradient at point \f$ x \f$ at the
 * extrema of the local neighborhood.
 *
 * This filter make use of the multi-threaded finite difference solver
 * hierarchy.  Updates are computed using a MinMaxCurvatureFlowFunction object.
 * A zero flux Neumann boundary condition is used when computing derivatives
 * near the data boundary.
 *
 * \warning This filter assumes that the input and output types have the
 * same dimensions. This filter also requires that the output image pixels
 * are of a real type. This filter works for any dimensional images,
 * however for dimensions greater than 3D, an expensive brute-force search
 * is used to compute the local threshold.
 *
 * Reference:
 * "Level Set Methods and Fast Marching Methods", J.A. Sethian,
 * Cambridge Press, Chapter 16, Second edition, 1999.
 *
 * \sa MinMaxCurvatureFlowFunction
 * \sa CurvatureFlowImageFilter
 * \sa BinaryMinMaxCurvatureFlowImageFilter
 *
 * \ingroup ImageEnhancement
 * \ingroup MultiThreaded
 *
 * \ingroup ITKCurvatureFlow
 *
 * \sphinx
 * \sphinxexample{Filtering/CurvatureFlow/SmoothImageUsingMinMaxCurvatureFlow,Smooth Image Using Min Max Curvature Flow}
 * \sphinxexample{Filtering/CurvatureFlow/SmoothRGBImageUsingMinMaxCurvatureFlow,SmoothRGBImageUsingMinMaxCurvatureFlow}
 * \endsphinx
 */
template <typename TInputImage, typename TOutputImage>
class ITK_TEMPLATE_EXPORT MinMaxCurvatureFlowImageFilter : public CurvatureFlowImageFilter<TInputImage, TOutputImage>
{
public:
  ITK_DISALLOW_COPY_AND_MOVE(MinMaxCurvatureFlowImageFilter);

  /** Standard class type aliases. */
  using Self = MinMaxCurvatureFlowImageFilter;
  using Superclass = CurvatureFlowImageFilter<TInputImage, TOutputImage>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;

  /** Method for creation through the object factory. */
  itkNewMacro(Self);

  /** \see LightObject::GetNameOfClass() */
  itkOverrideGetNameOfClassMacro(MinMaxCurvatureFlowImageFilter);

  /** Inherit type alias from Superclass. */
  using typename Superclass::FiniteDifferenceFunctionType;
  using typename Superclass::OutputImageType;

  /** MinMaxCurvatureFlowFunction type. */
  using MinMaxCurvatureFlowFunctionType = MinMaxCurvatureFlowFunction<OutputImageType>;

  /** Dimensionality of input and output data is assumed to be the same.
   * It is inherited from the superclass. */
  static constexpr unsigned int ImageDimension = Superclass::ImageDimension;

  /** Typedef support for the neighbour radius. */
  using RadiusType = typename FiniteDifferenceFunctionType::RadiusType;
  using RadiusValueType = typename RadiusType::SizeValueType;

  /** Set/Get the stencil radius. */
  itkSetMacro(StencilRadius, RadiusValueType);
  itkGetConstMacro(StencilRadius, RadiusValueType);

#ifdef ITK_USE_CONCEPT_CHECKING
  // Begin concept checking
  itkConceptMacro(UnsignedLongConvertibleToOutputCheck,
                  (Concept::Convertible<unsigned long, typename TOutputImage::PixelType>));
  itkConceptMacro(OutputLessThanComparableCheck, (Concept::LessThanComparable<typename TOutputImage::PixelType>));
  itkConceptMacro(LongConvertibleToOutputCheck, (Concept::Convertible<long, typename TOutputImage::PixelType>));
  itkConceptMacro(OutputDoubleComparableCheck, (Concept::Comparable<typename TOutputImage::PixelType, double>));
  itkConceptMacro(OutputDoubleMultiplyAndAssignOperatorCheck,
                  (Concept::MultiplyAndAssignOperator<typename TOutputImage::PixelType, double>));
  itkConceptMacro(OutputGreaterThanUnsignedLongCheck,
                  (Concept::GreaterThanComparable<typename TOutputImage::PixelType, unsigned long>));
  itkConceptMacro(UnsignedLongOutputAditiveOperatorsCheck,
                  (Concept::AdditiveOperators<unsigned long, typename TOutputImage::PixelType>));
  // End concept checking
#endif

protected:
  MinMaxCurvatureFlowImageFilter();
  ~MinMaxCurvatureFlowImageFilter() override = default;
  void
  PrintSelf(std::ostream & os, Indent indent) const override;

  /** Initialize the state of filter and equation before each iteration.
   * Progress feedback is implemented as part of this method. */
  void
  InitializeIteration() override;

private:
  RadiusValueType m_StencilRadius{};
};
} // namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkMinMaxCurvatureFlowImageFilter.hxx"
#endif

#endif