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 __itkHistogramMatchingImageFilter_h 00019 #define __itkHistogramMatchingImageFilter_h 00020 00021 #include "itkImageToImageFilter.h" 00022 #include "itkHistogram.h" 00023 #include "vnl/vnl_matrix.h" 00024 00025 namespace itk 00026 { 00060 /* THistogramMeasurement -- The precision level for which to do 00061 HistogramMeasurmenets */ 00062 template< class TInputImage, class TOutputImage, class THistogramMeasurement = typename TInputImage::PixelType > 00063 class ITK_EXPORT HistogramMatchingImageFilter: 00064 public ImageToImageFilter< TInputImage, TOutputImage > 00065 { 00066 public: 00068 typedef HistogramMatchingImageFilter Self; 00069 typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass; 00070 typedef SmartPointer< Self > Pointer; 00071 typedef SmartPointer< const Self > ConstPointer; 00072 00074 itkNewMacro(Self); 00075 00077 itkTypeMacro(HistogramMatchingImageFilter, ImageToImageFilter); 00078 00080 itkStaticConstMacro(ImageDimension, unsigned int, 00081 TInputImage::ImageDimension); 00082 itkStaticConstMacro(OutputImageDimension, unsigned int, 00083 TOutputImage::ImageDimension); 00085 00087 typedef typename TOutputImage::RegionType OutputImageRegionType; 00088 00090 typedef typename Superclass::InputImageType InputImageType; 00091 typedef typename Superclass::InputImagePointer InputImagePointer; 00092 typedef typename Superclass::InputImageConstPointer InputImageConstPointer; 00093 typedef typename Superclass::OutputImageType OutputImageType; 00094 typedef typename Superclass::OutputImagePointer OutputImagePointer; 00095 00097 typedef typename InputImageType::PixelType InputPixelType; 00098 typedef typename OutputImageType::PixelType OutputPixelType; 00099 00101 typedef Statistics::Histogram< THistogramMeasurement > HistogramType; 00102 typedef typename HistogramType::Pointer HistogramPointer; 00103 00105 void SetSourceImage(const InputImageType *source) 00106 { this->SetInput(source); } 00107 const InputImageType * GetSourceImage(void) 00108 { return this->GetInput(); } 00110 00112 void SetReferenceImage(const InputImageType *reference); 00113 00114 const InputImageType * GetReferenceImage(void); 00115 00117 itkSetMacro(NumberOfHistogramLevels, SizeValueType); 00118 itkGetConstMacro(NumberOfHistogramLevels, SizeValueType); 00120 00122 itkSetMacro(NumberOfMatchPoints, SizeValueType); 00123 itkGetConstMacro(NumberOfMatchPoints, SizeValueType); 00125 00131 itkSetMacro(ThresholdAtMeanIntensity, bool); 00132 itkGetConstMacro(ThresholdAtMeanIntensity, bool); 00133 itkBooleanMacro(ThresholdAtMeanIntensity); 00135 00137 virtual void GenerateInputRequestedRegion(); 00138 00142 itkGetObjectMacro(SourceHistogram, HistogramType); 00143 itkGetObjectMacro(ReferenceHistogram, HistogramType); 00144 itkGetObjectMacro(OutputHistogram, HistogramType); 00146 00147 #ifdef ITK_USE_CONCEPT_CHECKING 00148 00149 itkConceptMacro( IntConvertibleToInputCheck, 00150 ( Concept::Convertible< int, InputPixelType > ) ); 00151 itkConceptMacro( SameDimensionCheck, 00152 ( Concept::SameDimension< ImageDimension, OutputImageDimension > ) ); 00153 itkConceptMacro( DoubleConvertibleToInputCheck, 00154 ( Concept::Convertible< double, InputPixelType > ) ); 00155 itkConceptMacro( DoubleConvertibleToOutputCheck, 00156 ( Concept::Convertible< double, OutputPixelType > ) ); 00157 itkConceptMacro( InputConvertibleToDoubleCheck, 00158 ( Concept::Convertible< InputPixelType, double > ) ); 00159 itkConceptMacro( OutputConvertibleToDoubleCheck, 00160 ( Concept::Convertible< OutputPixelType, double > ) ); 00161 itkConceptMacro( SameTypeCheck, 00162 ( Concept::SameType< InputPixelType, OutputPixelType > ) ); 00163 00165 #endif 00166 protected: 00167 HistogramMatchingImageFilter(); 00168 ~HistogramMatchingImageFilter() {} 00169 void PrintSelf(std::ostream & os, Indent indent) const; 00171 00172 void BeforeThreadedGenerateData(); 00173 00174 void AfterThreadedGenerateData(); 00175 00176 void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread, 00177 ThreadIdType threadId); 00178 00184 virtual void VerifyInputInformation() {} 00185 00187 void ComputeMinMaxMean(const InputImageType *image, 00188 THistogramMeasurement & minValue, 00189 THistogramMeasurement & maxValue, 00190 THistogramMeasurement & meanValue); 00191 00193 void ConstructHistogram(const InputImageType *image, 00194 HistogramType *histogram, const THistogramMeasurement minValue, 00195 const THistogramMeasurement maxValue); 00196 00197 private: 00198 HistogramMatchingImageFilter(const Self &); //purposely not implemented 00199 void operator=(const Self &); //purposely not implemented 00200 00201 SizeValueType m_NumberOfHistogramLevels; 00202 SizeValueType m_NumberOfMatchPoints; 00203 bool m_ThresholdAtMeanIntensity; 00204 00205 InputPixelType m_SourceIntensityThreshold; 00206 InputPixelType m_ReferenceIntensityThreshold; 00207 OutputPixelType m_OutputIntensityThreshold; 00208 00209 THistogramMeasurement m_SourceMinValue; 00210 THistogramMeasurement m_SourceMaxValue; 00211 THistogramMeasurement m_SourceMeanValue; 00212 THistogramMeasurement m_ReferenceMinValue; 00213 THistogramMeasurement m_ReferenceMaxValue; 00214 THistogramMeasurement m_ReferenceMeanValue; 00215 THistogramMeasurement m_OutputMinValue; 00216 THistogramMeasurement m_OutputMaxValue; 00217 THistogramMeasurement m_OutputMeanValue; 00218 00219 HistogramPointer m_SourceHistogram; 00220 HistogramPointer m_ReferenceHistogram; 00221 HistogramPointer m_OutputHistogram; 00222 00223 typedef vnl_matrix< double > TableType; 00224 TableType m_QuantileTable; 00225 00226 typedef vnl_vector< double > GradientArrayType; 00227 GradientArrayType m_Gradients; 00228 double m_LowerGradient; 00229 double m_UpperGradient; 00230 }; 00231 } // end namespace itk 00232 00233 #ifndef ITK_MANUAL_INSTANTIATION 00234 #include "itkHistogramMatchingImageFilter.hxx" 00235 #endif 00236 00237 #endif 00238