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 __itkScalarImageToCooccurrenceMatrixFilter_h 00019 #define __itkScalarImageToCooccurrenceMatrixFilter_h 00020 00021 #include "itkImage.h" 00022 #include "itkHistogram.h" 00023 #include "itkVectorContainer.h" 00024 #include "itkNumericTraits.h" 00025 00026 namespace itk 00027 { 00028 namespace Statistics 00029 { 00092 template< class TImageType, 00093 class THistogramFrequencyContainer = DenseFrequencyContainer2 > 00094 class ITK_EXPORT ScalarImageToCooccurrenceMatrixFilter:public ProcessObject 00095 { 00096 public: 00098 typedef ScalarImageToCooccurrenceMatrixFilter Self; 00099 typedef ProcessObject Superclass; 00100 typedef SmartPointer< Self > Pointer; 00101 typedef SmartPointer< const Self > ConstPointer; 00102 00104 itkTypeMacro(ScalarImageToCooccurrenceMatrixFilter, ProcessObject); 00105 00107 itkNewMacro(Self); 00108 00109 typedef TImageType ImageType; 00110 typedef typename ImageType::Pointer ImagePointer; 00111 typedef typename ImageType::ConstPointer ImageConstPointer; 00112 typedef typename ImageType::PixelType PixelType; 00113 typedef typename ImageType::RegionType RegionType; 00114 typedef typename ImageType::SizeType RadiusType; 00115 typedef typename ImageType::OffsetType OffsetType; 00116 typedef VectorContainer< unsigned char, OffsetType > OffsetVector; 00117 typedef typename OffsetVector::Pointer OffsetVectorPointer; 00118 typedef typename OffsetVector::ConstPointer OffsetVectorConstPointer; 00119 00120 typedef typename NumericTraits< PixelType >::RealType MeasurementType; 00121 00122 typedef Histogram< MeasurementType, THistogramFrequencyContainer > HistogramType; 00123 typedef typename HistogramType::Pointer HistogramPointer; 00124 typedef typename HistogramType::ConstPointer HistogramConstPointer; 00125 typedef typename HistogramType::MeasurementVectorType MeasurementVectorType; 00126 00127 itkStaticConstMacro(DefaultBinsPerAxis, unsigned int, 256); 00128 00131 itkSetConstObjectMacro(Offsets, OffsetVector); 00132 itkGetConstObjectMacro(Offsets, OffsetVector); 00133 void SetOffset(const OffsetType offset); 00135 00137 itkSetMacro(NumberOfBinsPerAxis, unsigned int); 00138 itkGetConstMacro(NumberOfBinsPerAxis, unsigned int); 00140 00143 void SetPixelValueMinMax(PixelType min, PixelType max); 00144 00145 itkGetConstMacro(Min, PixelType); 00146 itkGetConstMacro(Max, PixelType); 00147 00150 itkSetMacro(Normalize, bool); 00151 itkGetConstMacro(Normalize, bool); 00152 itkBooleanMacro(Normalize); 00154 00156 using Superclass::SetInput; 00157 void SetInput(const ImageType *image); 00158 00159 const ImageType * GetInput() const; 00160 00162 void SetMaskImage(const ImageType *image); 00163 00164 const ImageType * GetMaskImage() const; 00165 00167 const HistogramType * GetOutput() const; 00168 00171 itkSetMacro(InsidePixelValue, PixelType); 00172 itkGetConstMacro(InsidePixelValue, PixelType); 00173 protected: 00174 ScalarImageToCooccurrenceMatrixFilter(); 00175 virtual ~ScalarImageToCooccurrenceMatrixFilter() {} 00176 void PrintSelf(std::ostream & os, Indent indent) const; 00178 00179 virtual void FillHistogram(RadiusType radius, RegionType region); 00180 00181 virtual void FillHistogramWithMask(RadiusType radius, RegionType region, const ImageType *maskImage); 00182 00184 typedef DataObject::Pointer DataObjectPointer; 00185 00186 typedef ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType; 00187 using Superclass::MakeOutput; 00188 virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx); 00189 00191 virtual void GenerateData(); 00192 00193 private: 00194 ScalarImageToCooccurrenceMatrixFilter(const Self &); //purposely not 00195 // implemented 00196 void operator=(const Self &); //purposely not 00197 00198 // implemented 00199 00200 void NormalizeHistogram(void); 00201 00202 OffsetVectorConstPointer m_Offsets; 00203 PixelType m_Min; 00204 PixelType m_Max; 00205 00206 unsigned int m_NumberOfBinsPerAxis; 00207 MeasurementVectorType m_LowerBound; 00208 MeasurementVectorType m_UpperBound; 00209 bool m_Normalize; 00210 00211 PixelType m_InsidePixelValue; 00212 }; 00213 } // end of namespace Statistics 00214 } // end of namespace itk 00215 00216 #ifndef ITK_MANUAL_INSTANTIATION 00217 #include "itkScalarImageToCooccurrenceMatrixFilter.hxx" 00218 #endif 00219 00220 #endif 00221