ITK
4.1.0
Insight Segmentation and Registration Toolkit
|
00001 /*========================================================================= 00002 * 00003 * Copyright Insight Software Consortium 00004 * 00005 * Licensed under the Apache License, Version 2.0 (the "License"); 00006 * you may not use this file except in compliance with the License. 00007 * You may obtain a copy of the License at 00008 * 00009 * http://www.apache.org/licenses/LICENSE-2.0.txt 00010 * 00011 * Unless required by applicable law or agreed to in writing, software 00012 * distributed under the License is distributed on an "AS IS" BASIS, 00013 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. 00014 * See the License for the specific language governing permissions and 00015 * limitations under the License. 00016 * 00017 *=========================================================================*/ 00018 #ifndef __itkHessianToObjectnessMeasureImageFilter_h 00019 #define __itkHessianToObjectnessMeasureImageFilter_h 00020 00021 #include "itkSymmetricSecondRankTensor.h" 00022 #include "itkImageToImageFilter.h" 00023 00024 namespace itk 00025 { 00061 template< typename TInputImage, typename TOutputImage > 00062 class ITK_EXPORT HessianToObjectnessMeasureImageFilter:public 00063 ImageToImageFilter< TInputImage, TOutputImage > 00064 { 00065 public: 00067 typedef HessianToObjectnessMeasureImageFilter Self; 00068 00069 typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass; 00070 00071 typedef SmartPointer< Self > Pointer; 00072 typedef SmartPointer< const Self > ConstPointer; 00073 00074 typedef typename Superclass::InputImageType InputImageType; 00075 typedef typename Superclass::OutputImageType OutputImageType; 00076 typedef typename InputImageType::PixelType InputPixelType; 00077 typedef typename OutputImageType::PixelType OutputPixelType; 00078 typedef typename OutputImageType::RegionType OutputImageRegionType; 00079 00081 itkStaticConstMacro(ImageDimension, unsigned int, ::itk::GetImageDimension< InputImageType >::ImageDimension); 00082 00083 typedef double EigenValueType; 00084 typedef itk::FixedArray< EigenValueType, itkGetStaticConstMacro(ImageDimension) > EigenValueArrayType; 00085 00087 itkNewMacro(Self); 00088 00090 itkTypeMacro(HessianToObjectnessMeasureImageFilter, ImageToImageFilter); 00091 00095 itkSetMacro(Alpha, double); 00096 itkGetConstMacro(Alpha, double); 00098 00102 itkSetMacro(Beta, double); 00103 itkGetConstMacro(Beta, double); 00105 00108 itkSetMacro(Gamma, double); 00109 itkGetConstMacro(Gamma, double); 00111 00114 itkSetMacro(ScaleObjectnessMeasure, bool); 00115 itkGetConstMacro(ScaleObjectnessMeasure, bool); 00116 itkBooleanMacro(ScaleObjectnessMeasure); 00118 00122 itkSetMacro(ObjectDimension, unsigned int); 00123 itkGetConstMacro(ObjectDimension, unsigned int); 00125 00128 itkSetMacro(BrightObject, bool); 00129 itkGetConstMacro(BrightObject, bool); 00130 itkBooleanMacro(BrightObject); 00132 00133 #ifdef ITK_USE_CONCEPT_CHECKING 00134 00135 itkConceptMacro( DoubleConvertibleToOutputCheck, ( Concept::Convertible< double, OutputPixelType > ) ); 00136 00138 #endif 00139 protected: 00140 HessianToObjectnessMeasureImageFilter(); 00141 ~HessianToObjectnessMeasureImageFilter() {} 00142 void PrintSelf(std::ostream & os, Indent indent) const; 00144 00145 void BeforeThreadedGenerateData(void); 00146 00147 void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, ThreadIdType threadId); 00148 00149 private: 00150 HessianToObjectnessMeasureImageFilter(const Self &); //purposely not 00151 // implemented 00152 void operator=(const Self &); //purposely not 00153 // implemented 00154 00155 // functor used to sort the eigenvalues are to be sorted 00156 // |e1|<=|e2|<=...<=|eN| 00161 struct AbsLessEqualCompare { 00162 bool operator()(EigenValueType a, EigenValueType b) 00163 { 00164 return vnl_math_abs(a) <= vnl_math_abs(b); 00165 } 00166 }; 00167 00168 double m_Alpha; 00169 double m_Beta; 00170 double m_Gamma; 00171 unsigned int m_ObjectDimension; 00172 bool m_BrightObject; 00173 bool m_ScaleObjectnessMeasure; 00174 }; 00175 } // end namespace itk 00176 00177 #ifndef ITK_MANUAL_INSTANTIATION 00178 #include "itkHessianToObjectnessMeasureImageFilter.hxx" 00179 #endif 00180 00181 #endif 00182