ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkMultiScaleHessianBasedMeasureImageFilter.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 itkMultiScaleHessianBasedMeasureImageFilter_h
19 #define itkMultiScaleHessianBasedMeasureImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
64 template< typename TInputImage,
65  typename THessianImage,
66  typename TOutputImage = TInputImage >
67 class ITK_TEMPLATE_EXPORT MultiScaleHessianBasedMeasureImageFilter:
68  public ImageToImageFilter< TInputImage, TOutputImage >
69 {
70 public:
71  ITK_DISALLOW_COPY_AND_ASSIGN(MultiScaleHessianBasedMeasureImageFilter);
72 
76 
79 
80  using InputImageType = TInputImage;
81  using OutputImageType = TOutputImage;
82  using HessianImageType = THessianImage;
83 
85 
86  using InputPixelType = typename TInputImage::PixelType;
87  using OutputPixelType = typename TOutputImage::PixelType;
89 
91  static constexpr unsigned int ImageDimension = InputImageType ::ImageDimension;
92 
94  using ScalesPixelType = float;
96 
99 
105 
107 
109  itkNewMacro(Self);
110 
114 
116  itkSetMacro(SigmaMinimum, double);
117  itkGetConstMacro(SigmaMinimum, double);
119 
121  itkSetMacro(SigmaMaximum, double);
122  itkGetConstMacro(SigmaMaximum, double);
124 
126  itkSetMacro(NumberOfSigmaSteps, unsigned int);
127  itkGetConstMacro(NumberOfSigmaSteps, unsigned int);
129 
133  itkSetObjectMacro(HessianToMeasureFilter, HessianToMeasureFilterType);
134  itkGetModifiableObjectMacro(HessianToMeasureFilter, HessianToMeasureFilterType);
136 
143  itkSetMacro(NonNegativeHessianBasedMeasure, bool);
144  itkGetConstMacro(NonNegativeHessianBasedMeasure, bool);
145  itkBooleanMacro(NonNegativeHessianBasedMeasure);
147 
148  typedef enum { EquispacedSigmaSteps = 0,
149  LogarithmicSigmaSteps = 1 } SigmaStepMethodType;
150 
153  itkSetMacro(SigmaStepMethod, SigmaStepMethodType);
154  itkGetConstMacro(SigmaStepMethod, SigmaStepMethodType);
156 
158  void SetSigmaStepMethodToEquispaced();
159 
161  void SetSigmaStepMethodToLogarithmic();
162 
165  const HessianImageType * GetHessianOutput() const;
166 
169  const ScalesImageType * GetScalesOutput() const;
170 
173  itkSetMacro(GenerateScalesOutput, bool);
174  itkGetConstMacro(GenerateScalesOutput, bool);
175  itkBooleanMacro(GenerateScalesOutput);
177 
180  itkSetMacro(GenerateHessianOutput, bool);
181  itkGetConstMacro(GenerateHessianOutput, bool);
182  itkBooleanMacro(GenerateHessianOutput);
184 
187 
188 protected:
190  ~MultiScaleHessianBasedMeasureImageFilter() override = default;
191  void PrintSelf(std::ostream & os, Indent indent) const override;
192 
194  void GenerateData() override;
195 
196  void EnlargeOutputRequestedRegion(DataObject *) override;
197 
198  using Superclass::MakeOutput;
199  DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override;
200 
201 private:
202  void UpdateMaximumResponse(double sigma);
203 
204  double ComputeSigmaValue(int scaleLevel);
205 
206  void AllocateUpdateBuffer();
207 
209 
212 
213  unsigned int m_NumberOfSigmaSteps;
215 
217 
219 
221 
224 };
225 } // end namespace itk
226 
227 #ifndef ITK_MANUAL_INSTANTIATION
228 #include "itkMultiScaleHessianBasedMeasureImageFilter.hxx"
229 #endif
230 
231 #endif
Superclass::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
TPixel ValueType
Definition: itkImage.h:98
A filter to enhance structures using Hessian eigensystem-based measures in a multiscale framework...
TOutputImage OutputImageType
Computes the Hessian matrix of an image by convolution with the Second and Cross derivatives of a Gau...
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
Base class for all data objects in ITK.
Templated n-dimensional image class.
Definition: itkImage.h:75