ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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{
72
73template <typename TInputImage, typename TOutputImage>
74class ITK_TEMPLATE_EXPORT BilateralImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
75{
76public:
77 ITK_DISALLOW_COPY_AND_MOVE(BilateralImageFilter);
78
84
86 itkNewMacro(Self);
87
89 itkOverrideGetNameOfClassMacro(BilateralImageFilter);
90
92 using InputImageType = TInputImage;
93 using OutputImageType = TOutputImage;
94
97
100 using OutputPixelType = typename TOutputImage::PixelType;
101 using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
103 using InputPixelType = typename TInputImage::PixelType;
104 using InputInternalPixelType = typename TInputImage::InternalPixelType;
105
108 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
109
112
115
120
124
127
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);
143 void
144 SetDomainSigma(const double v)
145 {
146 m_DomainSigma.Fill(v);
147 }
148
155 itkBooleanMacro(AutomaticKernelSize);
156 itkGetConstMacro(AutomaticKernelSize, bool);
157 itkSetMacro(AutomaticKernelSize, bool);
161 void
163
164 itkSetMacro(Radius, SizeType);
165 itkGetConstReferenceMacro(Radius, SizeType);
166
171 itkSetMacro(NumberOfRangeGaussianSamples, unsigned long);
172 itkGetConstMacro(NumberOfRangeGaussianSamples, unsigned long);
174 itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<OutputPixelType>));
175
176protected:
179
181 ~BilateralImageFilter() override = default;
182
184 void
185 PrintSelf(std::ostream & os, Indent indent) const override;
186
188 void
190
193 void
194 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
195
202 void
204
205private:
208 double m_RangeSigma{};
209
213
216 double m_DomainMu{};
217 double m_RangeMu{};
218
221
226
231 std::vector<double> m_RangeGaussianTable{};
232};
233} // end namespace itk
234
235#ifndef ITK_MANUAL_INSTANTIATION
236# include "itkBilateralImageFilter.hxx"
237#endif
238
239#endif
typename TInputImage::PixelType InputPixelType
FixedArray< double, Self::ImageDimension > ArrayType
void SetRadius(const SizeValueType)
typename OutputImageType::RegionType OutputImageRegionType
typename TOutputImage::PixelType OutputPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
Neighborhood< double, Self::ImageDimension > KernelType
typename KernelType::SizeType SizeType
SmartPointer< const Self > ConstPointer
~BilateralImageFilter() override=default
void BeforeThreadedGenerateData() override
Image< double, Self::ImageDimension > GaussianImageType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
typename KernelType::ConstIterator KernelConstIteratorType
std::vector< double > m_RangeGaussianTable
void GenerateInputRequestedRegion() override
typename KernelType::SizeValueType SizeValueType
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
typename TOutputImage::InternalPixelType OutputInternalPixelType
static constexpr unsigned int ImageDimension
typename NumericTraits< OutputPixelType >::RealType OutputPixelRealType
void PrintSelf(std::ostream &os, Indent indent) const override
ConstNeighborhoodIterator< TInputImage > NeighborhoodIteratorType
typename KernelType::Iterator KernelIteratorType
Const version of NeighborhoodIterator, defining iteration of a local N-dimensional neighborhood of pi...
Simulate a standard C array with copy semantics.
typename OutputImageType::RegionType OutputImageRegionType
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
A light-weight container object for storing an N-dimensional neighborhood of values.
typename AllocatorType::const_iterator ConstIterator
typename SizeType::SizeValueType SizeValueType
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....