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 
21 #include "itkImageToImageFilter.h"
22 #include "itkFixedArray.h"
24 #include "itkNeighborhood.h"
25 
26 namespace itk
27 {
74 template <typename TInputImage, typename TOutputImage>
75 class ITK_TEMPLATE_EXPORT BilateralImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
76 {
77 public:
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 
119  using SizeType = typename KernelType::SizeType;
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
163  SetRadius(const SizeValueType);
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 
181 protected:
184 
186  ~BilateralImageFilter() override = default;
187 
189  void
190  PrintSelf(std::ostream & os, Indent indent) const override;
191 
193  void
194  BeforeThreadedGenerateData() override;
195 
198  void
199  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
200 
207  void
208  GenerateInputRequestedRegion() override;
209 
210 private:
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
itk::BilateralImageFilter::OutputInternalPixelType
typename TOutputImage::InternalPixelType OutputInternalPixelType
Definition: itkBilateralImageFilter.h:102
itk::BilateralImageFilter
Blurs an image while preserving edges.
Definition: itkBilateralImageFilter.h:75
itk::BilateralImageFilter::KernelIteratorType
typename KernelType::Iterator KernelIteratorType
Definition: itkBilateralImageFilter.h:123
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:716
itkNeighborhoodIterator.h
itk::BilateralImageFilter::InputInternalPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
Definition: itkBilateralImageFilter.h:105
itk::Neighborhood< double, Self::ImageDimension >
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::BilateralImageFilter::KernelConstIteratorType
typename KernelType::ConstIterator KernelConstIteratorType
Definition: itkBilateralImageFilter.h:124
itkNeighborhood.h
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::BilateralImageFilter::SizeType
typename KernelType::SizeType SizeType
Definition: itkBilateralImageFilter.h:119
itk::BilateralImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkBilateralImageFilter.h:104
itk::Neighborhood< double, Self::ImageDimension >::Iterator
typename AllocatorType::iterator Iterator
Definition: itkNeighborhood.h:75
itkFixedArray.h
itk::BilateralImageFilter::SetDomainSigma
void SetDomainSigma(const double v)
Definition: itkBilateralImageFilter.h:145
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::BilateralImageFilter::SizeValueType
typename KernelType::SizeValueType SizeValueType
Definition: itkBilateralImageFilter.h:120
itk::FixedArray< double, Self::ImageDimension >
itk::Neighborhood< double, Self::ImageDimension >::ConstIterator
typename AllocatorType::const_iterator ConstIterator
Definition: itkNeighborhood.h:76
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::BilateralImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkBilateralImageFilter.h:101
itk::BilateralImageFilter::OutputPixelRealType
typename NumericTraits< OutputPixelType >::RealType OutputPixelRealType
Definition: itkBilateralImageFilter.h:103
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ConstNeighborhoodIterator< TInputImage >
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:88
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90