ITK  6.0.0
Insight Toolkit
itkHessianRecursiveGaussianImageFilter.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 itkHessianRecursiveGaussianImageFilter_h
19#define itkHessianRecursiveGaussianImageFilter_h
20
23#include "itkImage.h"
25#include "itkPixelTraits.h"
26
27namespace itk
28{
42template <typename TInputImage,
43 typename TOutputImage =
44 Image<SymmetricSecondRankTensor<typename NumericTraits<typename TInputImage::PixelType>::RealType,
45 TInputImage::ImageDimension>,
46 TInputImage::ImageDimension>>
47class ITK_TEMPLATE_EXPORT HessianRecursiveGaussianImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
48{
49public:
50 ITK_DISALLOW_COPY_AND_MOVE(HessianRecursiveGaussianImageFilter);
51
57
59 using InputImageType = TInputImage;
60 using PixelType = typename TInputImage::PixelType;
62
64 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
65
67 static constexpr unsigned int NumberOfSmoothingFilters =
68 (TInputImage::ImageDimension > 2) ? (TInputImage::ImageDimension - 2) : (0);
69
74 using InternalRealType = float;
76
82
84
87
90
92
95 using GaussianFiltersArray = std::vector<GaussianFilterPointer>;
96
100
103
105 using OutputImageType = TOutputImage;
106 using OutputPixelType = typename OutputImageType::PixelType;
108
110 itkOverrideGetNameOfClassMacro(HessianRecursiveGaussianImageFilter);
111
113 itkNewMacro(Self);
114
116 void
119 GetSigma() const;
125 void
126 SetNormalizeAcrossScale(bool normalize);
127 itkGetConstMacro(NormalizeAcrossScale, bool);
128 itkBooleanMacro(NormalizeAcrossScale);
136 void
138
139#ifdef ITK_USE_CONCEPT_CHECKING
140 // Begin concept checking
141 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelType>));
142 itkConceptMacro(OutputHasPixelTraitsCheck, (Concept::HasPixelTraits<OutputPixelType>));
143 // End concept checking
144#endif
145
146protected:
149 void
150 PrintSelf(std::ostream & os, Indent indent) const override;
151
153 void
154 GenerateData() override;
155
156 // Override since the filter produces the entire dataset
157 void
159
160private:
161 GaussianFiltersArray m_SmoothingFilters{};
162 DerivativeFilterAPointer m_DerivativeFilterA{};
163 DerivativeFilterBPointer m_DerivativeFilterB{};
165
167 bool m_NormalizeAcrossScale{};
168};
169} // end namespace itk
170
171#ifndef ITK_MANUAL_INSTANTIATION
172# include "itkHessianRecursiveGaussianImageFilter.hxx"
173#endif
174
175#endif
Base class for all data objects in ITK.
Computes the Hessian matrix of an image by convolution with the Second and Cross derivatives of a Gau...
typename NumericTraits< PixelType >::RealType RealType
void PrintSelf(std::ostream &os, Indent indent) const override
typename OutputImageAdaptorType::Pointer OutputImageAdaptorPointer
~HessianRecursiveGaussianImageFilter() override=default
void EnlargeOutputRequestedRegion(DataObject *output) override
void SetNormalizeAcrossScale(bool normalize)
typename DerivativeFilterAType::Pointer DerivativeFilterAPointer
typename DerivativeFilterBType::Pointer DerivativeFilterBPointer
typename PixelTraits< OutputPixelType >::ValueType OutputComponentType
Base class for filters that take an image as input and produce an image as output.
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Presents an image as being composed of the N-th element of its pixels.
Define additional traits for native types such as int or float.
typename TPixelType::ValueType ValueType
Base class for computing IIR convolution with an approximation of a Gaussian kernel.
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....