ITK  4.6.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 >
76  public ImageToImageFilter< TInputImage, TOutputImage >
77 {
78 public:
84 
86  itkNewMacro(Self);
87 
90 
92  typedef TInputImage InputImageType;
93  typedef TOutputImage OutputImageType;
94 
97 
100  typedef typename TOutputImage::PixelType OutputPixelType;
101  typedef typename TOutputImage::InternalPixelType OutputInternalPixelType;
103  typedef typename TInputImage::PixelType InputPixelType;
104  typedef typename TInputImage::InternalPixelType InputInternalPixelType;
105 
108  itkStaticConstMacro(ImageDimension, unsigned int,
109  TOutputImage::ImageDimension);
110 
113 
116 
118  typedef
120  typedef typename KernelType::SizeType SizeType;
122 
126 
128  typedef
130 
134  itkSetMacro(DomainSigma, ArrayType);
135  itkGetConstMacro(DomainSigma, const ArrayType);
136  itkSetMacro(DomainMu, double);
137  itkGetConstReferenceMacro(DomainMu, double);
138  itkSetMacro(RangeSigma, double);
139  itkGetConstMacro(RangeSigma, double);
140  itkGetConstMacro(FilterDimensionality, unsigned int);
141  itkSetMacro(FilterDimensionality, unsigned int);
143 
146  void SetDomainSigma(const double v)
147  {
148  m_DomainSigma.Fill(v);
149  }
150 
156  itkBooleanMacro(AutomaticKernelSize);
157  itkGetConstMacro(AutomaticKernelSize, bool);
158  itkSetMacro(AutomaticKernelSize, bool);
160 
163  void SetRadius(const SizeValueType);
164 
165  itkSetMacro(Radius, SizeType);
166  itkGetConstReferenceMacro(Radius, SizeType);
167 
171  itkSetMacro(NumberOfRangeGaussianSamples, unsigned long);
172  itkGetConstMacro(NumberOfRangeGaussianSamples, unsigned long);
174 
175 #ifdef ITK_USE_CONCEPT_CHECKING
176  // Begin concept checking
177  itkConceptMacro( OutputHasNumericTraitsCheck,
179  // End concept checking
180 #endif
181 
182 protected:
185 
188 
190  void PrintSelf(std::ostream & os, Indent indent) const;
191 
194 
197  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
198  ThreadIdType threadId);
199 
206  virtual void GenerateInputRequestedRegion();
207 
208 private:
209  BilateralImageFilter(const Self &); //purposely not implemented
210  void operator=(const Self &); //purposely not implemented
211 
214  double m_RangeSigma;
215 
219 
222  double m_DomainMu;
223  double m_RangeMu;
224 
227 
232 
237  std::vector< double > m_RangeGaussianTable;
238 };
239 } // end namespace itk
240 
241 #ifndef ITK_MANUAL_INSTANTIATION
242 #include "itkBilateralImageFilter.hxx"
243 #endif
244 
245 #endif
SmartPointer< const Self > ConstPointer
Neighborhood< double, itkGetStaticConstMacro(ImageDimension) > KernelType
ConstNeighborhoodIterator< TInputImage > NeighborhoodIteratorType
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId)
TOutputImage::PixelType OutputPixelType
void PrintSelf(std::ostream &os, Indent indent) const
void SetRadius(const SizeValueType)
Base class for all process objects that output image data.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
static const unsigned int ImageDimension
void Fill(const ValueType &)
KernelType::ConstIterator KernelConstIteratorType
void SetDomainSigma(const double v)
std::vector< double > m_RangeGaussianTable
FixedArray< double, itkGetStaticConstMacro(ImageDimension) > ArrayType
void operator=(const Self &)
TInputImage::PixelType InputPixelType
Blurs an image while preserving edges.
NumericTraits< OutputPixelType >::RealType OutputPixelRealType
TOutputImage::InternalPixelType OutputInternalPixelType
TInputImage::InternalPixelType InputInternalPixelType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
KernelType::Iterator KernelIteratorType
virtual void GenerateInputRequestedRegion()
Superclass::OutputImageRegionType OutputImageRegionType
KernelType::SizeValueType SizeValueType
Define additional traits for native types such as int or float.
#define itkConceptMacro(name, concept)
Image< double, itkGetStaticConstMacro(ImageDimension) > GaussianImageType
Templated n-dimensional image class.
Definition: itkImage.h:75
Superclass::OutputImageRegionType OutputImageRegionType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
ImageToImageFilter< TInputImage, TOutputImage > Superclass