ITK  4.6.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 >
68  public
69  ImageToImageFilter< TInputImage, TOutputImage >
70 {
71 public:
75 
78 
79  typedef TInputImage InputImageType;
80  typedef TOutputImage OutputImageType;
81  typedef THessianImage HessianImageType;
82 
84 
85  typedef typename TInputImage::PixelType InputPixelType;
86  typedef typename TOutputImage::PixelType OutputPixelType;
87  typedef typename TOutputImage::RegionType OutputRegionType;
88 
90  itkStaticConstMacro(ImageDimension, unsigned int, InputImageType ::ImageDimension);
91 
93  typedef float ScalesPixelType;
95 
98 
104 
106 
108  itkNewMacro(Self);
109 
113 
115  itkSetMacro(SigmaMinimum, double);
116  itkGetConstMacro(SigmaMinimum, double);
118 
120  itkSetMacro(SigmaMaximum, double);
121  itkGetConstMacro(SigmaMaximum, double);
123 
125  itkSetMacro(NumberOfSigmaSteps, unsigned int);
126  itkGetConstMacro(NumberOfSigmaSteps, unsigned int);
128 
132  itkSetObjectMacro(HessianToMeasureFilter, HessianToMeasureFilterType);
133  itkGetModifiableObjectMacro(HessianToMeasureFilter, HessianToMeasureFilterType);
135 
142  itkSetMacro(NonNegativeHessianBasedMeasure, bool);
143  itkGetConstMacro(NonNegativeHessianBasedMeasure, bool);
144  itkBooleanMacro(NonNegativeHessianBasedMeasure);
146 
147  typedef enum { EquispacedSigmaSteps = 0,
149 
152  itkSetMacro(SigmaStepMethod, SigmaStepMethodType);
153  itkGetConstMacro(SigmaStepMethod, SigmaStepMethodType);
155 
158 
161 
164  const HessianImageType * GetHessianOutput() const;
165 
168  const ScalesImageType * GetScalesOutput() const;
169 
171 
174  itkSetMacro(GenerateScalesOutput, bool);
175  itkGetConstMacro(GenerateScalesOutput, bool);
176  itkBooleanMacro(GenerateScalesOutput);
178 
181  itkSetMacro(GenerateHessianOutput, bool);
182  itkGetConstMacro(GenerateHessianOutput, bool);
183  itkBooleanMacro(GenerateHessianOutput);
185 
190 
191 protected:
194  void PrintSelf(std::ostream & os, Indent indent) const;
195 
197  void GenerateData(void);
198 
199 private:
200  void UpdateMaximumResponse(double sigma);
201 
202  double ComputeSigmaValue(int scaleLevel);
203 
204  void AllocateUpdateBuffer();
205 
206  //purposely not implemented
208  void operator=(const Self &); //purposely not implemented
209 
211 
214 
215  unsigned int m_NumberOfSigmaSteps;
217 
219 
221 
223 
226 };
227 } // end namespace itk
228 
229 #ifndef ITK_MANUAL_INSTANTIATION
230 #include "itkMultiScaleHessianBasedMeasureImageFilter.hxx"
231 #endif
232 
233 #endif
Superclass::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
Image< double, itkGetStaticConstMacro(ImageDimension) > UpdateBufferType
Image< ScalesPixelType, itkGetStaticConstMacro(ImageDimension) > ScalesImageType
virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx)
virtual ProcessObject::DataObjectPointer MakeOutput(ProcessObject::DataObjectPointerArraySizeType idx) ITK_OVERRIDE
const HessianImageType * GetHessianOutput() const
Base class for all process objects that output image data.
TPixel ValueType
Definition: itkImage.h:96
HessianRecursiveGaussianImageFilter< InputImageType, HessianImageType > HessianFilterType
A filter to enhance structures using Hessian eigensystem-based measures in a multiscale framework...
ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
const ScalesImageType * GetScalesOutput() const
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
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
ImageToImageFilter< HessianImageType, OutputImageType > HessianToMeasureFilterType
Base class for all data objects in ITK.
void PrintSelf(std::ostream &os, Indent indent) const
Templated n-dimensional image class.
Definition: itkImage.h:75
ImageToImageFilter< TInputImage, TOutputImage > Superclass