ITK  5.4.0
Insight Toolkit
itkHessianToObjectnessMeasureImageFilter.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 itkHessianToObjectnessMeasureImageFilter_h
19#define itkHessianToObjectnessMeasureImageFilter_h
20
23
24namespace itk
25{
60template <typename TInputImage, typename TOutputImage>
61class ITK_TEMPLATE_EXPORT HessianToObjectnessMeasureImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
62{
63public:
64 ITK_DISALLOW_COPY_AND_MOVE(HessianToObjectnessMeasureImageFilter);
65
71
72 using typename Superclass::InputImageType;
73 using typename Superclass::OutputImageType;
74 using InputPixelType = typename InputImageType::PixelType;
75 using OutputPixelType = typename OutputImageType::PixelType;
77
79 static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;
80
81 using EigenValueType = double;
83
85 itkNewMacro(Self);
86
88 itkOverrideGetNameOfClassMacro(HessianToObjectnessMeasureImageFilter);
89
93 itkSetMacro(Alpha, double);
94 itkGetConstMacro(Alpha, double);
100 itkSetMacro(Beta, double);
101 itkGetConstMacro(Beta, double);
106 itkSetMacro(Gamma, double);
107 itkGetConstMacro(Gamma, double);
112 itkSetMacro(ScaleObjectnessMeasure, bool);
113 itkGetConstMacro(ScaleObjectnessMeasure, bool);
114 itkBooleanMacro(ScaleObjectnessMeasure);
120 itkSetMacro(ObjectDimension, unsigned int);
121 itkGetConstMacro(ObjectDimension, unsigned int);
126 itkSetMacro(BrightObject, bool);
127 itkGetConstMacro(BrightObject, bool);
128 itkBooleanMacro(BrightObject);
131#ifdef ITK_USE_CONCEPT_CHECKING
132 // Begin concept checking
133 itkConceptMacro(DoubleConvertibleToOutputCheck, (Concept::Convertible<double, OutputPixelType>));
134 // End concept checking
135#endif
136
137protected:
140 void
141 PrintSelf(std::ostream & os, Indent indent) const override;
142
143 void
144 VerifyPreconditions() ITKv5_CONST override;
145
146 void
147 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
148
149
150private:
151 // functor used to sort the eigenvalues are to be sorted
152 // |e1|<=|e2|<=...<=|eN|
153 //
154 // Returns ( itk::Math::abs(a) < itk::Math::abs(b) )
156 {
157 bool
159 {
160 return itk::Math::abs(a) < itk::Math::abs(b);
161 }
162 };
163
164 double m_Alpha{ 0.5 };
165 double m_Beta{ 0.5 };
166 double m_Gamma{ 5.0 };
167 unsigned int m_ObjectDimension{ 1 };
168 bool m_BrightObject{ true };
169 bool m_ScaleObjectnessMeasure{ true };
170};
171} // end namespace itk
172
173#ifndef ITK_MANUAL_INSTANTIATION
174# include "itkHessianToObjectnessMeasureImageFilter.hxx"
175#endif
176
177#endif
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:54
A filter to enhance M-dimensional objects in N-dimensional images.
void PrintSelf(std::ostream &os, Indent indent) const override
void VerifyPreconditions() ITKv5_CONST override
Verifies that the process object has been configured correctly, that all required inputs are set,...
~HessianToObjectnessMeasureImageFilter() override=default
Base class for all process objects that output image data.
typename OutputImageType::RegionType OutputImageRegionType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
#define itkConceptMacro(name, concept)
bool abs(bool x)
Definition: itkMath.h:844
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....