ITK  5.4.0
Insight Toolkit
itkDiscreteGaussianImageFilter.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 itkDiscreteGaussianImageFilter_h
19 #define itkDiscreteGaussianImageFilter_h
20 
21 #include "itkGaussianOperator.h"
22 #include "itkImageToImageFilter.h"
23 #include "itkImage.h"
25 
26 namespace itk
27 {
63 template <typename TInputImage, typename TOutputImage = TInputImage>
64 class ITK_TEMPLATE_EXPORT DiscreteGaussianImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
65 {
66 public:
67  ITK_DISALLOW_COPY_AND_MOVE(DiscreteGaussianImageFilter);
68 
74 
76  itkNewMacro(Self);
77 
79  itkOverrideGetNameOfClassMacro(DiscreteGaussianImageFilter);
80 
82  using InputImageType = TInputImage;
83  using OutputImageType = TOutputImage;
84 
87  using OutputPixelType = typename TOutputImage::PixelType;
88  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
89  using InputPixelType = typename TInputImage::PixelType;
90  using InputInternalPixelType = typename TInputImage::InternalPixelType;
91 
95 
98  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
99 
104 
111 
115  using ScalarRealType = double;
116 
120 
127  itkSetMacro(Variance, ArrayType);
128  itkGetConstMacro(Variance, const ArrayType);
134  itkSetMacro(MaximumError, ArrayType);
135  itkGetConstMacro(MaximumError, const ArrayType);
140  itkGetConstMacro(MaximumKernelWidth, unsigned int);
141  itkSetMacro(MaximumKernelWidth, unsigned int);
149  itkGetConstMacro(FilterDimensionality, unsigned int);
150  itkSetMacro(FilterDimensionality, unsigned int);
154  itkSetMacro(InputBoundaryCondition, InputBoundaryConditionPointerType);
155  itkGetConstMacro(InputBoundaryCondition, InputBoundaryConditionPointerType);
156  itkSetMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
157  itkGetConstMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
162  void
164  {
165  m_Variance.Fill(v);
166  this->Modified();
167  }
170  void
172  {
173  m_MaximumError.Fill(v);
174  this->Modified();
175  }
176 
177  void
178  SetVariance(const double * v)
179  {
180  ArrayType dv;
181 
182  for (unsigned int i = 0; i < ImageDimension; ++i)
183  {
184  dv[i] = v[i];
185  }
186  this->SetVariance(dv);
187  }
188 
189  void
190  SetVariance(const float * v)
191  {
192  ArrayType dv;
193 
194  for (unsigned int i = 0; i < ImageDimension; ++i)
195  {
196  dv[i] = v[i];
197  }
198  this->SetVariance(dv);
199  }
200 
203  void
204  SetSigma(const ArrayType & sigma)
205  {
206  ArrayType variance;
207  for (unsigned int i = 0; i < ImageDimension; ++i)
208  {
209  variance[i] = sigma[i] * sigma[i];
210  }
211  this->SetVariance(variance);
212  }
213 
216  void
217  SetSigmaArray(const ArrayType & sigmas)
218  {
219  this->SetSigma(sigmas);
220  }
221  void
222  SetSigma(double sigma)
223  {
224  this->SetVariance(sigma * sigma);
225  }
229  ArrayType
231  {
232  ArrayType sigmas;
233  for (unsigned int i = 0; i < ImageDimension; ++i)
234  {
235  sigmas[i] = std::sqrt(m_Variance[i]);
236  }
237  return sigmas;
238  }
243  double
244  GetSigma() const
245  {
246  return std::sqrt(m_Variance[0]);
247  }
248 
249  void
250  SetMaximumError(const double * v)
251  {
252  ArrayType dv;
253 
254  for (unsigned int i = 0; i < ImageDimension; ++i)
255  {
256  dv[i] = v[i];
257  }
258  this->SetMaximumError(dv);
259  }
260 
261  void
262  SetMaximumError(const float * v)
263  {
264  ArrayType dv;
265 
266  for (unsigned int i = 0; i < ImageDimension; ++i)
267  {
268  dv[i] = v[i];
269  }
270  this->SetMaximumError(dv);
271  }
272 
274  unsigned int
275  GetKernelRadius(const unsigned int dimension) const;
276 
278  ArrayType
279  GetKernelRadius() const;
280 
283  ArrayType
284  GetKernelSize() const;
285 
291  itkSetMacro(UseImageSpacing, bool);
292  itkGetConstMacro(UseImageSpacing, bool);
293  itkBooleanMacro(UseImageSpacing);
296 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
297 
301  void
302  SetUseImageSpacingOn()
303  {
304  this->SetUseImageSpacing(true);
305  }
306 
310  void
311  SetUseImageSpacingOff()
312  {
313  this->SetUseImageSpacing(false);
314  }
315 #endif
316 
326  itkLegacyMacro(unsigned int GetInternalNumberOfStreamDivisions() const);
327  itkLegacyMacro(void SetInternalNumberOfStreamDivisions(unsigned int));
328 
329 #ifdef ITK_USE_CONCEPT_CHECKING
330  // Begin concept checking
331 
332  itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<OutputPixelValueType>));
333 
334  // End concept checking
335 #endif
336 
337 protected:
339  {
340  m_Variance.Fill(0.0);
341  m_MaximumError.Fill(0.01);
342  m_MaximumKernelWidth = 32;
343  m_UseImageSpacing = true;
344  m_FilterDimensionality = ImageDimension;
345  m_InputBoundaryCondition = &m_InputDefaultBoundaryCondition;
346  m_RealBoundaryCondition = &m_RealDefaultBoundaryCondition;
347  }
348 
349  ~DiscreteGaussianImageFilter() override = default;
350  void
351  PrintSelf(std::ostream & os, Indent indent) const override;
352 
359  void
360  GenerateInputRequestedRegion() override;
361 
367  void
368  GenerateData() override;
369 
371  void
372  GenerateKernel(const unsigned int dimension, KernelType & oper) const;
373 
375  ArrayType
376  GetKernelVarianceArray() const;
377 
378 private:
381  ArrayType m_Variance{};
382 
386  ArrayType m_MaximumError{};
387 
390  unsigned int m_MaximumKernelWidth{};
391 
393  unsigned int m_FilterDimensionality{};
394 
396  bool m_UseImageSpacing{};
397 
400  InputBoundaryConditionPointerType m_InputBoundaryCondition{};
401 
403  InputDefaultBoundaryConditionType m_InputDefaultBoundaryCondition{};
404 
406  RealBoundaryConditionPointerType m_RealBoundaryCondition{};
407 
409  RealDefaultBoundaryConditionType m_RealDefaultBoundaryCondition{};
410 };
411 } // end namespace itk
412 
413 #ifndef ITK_MANUAL_INSTANTIATION
414 # include "itkDiscreteGaussianImageFilter.hxx"
415 #endif
416 
417 #endif
itk::DiscreteGaussianImageFilter::SetVariance
void SetVariance(const float *v)
Definition: itkDiscreteGaussianImageFilter.h:190
itk::DiscreteGaussianImageFilter
Blurs an image by separable convolution with discrete gaussian kernels. This filter performs Gaussian...
Definition: itkDiscreteGaussianImageFilter.h:64
itk::DiscreteGaussianImageFilter::SetMaximumError
void SetMaximumError(const typename ArrayType::ValueType v)
Definition: itkDiscreteGaussianImageFilter.h:171
itk::DiscreteGaussianImageFilter::InputInternalPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
Definition: itkDiscreteGaussianImageFilter.h:90
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:71
itk::NumericTraits::ValueType
T ValueType
Definition: itkNumericTraits.h:66
itk::DiscreteGaussianImageFilter::DiscreteGaussianImageFilter
DiscreteGaussianImageFilter()
Definition: itkDiscreteGaussianImageFilter.h:338
itkImage.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::DiscreteGaussianImageFilter::SetVariance
void SetVariance(const typename ArrayType::ValueType v)
Definition: itkDiscreteGaussianImageFilter.h:163
itk::DiscreteGaussianImageFilter::SetMaximumError
void SetMaximumError(const float *v)
Definition: itkDiscreteGaussianImageFilter.h:262
itk::DiscreteGaussianImageFilter::SetSigmaArray
void SetSigmaArray(const ArrayType &sigmas)
Definition: itkDiscreteGaussianImageFilter.h:217
itk::DiscreteGaussianImageFilter::SetSigma
void SetSigma(const ArrayType &sigma)
Definition: itkDiscreteGaussianImageFilter.h:204
itk::DiscreteGaussianImageFilter::OutputPixelValueType
typename NumericTraits< OutputPixelType >::ValueType OutputPixelValueType
Definition: itkDiscreteGaussianImageFilter.h:94
itk::ImageBoundaryCondition
A virtual base object that defines an interface to a class of boundary condition objects for use by n...
Definition: itkImageBoundaryCondition.h:52
itk::DiscreteGaussianImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkDiscreteGaussianImageFilter.h:87
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::DiscreteGaussianImageFilter::RadiusType
typename KernelType::RadiusType RadiusType
Definition: itkDiscreteGaussianImageFilter.h:119
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::DiscreteGaussianImageFilter::GetSigmaArray
ArrayType GetSigmaArray() const
Definition: itkDiscreteGaussianImageFilter.h:230
itk::DiscreteGaussianImageFilter::SetMaximumError
void SetMaximumError(const double *v)
Definition: itkDiscreteGaussianImageFilter.h:250
itk::DiscreteGaussianImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkDiscreteGaussianImageFilter.h:89
itk::GaussianOperator
A NeighborhoodOperator whose coefficients are a one dimensional, discrete Gaussian kernel.
Definition: itkGaussianOperator.h:69
itk::DiscreteGaussianImageFilter::SetSigma
void SetSigma(double sigma)
Definition: itkDiscreteGaussianImageFilter.h:222
itk::DiscreteGaussianImageFilter::ScalarRealType
double ScalarRealType
Definition: itkDiscreteGaussianImageFilter.h:115
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::FixedArray< double, Self::ImageDimension >
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:59
itk::DiscreteGaussianImageFilter::RealOutputPixelType
typename NumericTraits< OutputPixelType >::RealType RealOutputPixelType
Definition: itkDiscreteGaussianImageFilter.h:101
itk::DiscreteGaussianImageFilter::OutputInternalPixelType
typename TOutputImage::InternalPixelType OutputInternalPixelType
Definition: itkDiscreteGaussianImageFilter.h:88
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
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::DiscreteGaussianImageFilter::RealOutputPixelValueType
typename NumericTraits< RealOutputPixelType >::ValueType RealOutputPixelValueType
Definition: itkDiscreteGaussianImageFilter.h:103
itk::ZeroFluxNeumannBoundaryCondition< TInputImage >
itk::DiscreteGaussianImageFilter::GetSigma
double GetSigma() const
Definition: itkDiscreteGaussianImageFilter.h:244
itkGaussianOperator.h
itk::DiscreteGaussianImageFilter::SetVariance
void SetVariance(const double *v)
Definition: itkDiscreteGaussianImageFilter.h:178
itk::DiscreteGaussianImageFilter::InputPixelValueType
typename NumericTraits< InputPixelType >::ValueType InputPixelValueType
Definition: itkDiscreteGaussianImageFilter.h:93
itkZeroFluxNeumannBoundaryCondition.h
itk::FixedArray< double, Self::ImageDimension >::ValueType
double ValueType
Definition: itkFixedArray.h:63
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90