ITK  5.2.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  * 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 itkDiscreteGaussianImageFilter_h
19 #define itkDiscreteGaussianImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkImage.h"
24 
25 namespace itk
26 {
62 template <typename TInputImage, typename TOutputImage = TInputImage>
63 class ITK_TEMPLATE_EXPORT DiscreteGaussianImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
64 {
65 public:
66  ITK_DISALLOW_COPY_AND_MOVE(DiscreteGaussianImageFilter);
67 
73 
75  itkNewMacro(Self);
76 
79 
81  using InputImageType = TInputImage;
82  using OutputImageType = TOutputImage;
83 
86  using OutputPixelType = typename TOutputImage::PixelType;
87  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
88  using InputPixelType = typename TInputImage::PixelType;
89  using InputInternalPixelType = typename TInputImage::InternalPixelType;
90 
94 
97  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
98 
103 
110 
114  using ScalarRealType = double;
115 
122  itkSetMacro(Variance, ArrayType);
123  itkGetConstMacro(Variance, const ArrayType);
125 
129  itkSetMacro(MaximumError, ArrayType);
130  itkGetConstMacro(MaximumError, const ArrayType);
132 
135  itkGetConstMacro(MaximumKernelWidth, int);
136  itkSetMacro(MaximumKernelWidth, int);
138 
144  itkGetConstMacro(FilterDimensionality, unsigned int);
145  itkSetMacro(FilterDimensionality, unsigned int);
147 
149  itkSetMacro(InputBoundaryCondition, InputBoundaryConditionPointerType);
150  itkGetConstMacro(InputBoundaryCondition, InputBoundaryConditionPointerType);
151  itkSetMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
152  itkGetConstMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
154 
157  void
159  {
160  m_Variance.Fill(v);
161  this->Modified();
162  }
164 
165  void
167  {
168  m_MaximumError.Fill(v);
169  this->Modified();
170  }
171 
172  void
173  SetVariance(const double * v)
174  {
175  ArrayType dv;
176 
177  for (unsigned int i = 0; i < ImageDimension; i++)
178  {
179  dv[i] = v[i];
180  }
181  this->SetVariance(dv);
182  }
183 
184  void
185  SetVariance(const float * v)
186  {
187  ArrayType dv;
188 
189  for (unsigned int i = 0; i < ImageDimension; i++)
190  {
191  dv[i] = v[i];
192  }
193  this->SetVariance(dv);
194  }
195 
200  void
201  SetSigmaArray(const ArrayType & sigmas)
202  {
203  ArrayType variance;
204  for (unsigned int i = 0; i < ImageDimension; i++)
205  {
206  variance[i] = sigmas[i] * sigmas[i];
207  }
208  this->SetVariance(variance);
209  }
210  void
211  SetSigma(double sigma)
212  {
213  this->SetVariance(sigma * sigma);
214  }
216 
218  ArrayType
220  {
221  ArrayType sigmas;
222  for (unsigned int i = 0; i < ImageDimension; i++)
223  {
224  sigmas[i] = std::sqrt(m_Variance[i]);
225  }
226  return sigmas;
227  }
229 
232  double
233  GetSigma() const
234  {
235  return std::sqrt(m_Variance[0]);
236  }
237 
238  void
239  SetMaximumError(const double * v)
240  {
241  ArrayType dv;
242 
243  for (unsigned int i = 0; i < ImageDimension; i++)
244  {
245  dv[i] = v[i];
246  }
247  this->SetMaximumError(dv);
248  }
249 
250  void
251  SetMaximumError(const float * v)
252  {
253  ArrayType dv;
254 
255  for (unsigned int i = 0; i < ImageDimension; i++)
256  {
257  dv[i] = v[i];
258  }
259  this->SetMaximumError(dv);
260  }
261 
267  itkSetMacro(UseImageSpacing, bool);
268  itkGetConstMacro(UseImageSpacing, bool);
269  itkBooleanMacro(UseImageSpacing);
271 
272 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
273 
277  void
278  SetUseImageSpacingOn()
279  {
280  this->SetUseImageSpacing(true);
281  }
282 
286  void
287  SetUseImageSpacingOff()
288  {
289  this->SetUseImageSpacing(false);
290  }
291 #endif
292 
302  itkLegacyMacro(unsigned int GetInternalNumberOfStreamDivisions() const);
303  itkLegacyMacro(void SetInternalNumberOfStreamDivisions(unsigned int));
304 
311  void
312  GenerateInputRequestedRegion() override;
313 
314 #ifdef ITK_USE_CONCEPT_CHECKING
315  // Begin concept checking
316 
317  itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<OutputPixelValueType>));
318 
319  // End concept checking
320 #endif
321 
322 protected:
324  {
325  m_Variance.Fill(0.0);
326  m_MaximumError.Fill(0.01);
327  m_MaximumKernelWidth = 32;
328  m_UseImageSpacing = true;
329  m_FilterDimensionality = ImageDimension;
330  m_InputBoundaryCondition = &m_InputDefaultBoundaryCondition;
331  m_RealBoundaryCondition = &m_RealDefaultBoundaryCondition;
332  }
333 
334  ~DiscreteGaussianImageFilter() override = default;
335  void
336  PrintSelf(std::ostream & os, Indent indent) const override;
337 
343  void
344  GenerateData() override;
345 
346 private:
350 
355 
359 
362 
365 
369 
372 
375 
378 };
379 } // end namespace itk
380 
381 #ifndef ITK_MANUAL_INSTANTIATION
382 # include "itkDiscreteGaussianImageFilter.hxx"
383 #endif
384 
385 #endif
itk::DiscreteGaussianImageFilter::SetVariance
void SetVariance(const float *v)
Definition: itkDiscreteGaussianImageFilter.h:185
itk::DiscreteGaussianImageFilter
Blurs an image by separable convolution with discrete gaussian kernels. This filter performs Gaussian...
Definition: itkDiscreteGaussianImageFilter.h:63
itk::DiscreteGaussianImageFilter::SetMaximumError
void SetMaximumError(const typename ArrayType::ValueType v)
Definition: itkDiscreteGaussianImageFilter.h:166
itk::DiscreteGaussianImageFilter::InputInternalPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
Definition: itkDiscreteGaussianImageFilter.h:89
itk::NumericTraits::ValueType
T ValueType
Definition: itkNumericTraits.h:65
itk::DiscreteGaussianImageFilter::DiscreteGaussianImageFilter
DiscreteGaussianImageFilter()
Definition: itkDiscreteGaussianImageFilter.h:323
itk::DiscreteGaussianImageFilter::m_InputDefaultBoundaryCondition
InputDefaultBoundaryConditionType m_InputDefaultBoundaryCondition
Definition: itkDiscreteGaussianImageFilter.h:371
itk::DiscreteGaussianImageFilter::m_MaximumKernelWidth
int m_MaximumKernelWidth
Definition: itkDiscreteGaussianImageFilter.h:358
itk::DiscreteGaussianImageFilter::m_RealBoundaryCondition
RealBoundaryConditionPointerType m_RealBoundaryCondition
Definition: itkDiscreteGaussianImageFilter.h:374
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:158
itk::DiscreteGaussianImageFilter::m_FilterDimensionality
unsigned int m_FilterDimensionality
Definition: itkDiscreteGaussianImageFilter.h:361
itk::DiscreteGaussianImageFilter::SetMaximumError
void SetMaximumError(const float *v)
Definition: itkDiscreteGaussianImageFilter.h:251
itk::DiscreteGaussianImageFilter::SetSigmaArray
void SetSigmaArray(const ArrayType &sigmas)
Definition: itkDiscreteGaussianImageFilter.h:201
itk::DiscreteGaussianImageFilter::OutputPixelValueType
typename NumericTraits< OutputPixelType >::ValueType OutputPixelValueType
Definition: itkDiscreteGaussianImageFilter.h:93
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:86
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::DiscreteGaussianImageFilter::GetSigmaArray
ArrayType GetSigmaArray() const
Definition: itkDiscreteGaussianImageFilter.h:219
itk::DiscreteGaussianImageFilter::SetMaximumError
void SetMaximumError(const double *v)
Definition: itkDiscreteGaussianImageFilter.h:239
itk::DiscreteGaussianImageFilter::m_MaximumError
ArrayType m_MaximumError
Definition: itkDiscreteGaussianImageFilter.h:354
itk::DiscreteGaussianImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkDiscreteGaussianImageFilter.h:88
itk::DiscreteGaussianImageFilter::SetSigma
void SetSigma(double sigma)
Definition: itkDiscreteGaussianImageFilter.h:211
itk::DiscreteGaussianImageFilter::ScalarRealType
double ScalarRealType
Definition: itkDiscreteGaussianImageFilter.h:114
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:58
itk::DiscreteGaussianImageFilter::RealOutputPixelType
typename NumericTraits< OutputPixelType >::RealType RealOutputPixelType
Definition: itkDiscreteGaussianImageFilter.h:100
itk::DiscreteGaussianImageFilter::m_InputBoundaryCondition
InputBoundaryConditionPointerType m_InputBoundaryCondition
Definition: itkDiscreteGaussianImageFilter.h:368
itk::DiscreteGaussianImageFilter::OutputInternalPixelType
typename TOutputImage::InternalPixelType OutputInternalPixelType
Definition: itkDiscreteGaussianImageFilter.h:87
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::DiscreteGaussianImageFilter::m_UseImageSpacing
bool m_UseImageSpacing
Definition: itkDiscreteGaussianImageFilter.h:364
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::DiscreteGaussianImageFilter::RealOutputPixelValueType
typename NumericTraits< RealOutputPixelType >::ValueType RealOutputPixelValueType
Definition: itkDiscreteGaussianImageFilter.h:102
itk::DiscreteGaussianImageFilter::m_RealDefaultBoundaryCondition
RealDefaultBoundaryConditionType m_RealDefaultBoundaryCondition
Definition: itkDiscreteGaussianImageFilter.h:377
itk::ZeroFluxNeumannBoundaryCondition< TInputImage >
itk::DiscreteGaussianImageFilter::GetSigma
double GetSigma() const
Definition: itkDiscreteGaussianImageFilter.h:233
itk::DiscreteGaussianImageFilter::SetVariance
void SetVariance(const double *v)
Definition: itkDiscreteGaussianImageFilter.h:173
itk::DiscreteGaussianImageFilter::InputPixelValueType
typename NumericTraits< InputPixelType >::ValueType InputPixelValueType
Definition: itkDiscreteGaussianImageFilter.h:92
itkZeroFluxNeumannBoundaryCondition.h
itk::DiscreteGaussianImageFilter::m_Variance
ArrayType m_Variance
Definition: itkDiscreteGaussianImageFilter.h:349
itk::FixedArray< double, Self::ImageDimension >::ValueType
double ValueType
Definition: itkFixedArray.h:62
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90