ITK  5.4.0
Insight Toolkit
itkDiscreteGaussianDerivativeImageFilter.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 itkDiscreteGaussianDerivativeImageFilter_h
19 #define itkDiscreteGaussianDerivativeImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkImage.h"
23 
24 namespace itk
25 {
58 template <typename TInputImage, typename TOutputImage>
59 class ITK_TEMPLATE_EXPORT DiscreteGaussianDerivativeImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
60 {
61 public:
62  ITK_DISALLOW_COPY_AND_MOVE(DiscreteGaussianDerivativeImageFilter);
63 
69 
71  itkNewMacro(Self);
72 
74  itkOverrideGetNameOfClassMacro(DiscreteGaussianDerivativeImageFilter);
75 
77  using InputImageType = TInputImage;
78  using OutputImageType = TOutputImage;
79 
82  using OutputPixelType = typename TOutputImage::PixelType;
83  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
84  using InputPixelType = typename TInputImage::PixelType;
85  using InputInternalPixelType = typename TInputImage::InternalPixelType;
86 
89  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
90 
93 
96 
100  itkSetMacro(Order, OrderArrayType);
101  itkGetConstMacro(Order, const OrderArrayType);
110  itkSetMacro(Variance, ArrayType);
111  itkGetConstMacro(Variance, const ArrayType);
117  itkSetMacro(MaximumError, ArrayType);
118  itkGetConstMacro(MaximumError, const ArrayType);
123  itkGetConstMacro(MaximumKernelWidth, int);
124  itkSetMacro(MaximumKernelWidth, int);
136  itkSetMacro(InternalNumberOfStreamDivisions, unsigned int);
137  itkGetConstMacro(InternalNumberOfStreamDivisions, unsigned int);
138 
143  void
145  {
146  OrderArrayType a;
147 
148  a.Fill(v);
149  this->SetOrder(a);
150  }
151 
152  void
154  {
155  ArrayType a;
156 
157  a.Fill(v);
158  this->SetVariance(a);
159  }
160 
161  void
163  {
164  ArrayType a;
165 
166  a.Fill(v);
167  this->SetMaximumError(a);
168  }
169 
174  itkSetMacro(UseImageSpacing, bool);
175  itkGetConstMacro(UseImageSpacing, bool);
176  itkBooleanMacro(UseImageSpacing);
182  itkSetMacro(NormalizeAcrossScale, bool);
183  itkGetConstMacro(NormalizeAcrossScale, bool);
184  itkBooleanMacro(NormalizeAcrossScale);
187 #ifdef ITK_USE_CONCEPT_CHECKING
188  // Begin concept checking
189  itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<OutputPixelType>));
190  // End concept checking
191 #endif
192 
193 protected:
195  {
196  m_Order.Fill(1);
197  m_Variance.Fill(0.0);
198  m_MaximumError.Fill(0.01);
199  m_MaximumKernelWidth = 32;
200  m_UseImageSpacing = true;
201  m_NormalizeAcrossScale = false;
202  m_InternalNumberOfStreamDivisions = ImageDimension * ImageDimension;
203  }
204 
205  ~DiscreteGaussianDerivativeImageFilter() override = default;
206  void
207  PrintSelf(std::ostream & os, Indent indent) const override;
208 
215  void
216  GenerateInputRequestedRegion() override;
217 
223  void
224  GenerateData() override;
225 
226 private:
228  OrderArrayType m_Order{};
229 
232  ArrayType m_Variance{};
233 
237  ArrayType m_MaximumError{};
238 
241  int m_MaximumKernelWidth{};
242 
244  bool m_UseImageSpacing{};
245 
247  bool m_NormalizeAcrossScale{};
248 
251  unsigned int m_InternalNumberOfStreamDivisions{};
252 };
253 } // end namespace itk
254 
255 #ifndef ITK_MANUAL_INSTANTIATION
256 # include "itkDiscreteGaussianDerivativeImageFilter.hxx"
257 #endif
258 
259 #endif
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:714
itk::DiscreteGaussianDerivativeImageFilter::OutputInternalPixelType
typename TOutputImage::InternalPixelType OutputInternalPixelType
Definition: itkDiscreteGaussianDerivativeImageFilter.h:83
itk::DiscreteGaussianDerivativeImageFilter::SetMaximumError
void SetMaximumError(const typename ArrayType::ValueType v)
Definition: itkDiscreteGaussianDerivativeImageFilter.h:162
itkImage.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::DiscreteGaussianDerivativeImageFilter
Calculates image derivatives using discrete derivative gaussian kernels. This filter calculates Gauss...
Definition: itkDiscreteGaussianDerivativeImageFilter.h:59
itk::DiscreteGaussianDerivativeImageFilter::SetVariance
void SetVariance(const typename ArrayType::ValueType v)
Definition: itkDiscreteGaussianDerivativeImageFilter.h:153
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::DiscreteGaussianDerivativeImageFilter::InputInternalPixelType
typename TInputImage::InternalPixelType InputInternalPixelType
Definition: itkDiscreteGaussianDerivativeImageFilter.h:85
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::DiscreteGaussianDerivativeImageFilter::OutputPixelType
typename TOutputImage::PixelType OutputPixelType
Definition: itkDiscreteGaussianDerivativeImageFilter.h:82
itkImageToImageFilter.h
itk::FixedArray< double, Self::ImageDimension >
itk::DiscreteGaussianDerivativeImageFilter::DiscreteGaussianDerivativeImageFilter
DiscreteGaussianDerivativeImageFilter()
Definition: itkDiscreteGaussianDerivativeImageFilter.h:194
itk::DiscreteGaussianDerivativeImageFilter::SetOrder
void SetOrder(const typename OrderArrayType::ValueType v)
Definition: itkDiscreteGaussianDerivativeImageFilter.h:144
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::DiscreteGaussianDerivativeImageFilter::InputPixelType
typename TInputImage::PixelType InputPixelType
Definition: itkDiscreteGaussianDerivativeImageFilter.h:84
itk::FixedArray< unsigned int, Self::ImageDimension >::ValueType
unsigned int ValueType
Definition: itkFixedArray.h:63
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90
itk::FixedArray::Fill
void Fill(const ValueType &)