ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkMRFImageFilter.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 itkMRFImageFilter_h
19 #define itkMRFImageFilter_h
20 
21 #include "vnl/vnl_vector.h"
22 #include "vnl/vnl_matrix.h"
23 
24 #include "itkImageClassifierBase.h"
25 
26 #include "itkImageToImageFilter.h"
27 
30 #include "itkSize.h"
31 
32 namespace itk
33 {
124 template< typename TInputImage, typename TClassifiedImage >
126  public ImageToImageFilter< TInputImage, TClassifiedImage >
127 {
128 public:
135 
137  itkNewMacro(Self);
138 
140  itkTypeMacro(MRFImageFilter, Object);
141 
143  typedef TInputImage InputImageType;
144  typedef typename TInputImage::Pointer InputImagePointer;
145  typedef typename TInputImage::ConstPointer InputImageConstPointer;
146 
148  typedef typename TInputImage::PixelType InputImagePixelType;
149 
151  typedef typename TInputImage::RegionType InputImageRegionType;
152 
156 
158  itkStaticConstMacro(InputImageDimension, unsigned int,
159  TInputImage::ImageDimension);
160 
162  typedef typename TClassifiedImage::Pointer TrainingImagePointer;
163 
165  typedef typename TClassifiedImage::PixelType TrainingImagePixelType;
166 
169  typedef typename TClassifiedImage::Pointer LabelledImagePointer;
170 
173  typedef typename TClassifiedImage::PixelType LabelledImagePixelType;
174 
177  typedef typename TClassifiedImage::RegionType LabelledImageRegionType;
178 
180  typedef typename TClassifiedImage::IndexType LabelledImageIndexType;
182 
184  typedef typename TClassifiedImage::OffsetType LabelledImageOffsetType;
185 
189 
191  itkStaticConstMacro(ClassifiedImageDimension, unsigned int,
192  TClassifiedImage::ImageDimension);
193 
196 
198  typedef typename TInputImage::SizeType SizeType;
199 
201  typedef typename TInputImage::SizeType NeighborhoodRadiusType;
202 
206 
209 
212 
215 
216  typedef typename InputImageFaceListType::iterator
218 
222 
225 
228 
231 
232  typedef typename LabelledImageFaceListType::iterator
234 
236  void SetClassifier(typename ClassifierType::Pointer ptrToClassifier);
237 
239  itkSetMacro(NumberOfClasses, unsigned int);
240  itkGetConstMacro(NumberOfClasses, unsigned int);
242 
245  itkSetMacro(MaximumNumberOfIterations, unsigned int);
246  itkGetConstMacro(MaximumNumberOfIterations, unsigned int);
248 
251  itkSetMacro(ErrorTolerance, double);
252  itkGetConstMacro(ErrorTolerance, double);
254 
257  itkSetMacro(SmoothingFactor, double);
258  itkGetConstMacro(SmoothingFactor, double);
260 
263 
268 
269  void SetNeighborhoodRadius(const SizeValueType *radiusArray);
270 
273  {
274  NeighborhoodRadiusType radius;
275 
276  for ( int i = 0; i < InputImageDimension; ++i )
277  {
278  radius[i] = m_InputImageNeighborhoodRadius[i];
279  }
280  return radius;
281  }
282 
288  virtual void SetMRFNeighborhoodWeight(std::vector< double > BetaMatrix);
289 
290  virtual std::vector< double > GetMRFNeighborhoodWeight()
291  {
293  }
294 
295 //Enum to get the stopping condition of the MRF filter
296  typedef enum {
300 
303  itkGetConstReferenceMacro(StopCondition, StopConditionType);
304 
305  /* Get macro for number of iterations */
306  itkGetConstReferenceMacro(NumberOfIterations, unsigned int);
307 
308 #ifdef ITK_USE_CONCEPT_CHECKING
309  // Begin concept checking
310  itkConceptMacro( UnsignedIntConvertibleToClassifiedCheck,
312  itkConceptMacro( ClassifiedConvertibleToUnsignedIntCheck,
314  itkConceptMacro( ClassifiedConvertibleToIntCheck,
316  itkConceptMacro( IntConvertibleToClassifiedCheck,
318  itkConceptMacro( SameDimensionCheck,
320  // End concept checking
321 #endif
322 
323 protected:
324  MRFImageFilter();
325  ~MRFImageFilter();
326  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
327 
329  void Allocate();
330 
335  virtual void ApplyMRFImageFilter();
336 
338  virtual void MinimizeFunctional();
339 
345 
349  //Function implementing the neighborhood operation
350 
351  virtual void DoNeighborhoodOperation(const InputImageNeighborhoodIterator & imageIter,
352  LabelledImageNeighborhoodIterator & labelledIter,
353  LabelStatusImageNeighborhoodIterator & labelStatusIter);
354 
355  virtual void GenerateData() ITK_OVERRIDE;
356 
357  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
358 
359  virtual void EnlargeOutputRequestedRegion(DataObject *) ITK_OVERRIDE;
360 
361  virtual void GenerateOutputInformation() ITK_OVERRIDE;
362 
363 private:
364  MRFImageFilter(const Self &); //purposely not implemented
365  void operator=(const Self &); //purposely not implemented
366 
367  typedef typename TInputImage::SizeType InputImageSizeType;
368 
369  typedef typename LabelStatusImageNeighborhoodIterator::RadiusType
371 
372  typedef NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< LabelStatusImageType >
374 
375  typedef typename LabelStatusImageFacesCalculator::FaceListType
377 
378  typedef typename LabelStatusImageFaceListType::iterator
380 
384 
385  unsigned int m_NumberOfClasses;
387  unsigned int m_KernelSize;
388 
395  double * m_ClassProbability; //Class liklihood
396  unsigned int m_NumberOfIterations;
398 
400 
401  std::vector< double > m_MRFNeighborhoodWeight;
402  std::vector< double > m_NeighborInfluence;
403  std::vector< double > m_MahalanobisDistance;
404  std::vector< double > m_DummyVector;
405 
408 
412  virtual void SetDefaultMRFNeighborhoodWeight();
413 
414  //Function implementing the ICM algorithm to label the images
415  void ApplyICMLabeller();
416 }; // class MRFImageFilter
417 } // namespace itk
418 
419 #ifndef ITK_MANUAL_INSTANTIATION
420 #include "itkMRFImageFilter.hxx"
421 #endif
422 
423 #endif
unsigned int m_NumberOfClasses
static const unsigned int ClassifiedImageDimension
Superclass::RegionType RegionType
Definition: itkImage.h:140
Image< int, itkGetStaticConstMacro(InputImageDimension) > LabelStatusImageType
TInputImage::RegionType InputImageRegionType
Light weight base class for most itk classes.
LabelledImageIndexType::IndexValueType IndexValueType
void PrintSelf(std::ostream &os, Indent indent) const override
LabelledImageFacesCalculator::FaceListType LabelledImageFaceListType
InputImageFaceListType::iterator InputImageFaceListIterator
NeighborhoodIterator< LabelStatusImageType > LabelStatusImageNeighborhoodIterator
TClassifiedImage::PixelType LabelledImagePixelType
LabelStatusImageType::RegionType LabelStatusRegionType
TClassifiedImage::IndexType LabelledImageIndexType
virtual void GenerateOutputInformation() override
std::vector< double > m_MahalanobisDistance
LabelStatusImageType::IndexType LabelStatusIndexType
signed long IndexValueType
Definition: itkIntTypes.h:150
SmartPointer< const Self > ConstPointer
InputImageFacesCalculator::FaceListType InputImageFaceListType
TInputImage::PixelType InputImagePixelType
ImageRegionIterator< TInputImage > InputImageRegionIterator
LabelStatusImageNeighborhoodIterator::RadiusType LabelStatusImageNeighborhoodRadiusType
void SetClassifier(typename ClassifierType::Pointer ptrToClassifier)
InputImageNeighborhoodIterator::RadiusType InputImageNeighborhoodRadiusType
NeighborhoodIterator< TClassifiedImage > LabelledImageNeighborhoodIterator
virtual void EnlargeOutputRequestedRegion(DataObject *) override
std::vector< double > m_MRFNeighborhoodWeight
TClassifiedImage::OffsetType LabelledImageOffsetType
virtual void GenerateData() override
LabelStatusImagePointer m_LabelStatusImage
LabelledImageNeighborhoodIterator::RadiusType LabelledImageNeighborhoodRadiusType
Base class for the ImageClassifierBase object.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
TClassifiedImage::Pointer TrainingImagePointer
TClassifiedImage::RegionType LabelledImageRegionType
std::vector< double > m_DummyVector
LabelStatusImageFaceListType::iterator LabelStatusImageFaceListIterator
Implementation of a labeller object that uses Markov Random Fields to classify pixels in an image dat...
ImageToImageFilter< TInputImage, TClassifiedImage > Superclass
unsigned int m_NumberOfIterations
virtual void SetMRFNeighborhoodWeight(std::vector< double > BetaMatrix)
unsigned int m_MaximumNumberOfIterations
ClassifierType::Pointer m_ClassifierPtr
ImageRegionIterator< TClassifiedImage > LabelledImageRegionIterator
Splits an image into a main region and several &quot;face&quot; regions which are used to handle computations o...
virtual void ApplyMRFImageFilter()
typedef(Concept::Convertible< unsigned int, LabelledImagePixelType >) UnsignedIntConvertibleToClassifiedCheck
TClassifiedImage::PixelType TrainingImagePixelType
OutputImageType::Pointer OutputImagePointer
ImageRegionConstIterator< TInputImage > InputImageRegionConstIterator
virtual void SetDefaultMRFNeighborhoodWeight()
ConstNeighborhoodIterator< TInputImage > InputImageNeighborhoodIterator
TClassifiedImage::Pointer LabelledImagePointer
A multi-dimensional iterator templated over image type that walks a region of pixels.
Superclass::IndexType IndexType
Definition: itkImage.h:122
virtual void DoNeighborhoodOperation(const InputImageNeighborhoodIterator &imageIter, LabelledImageNeighborhoodIterator &labelledIter, LabelStatusImageNeighborhoodIterator &labelStatusIter)
TInputImage::SizeType NeighborhoodRadiusType
LabelStatusImageNeighborhoodRadiusType m_LabelStatusImageNeighborhoodRadius
Superclass::OutputImagePointer OutputImagePointer
LabelledImageFaceListType::iterator LabelledImageFaceListIterator
SmartPointer< Self > Pointer
InputImageNeighborhoodRadiusType m_InputImageNeighborhoodRadius
TInputImage::SizeType SizeType
virtual std::vector< double > GetMRFNeighborhoodWeight()
NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TInputImage > InputImageFacesCalculator
Superclass::RadiusType RadiusType
const NeighborhoodRadiusType GetNeighborhoodRadius() const
virtual void GenerateInputRequestedRegion() override
TInputImage::Pointer InputImagePointer
LabelledImageNeighborhoodRadiusType m_LabelledImageNeighborhoodRadius
ImageClassifierBase< TInputImage, TClassifiedImage > ClassifierType
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
TInputImage::ConstPointer InputImageConstPointer
static const unsigned int InputImageDimension
std::vector< double > m_NeighborInfluence
LabelStatusImageFacesCalculator::FaceListType LabelStatusImageFaceListType
virtual void MinimizeFunctional()
StopConditionType m_StopCondition
Base class for most ITK classes.
Definition: itkObject.h:57
ImageRegionIterator< LabelStatusImageType > LabelStatusImageIterator
#define itkConceptMacro(name, concept)
void SetNeighborhoodRadius(const NeighborhoodRadiusType &)
NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TClassifiedImage > LabelledImageFacesCalculator
Base class for all data objects in ITK.
Defines iteration of a local N-dimensional neighborhood of pixels across an itk::Image.
Templated n-dimensional image class.
Definition: itkImage.h:75
A multi-dimensional iterator templated over image type that walks a region of pixels.
LabelStatusImageType::Pointer LabelStatusImagePointer
TInputImage::SizeType InputImageSizeType