OTB  10.0.0
Orfeo Toolbox
otbMarkovRandomFieldFilter.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbMarkovRandomFieldFilter_h
22 #define otbMarkovRandomFieldFilter_h
23 
24 #include "otbMacro.h"
25 
26 #include "vnl/vnl_vector.h"
27 #include "vnl/vnl_matrix.h"
28 #include "itkMersenneTwisterRandomVariateGenerator.h"
29 
30 #include "itkImageToImageFilter.h"
31 #include "itkImageRegionIterator.h"
32 
33 #include "itkNeighborhoodAlgorithm.h"
34 #include "itkNeighborhood.h"
35 #include "itkSize.h"
36 #include "otbMRFOptimizer.h"
37 #include "otbMRFSampler.h"
38 
39 namespace otb
40 {
84 
85 template <class TInputImage, class TClassifiedImage>
86 class ITK_EXPORT MarkovRandomFieldFilter : public itk::ImageToImageFilter<TInputImage, TClassifiedImage>
87 {
88 public:
91  typedef itk::ImageToImageFilter<TInputImage, TClassifiedImage> Superclass;
92  typedef itk::SmartPointer<Self> Pointer;
93  typedef itk::SmartPointer<const Self> ConstPointer;
94  typedef typename Superclass::OutputImagePointer OutputImagePointer;
95 
97  itkTypeMacro(MarkovRandomFieldFilter, itk::ImageToImageFilter);
98 
100  itkNewMacro(Self);
101 
103  typedef TInputImage InputImageType;
104  typedef typename TInputImage::Pointer InputImagePointer;
105  typedef typename TInputImage::ConstPointer InputImageConstPointer;
106 
108  typedef typename TInputImage::PixelType InputImagePixelType;
109 
111  typedef typename TInputImage::RegionType InputImageRegionType;
112 
114  typedef itk::ImageRegionIterator<TInputImage> InputImageRegionIterator;
115  typedef itk::ImageRegionConstIterator<TInputImage> InputImageRegionConstIterator;
116 
118  itkStaticConstMacro(InputImageDimension, unsigned int, TInputImage::ImageDimension);
119 
121  typedef TClassifiedImage TrainingImageType;
122  typedef typename TClassifiedImage::Pointer TrainingImagePointer;
123 
125  typedef typename TClassifiedImage::PixelType TrainingImagePixelType;
126 
129  typedef TClassifiedImage LabelledImageType;
130  typedef typename TClassifiedImage::Pointer LabelledImagePointer;
131 
134  typedef typename TClassifiedImage::PixelType LabelledImagePixelType;
135 
138  typedef typename TClassifiedImage::RegionType LabelledImageRegionType;
139 
141  typedef typename TClassifiedImage::IndexType LabelledImageIndexType;
142  typedef typename LabelledImageIndexType::IndexValueType IndexValueType;
143 
145  typedef typename TClassifiedImage::OffsetType LabelledImageOffsetType;
146 
148  typedef itk::ImageRegionIterator<TClassifiedImage> LabelledImageRegionIterator;
149 
150  typedef itk::ImageRegionConstIterator<TClassifiedImage> LabelledImageRegionConstIterator;
151 
153  itkStaticConstMacro(ClassifiedImageDimension, unsigned int, TClassifiedImage::ImageDimension);
154 
156  typedef typename TInputImage::SizeType SizeType;
157 
159  typedef typename TInputImage::SizeType NeighborhoodRadiusType;
160 
162  typedef itk::ConstNeighborhoodIterator<TInputImage> InputImageNeighborhoodIterator;
163 
164  typedef typename InputImageNeighborhoodIterator::RadiusType InputImageNeighborhoodRadiusType;
165 
166  typedef itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TInputImage> InputImageFacesCalculator;
167 
168  typedef typename InputImageFacesCalculator::FaceListType InputImageFaceListType;
169 
170  typedef typename InputImageFaceListType::iterator InputImageFaceListIterator;
171 
173  typedef itk::NeighborhoodIterator<TClassifiedImage> LabelledImageNeighborhoodIterator;
174 
175  typedef typename LabelledImageNeighborhoodIterator::RadiusType LabelledImageNeighborhoodRadiusType;
176 
177  typedef itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<TClassifiedImage> LabelledImageFacesCalculator;
178 
179  typedef typename LabelledImageFacesCalculator::FaceListType LabelledImageFaceListType;
180 
181  typedef typename LabelledImageFaceListType::iterator LabelledImageFaceListIterator;
182 
184  typedef itk::Statistics::MersenneTwisterRandomVariateGenerator RandomGeneratorType;
185 
189 
192 
195 
198 
202  itkSetObjectMacro(EnergyRegularization, EnergyRegularizationType);
205 
206  itkSetObjectMacro(EnergyFidelity, EnergyFidelityType);
208 
209  itkSetObjectMacro(Sampler, SamplerType);
211 
212  itkSetObjectMacro(Optimizer, OptimizerType);
214 
216  itkSetMacro(NumberOfClasses, unsigned int);
217  itkGetMacro(NumberOfClasses, unsigned int);
219 
222  itkSetMacro(MaximumNumberOfIterations, unsigned int);
223  itkGetMacro(MaximumNumberOfIterations, unsigned int);
225 
228  itkSetMacro(ErrorTolerance, double);
229  itkGetMacro(ErrorTolerance, double);
231 
234  itkSetMacro(SmoothingFactor, double);
235  itkGetMacro(SmoothingFactor, double);
237 
240  itkSetMacro(Lambda, double);
241  itkGetMacro(Lambda, double);
243 
245  void SetNeighborhoodRadius(const NeighborhoodRadiusType&);
246 
250  void SetNeighborhoodRadius(const unsigned long);
251  void SetNeighborhoodRadius(const unsigned long* radiusArray);
252 
255  {
256  NeighborhoodRadiusType neighborhoodRadius;
257 
258  for (int i = 0; i < InputImageDimension; ++i)
259  neighborhoodRadius[i] = m_InputImageNeighborhoodRadius[i];
260 
261  return neighborhoodRadius;
262  }
263 
271  virtual void SetTrainingInput(const TrainingImageType* trainingImage);
272  const TrainingImageType* GetTrainingInput(void);
274 
275  // Enum to get the stopping condition of the MRF filter
276  typedef enum { MaximumNumberOfIterations = 1, ErrorTolerance } StopConditionType;
277 
280  itkGetConstReferenceMacro(StopCondition, StopConditionType);
281 
283  itkGetConstReferenceMacro(NumberOfIterations, unsigned int);
284 
285 #ifdef ITK_USE_CONCEPT_CHECKING
287  itkConceptMacro(UnsignedIntConvertibleToClassifiedCheck, (itk::Concept::Convertible<unsigned int, LabelledImagePixelType>));
288  itkConceptMacro(ClassifiedConvertibleToUnsignedIntCheck, (itk::Concept::Convertible<LabelledImagePixelType, unsigned int>));
289  itkConceptMacro(ClassifiedConvertibleToIntCheck, (itk::Concept::Convertible<LabelledImagePixelType, int>));
290  itkConceptMacro(IntConvertibleToClassifiedCheck, (itk::Concept::Convertible<int, LabelledImagePixelType>));
291  itkConceptMacro(SameDimensionCheck, (itk::Concept::SameDimension<InputImageDimension, ClassifiedImageDimension>));
292 
294 #endif
295 
297  void InitializeSeed(int seed)
298  {
299  m_Generator->SetSeed(seed);
300  }
302  {
303  m_Generator->SetSeed();
304  }
306 
307 protected:
310  {
311  }
312  void PrintSelf(std::ostream& os, itk::Indent indent) const override;
313 
317  void Allocate();
318 
322  void Initialize();
323 
324  virtual void ApplyMarkovRandomFieldFilter();
325 
326  void GenerateData() override;
327  void GenerateInputRequestedRegion() override;
328  void EnlargeOutputRequestedRegion(itk::DataObject*) override;
329  void GenerateOutputInformation() override;
330 
331  MarkovRandomFieldFilter(const Self&) = delete;
332  void operator=(const Self&) = delete;
333 
334  typedef typename TInputImage::SizeType InputImageSizeType;
335 
338 
339  unsigned int m_NumberOfClasses;
341 
344 
350 
351  unsigned int m_NumberOfIterations;
352 
353  double m_Lambda;
356 
358 
359  std::vector<double> m_MRFNeighborhoodWeight;
360  std::vector<double> m_NeighborInfluence;
361  std::vector<double> m_DummyVector;
362 
363  RandomGeneratorType::Pointer m_Generator;
364 
371 
372  virtual void MinimizeOnce();
373 
374 private:
375 }; // class MarkovRandomFieldFilter
376 
377 } // namespace otb
378 
379 #ifndef OTB_MANUAL_INSTANTIATION
381 #endif
382 
383 #endif
This is the base class for energy function used in the MRF framework.
Definition: otbMRFEnergy.h:44
itk::SmartPointer< Self > Pointer
Definition: otbMRFEnergy.h:48
This is the base class for optimizer used in the MRF framework.
itk::SmartPointer< Self > Pointer
This is the base class for sampler methods used in the MRF framework.
Definition: otbMRFSampler.h:45
itk::SmartPointer< Self > Pointer
Definition: otbMRFSampler.h:49
This is the class to use the Markov Random Field framework in OTB.
LabelledImageNeighborhoodIterator::RadiusType LabelledImageNeighborhoodRadiusType
std::vector< double > m_MRFNeighborhoodWeight
EnergyRegularizationPointer m_EnergyRegularization
itkGetObjectMacro(EnergyFidelity, EnergyFidelityType)
InputImageFaceListType::iterator InputImageFaceListIterator
TInputImage::RegionType InputImageRegionType
LabelledImageNeighborhoodRadiusType m_LabelledImageNeighborhoodRadius
LabelledImageIndexType::IndexValueType IndexValueType
TInputImage::PixelType InputImagePixelType
itkGetObjectMacro(Sampler, SamplerType)
itk::ConstNeighborhoodIterator< TInputImage > InputImageNeighborhoodIterator
itk::ImageRegionIterator< TInputImage > InputImageRegionIterator
itkGetObjectMacro(Optimizer, OptimizerType)
LabelledImageFaceListType::iterator LabelledImageFaceListIterator
itk::Statistics::MersenneTwisterRandomVariateGenerator RandomGeneratorType
itkGetObjectMacro(EnergyRegularization, EnergyRegularizationType)
MRFSampler< TInputImage, TClassifiedImage > SamplerType
TClassifiedImage::Pointer TrainingImagePointer
EnergyRegularizationType::Pointer EnergyRegularizationPointer
TClassifiedImage::PixelType TrainingImagePixelType
TInputImage::ConstPointer InputImageConstPointer
itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TInputImage > InputImageFacesCalculator
InputImageNeighborhoodRadiusType m_InputImageNeighborhoodRadius
EnergyFidelityType::Pointer EnergyFidelityPointer
MRFEnergy< TClassifiedImage, TClassifiedImage > EnergyRegularizationType
RandomGeneratorType::Pointer m_Generator
itk::ImageToImageFilter< TInputImage, TClassifiedImage > Superclass
TClassifiedImage::PixelType LabelledImagePixelType
TClassifiedImage::Pointer LabelledImagePointer
LabelledImageFacesCalculator::FaceListType LabelledImageFaceListType
const NeighborhoodRadiusType GetNeighborhoodRadius() const
MRFEnergy< TInputImage, TClassifiedImage > EnergyFidelityType
Superclass::OutputImagePointer OutputImagePointer
void operator=(const Self &)=delete
TClassifiedImage::IndexType LabelledImageIndexType
itk::SmartPointer< const Self > ConstPointer
InputImageFacesCalculator::FaceListType InputImageFaceListType
TInputImage::SizeType NeighborhoodRadiusType
TClassifiedImage::RegionType LabelledImageRegionType
itk::ImageRegionConstIterator< TInputImage > InputImageRegionConstIterator
TClassifiedImage::OffsetType LabelledImageOffsetType
itk::NeighborhoodIterator< TClassifiedImage > LabelledImageNeighborhoodIterator
MarkovRandomFieldFilter(const Self &)=delete
InputImageNeighborhoodIterator::RadiusType InputImageNeighborhoodRadiusType
itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator< TClassifiedImage > LabelledImageFacesCalculator
itk::ImageRegionIterator< TClassifiedImage > LabelledImageRegionIterator
itk::ImageRegionConstIterator< TClassifiedImage > LabelledImageRegionConstIterator
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.