ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkAdaptiveNonLocalMeansDenoisingImageFilter.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 itkAdaptiveNonLocalMeansDenoisingImageFilter_h
19#define itkAdaptiveNonLocalMeansDenoisingImageFilter_h
20
22
24#include "itkGaussianOperator.h"
25
26namespace itk
27{
28
46
47template <typename TInputImage,
48 typename TOutputImage = TInputImage,
50class ITK_TEMPLATE_EXPORT AdaptiveNonLocalMeansDenoisingImageFilter final
51 : public NonLocalPatchBasedImageFilter<TInputImage, TOutputImage>
52{
53public:
54 ITK_DISALLOW_COPY_AND_MOVE(AdaptiveNonLocalMeansDenoisingImageFilter);
55
61
63 itkOverrideGetNameOfClassMacro(AdaptiveNonLocalMeansDenoisingImageFilter);
64
66 itkNewMacro(Self);
67
69 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
70
72 using InputImageType = TInputImage;
73 using InputPixelType = typename InputImageType::PixelType;
74 using OutputImageType = TOutputImage;
76
77 using MaskImageType = TMaskImage;
78 using MaskPixelType = typename MaskImageType::PixelType;
79 using LabelType = typename MaskImageType::PixelType;
80
85
90
92
96 void
98 {
99 this->SetInput(image);
100 }
101
106 void
108 {
109 this->SetNthInput(1, const_cast<MaskImageType *>(mask));
110 }
111 void
113 {
114 this->SetMaskImage(mask);
115 }
116
121 const MaskImageType *
123 {
124 return static_cast<const MaskImageType *>(this->ProcessObject::GetInput(1));
125 }
126
131 itkNonVirtualSetMacro(UseRicianNoiseModel, bool);
132 itkNonVirtualGetConstMacro(UseRicianNoiseModel, bool);
133 itkNonVirtualBooleanMacro(UseRicianNoiseModel);
134
140
144 itkNonVirtualSetMacro(SmoothingVariance, RealType);
146
153
160
165 itkNonVirtualSetMacro(VarianceThreshold, RealType);
167
172 itkNonVirtualSetMacro(NeighborhoodRadiusForLocalMeanAndVariance, NeighborhoodRadiusType);
173 itkNonVirtualGetConstMacro(NeighborhoodRadiusForLocalMeanAndVariance, NeighborhoodRadiusType);
174
175protected:
178
179 void
180 PrintSelf(std::ostream & os, Indent indent) const override;
181
182 void
184
185 void
187
188 void
190
191private:
193
195
197
203
206
212
214};
215
216} // end namespace itk
217
218#ifndef ITK_MANUAL_INSTANTIATION
219# include "itkAdaptiveNonLocalMeansDenoisingImageFilter.hxx"
220#endif
221
222#endif
itkNonVirtualSetMacro(UseRicianNoiseModel, bool)
itkNonVirtualGetConstMacro(NeighborhoodRadiusForLocalMeanAndVariance, NeighborhoodRadiusType)
itkNonVirtualSetMacro(VarianceThreshold, RealType)
itkNonVirtualSetMacro(SmoothingVariance, RealType)
itkNonVirtualGetConstMacro(VarianceThreshold, RealType)
void ThreadedGenerateData(const RegionType &, ThreadIdType) override
void PrintSelf(std::ostream &os, Indent indent) const override
itkNonVirtualSetMacro(NeighborhoodRadiusForLocalMeanAndVariance, NeighborhoodRadiusType)
itkNonVirtualSetMacro(MeanThreshold, RealType)
typename Superclass::ConstNeighborhoodIteratorType ConstNeighborhoodIteratorType
itkNonVirtualGetConstMacro(SmoothingVariance, RealType)
itkNonVirtualSetMacro(SmoothingFactor, RealType)
NonLocalPatchBasedImageFilter< TInputImage, TOutputImage > Superclass
itkNonVirtualGetConstMacro(SmoothingFactor, RealType)
itkNonVirtualGetConstMacro(UseRicianNoiseModel, bool)
itkNonVirtualGetConstMacro(MeanThreshold, RealType)
~AdaptiveNonLocalMeansDenoisingImageFilter() override=default
typename Superclass::NeighborhoodOffsetListType NeighborhoodOffsetListType
A NeighborhoodOperator whose coefficients are a one dimensional, discrete Gaussian kernel.
virtual void SetInput(const InputImageType *input)
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:51
typename ConstNeighborhoodIteratorType::RadiusType NeighborhoodRadiusType
typename ConstNeighborhoodIteratorType::OffsetType NeighborhoodOffsetType
virtual void SetNthInput(DataObjectPointerArraySizeType idx, DataObject *input)
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType