ITK  6.0.0
Insight Toolkit
itkBilateralImageFilter.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 itkBilateralImageFilter_h
19#define itkBilateralImageFilter_h
20
22#include "itkFixedArray.h"
24#include "itkNeighborhood.h"
25
26namespace itk
27{
74template <typename TInputImage, typename TOutputImage>
75class ITK_TEMPLATE_EXPORT BilateralImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
76{
77public:
78 ITK_DISALLOW_COPY_AND_MOVE(BilateralImageFilter);
79
85
87 itkNewMacro(Self);
88
90 itkOverrideGetNameOfClassMacro(BilateralImageFilter);
91
93 using InputImageType = TInputImage;
94 using OutputImageType = TOutputImage;
95
97 using typename Superclass::OutputImageRegionType;
98
101 using OutputPixelType = typename TOutputImage::PixelType;
102 using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
104 using InputPixelType = typename TInputImage::PixelType;
105 using InputInternalPixelType = typename TInputImage::InternalPixelType;
106
109 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
110
113
116
121
125
128
132 itkSetMacro(DomainSigma, ArrayType);
133 itkGetConstMacro(DomainSigma, const ArrayType);
134 itkSetMacro(DomainMu, double);
135 itkGetConstReferenceMacro(DomainMu, double);
136 itkSetMacro(RangeSigma, double);
137 itkGetConstMacro(RangeSigma, double);
138 itkGetConstMacro(FilterDimensionality, unsigned int);
139 itkSetMacro(FilterDimensionality, unsigned int);
144 void
145 SetDomainSigma(const double v)
146 {
147 m_DomainSigma.Fill(v);
148 }
149
155 itkBooleanMacro(AutomaticKernelSize);
156 itkGetConstMacro(AutomaticKernelSize, bool);
157 itkSetMacro(AutomaticKernelSize, bool);
162 void
164
165 itkSetMacro(Radius, SizeType);
166 itkGetConstReferenceMacro(Radius, SizeType);
167
171 itkSetMacro(NumberOfRangeGaussianSamples, unsigned long);
172 itkGetConstMacro(NumberOfRangeGaussianSamples, unsigned long);
175#ifdef ITK_USE_CONCEPT_CHECKING
176 // Begin concept checking
177 itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<OutputPixelType>));
178 // End concept checking
179#endif
180
181protected:
184
186 ~BilateralImageFilter() override = default;
187
189 void
190 PrintSelf(std::ostream & os, Indent indent) const override;
191
193 void
195
198 void
199 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
200
207 void
209
210private:
213 double m_RangeSigma{};
214
217 ArrayType m_DomainSigma{};
218
221 double m_DomainMu{};
222 double m_RangeMu{};
223
225 unsigned int m_FilterDimensionality{};
226
228 KernelType m_GaussianKernel{};
229 SizeType m_Radius{};
230 bool m_AutomaticKernelSize{};
231
233 unsigned long m_NumberOfRangeGaussianSamples{};
234 double m_DynamicRange{};
235 double m_DynamicRangeUsed{};
236 std::vector<double> m_RangeGaussianTable{};
237};
238} // end namespace itk
239
240#ifndef ITK_MANUAL_INSTANTIATION
241# include "itkBilateralImageFilter.hxx"
242#endif
243
244#endif
Blurs an image while preserving edges.
typename TInputImage::PixelType InputPixelType
void SetDomainSigma(const double v)
void SetRadius(const SizeValueType)
typename TOutputImage::PixelType OutputPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
typename KernelType::SizeType SizeType
~BilateralImageFilter() override=default
void BeforeThreadedGenerateData() override
typename KernelType::ConstIterator KernelConstIteratorType
void GenerateInputRequestedRegion() override
typename KernelType::SizeValueType SizeValueType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
typename TOutputImage::InternalPixelType OutputInternalPixelType
typename NumericTraits< OutputPixelType >::RealType OutputPixelRealType
void PrintSelf(std::ostream &os, Indent indent) const override
typename KernelType::Iterator KernelIteratorType
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
Base class for filters that take an image as input and produce an image as output.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
typename AllocatorType::iterator Iterator
typename AllocatorType::const_iterator ConstIterator
Define additional traits for native types such as int or float.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86