ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkBilateralImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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:
76  public ImageToImageFilter< TInputImage, TOutputImage >
77 {
78 public:
79  ITK_DISALLOW_COPY_AND_ASSIGN(BilateralImageFilter);
80 
86 
88  itkNewMacro(Self);
89 
92 
94  using InputImageType = TInputImage;
95  using OutputImageType = TOutputImage;
96 
98  using OutputImageRegionType = typename Superclass::OutputImageRegionType;
99 
102  using OutputPixelType = typename TOutputImage::PixelType;
103  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
105  using InputPixelType = typename TInputImage::PixelType;
106  using InputInternalPixelType = typename TInputImage::InternalPixelType;
107 
110  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
111 
114 
117 
120  using SizeType = typename KernelType::SizeType;
122 
126 
129 
133  itkSetMacro(DomainSigma, ArrayType);
134  itkGetConstMacro(DomainSigma, const ArrayType);
135  itkSetMacro(DomainMu, double);
136  itkGetConstReferenceMacro(DomainMu, double);
137  itkSetMacro(RangeSigma, double);
138  itkGetConstMacro(RangeSigma, double);
139  itkGetConstMacro(FilterDimensionality, unsigned int);
140  itkSetMacro(FilterDimensionality, unsigned int);
142 
145  void SetDomainSigma(const double v)
146  {
147  m_DomainSigma.Fill(v);
148  }
149 
155  itkBooleanMacro(AutomaticKernelSize);
156  itkGetConstMacro(AutomaticKernelSize, bool);
157  itkSetMacro(AutomaticKernelSize, bool);
159 
162  void SetRadius(const SizeValueType);
163 
164  itkSetMacro(Radius, SizeType);
165  itkGetConstReferenceMacro(Radius, SizeType);
166 
170  itkSetMacro(NumberOfRangeGaussianSamples, unsigned long);
171  itkGetConstMacro(NumberOfRangeGaussianSamples, unsigned long);
173 
174 #ifdef ITK_USE_CONCEPT_CHECKING
175  // Begin concept checking
176  itkConceptMacro( OutputHasNumericTraitsCheck,
178  // End concept checking
179 #endif
180 
181 protected:
184 
186  ~BilateralImageFilter() override = default;
187 
189  void PrintSelf(std::ostream & os, Indent indent) const override;
190 
192  void BeforeThreadedGenerateData() override;
193 
196  void DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
197 
204  void GenerateInputRequestedRegion() override;
205 
206 private:
209  double m_RangeSigma;
210 
214 
217  double m_DomainMu;
218  double m_RangeMu;
219 
222 
227 
232  std::vector< double > m_RangeGaussianTable;
233 };
234 } // end namespace itk
235 
236 #ifndef ITK_MANUAL_INSTANTIATION
237 #include "itkBilateralImageFilter.hxx"
238 #endif
239 
240 #endif
Define numeric traits for std::vector.
unsigned long SizeValueType
Definition: itkIntTypes.h:83
typename TInputImage::InternalPixelType InputInternalPixelType
typename AllocatorType::iterator Iterator
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename TOutputImage::PixelType OutputPixelType
typename AllocatorType::const_iterator ConstIterator
typename KernelType::ConstIterator KernelConstIteratorType
void SetDomainSigma(const double v)
std::vector< double > m_RangeGaussianTable
typename TOutputImage::InternalPixelType OutputInternalPixelType
typename OutputImageType::RegionType OutputImageRegionType
typename KernelType::Iterator KernelIteratorType
TOutputImage OutputImageType
Blurs an image while preserving edges.
typename NumericTraits< OutputPixelType >::RealType OutputPixelRealType
typename TInputImage::PixelType InputPixelType
typename KernelType::SizeType SizeType
Base class for filters that take an image as input and produce an image as output.
typename KernelType::SizeValueType SizeValueType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
#define itkConceptMacro(name, concept)
Templated n-dimensional image class.
Definition: itkImage.h:75