OTB  10.0.0
Orfeo Toolbox
otbMRFSamplerRandomMAP.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 
22 #ifndef otbMRFSamplerRandomMAP_h
23 #define otbMRFSamplerRandomMAP_h
24 
25 #include "itkMersenneTwisterRandomVariateGenerator.h"
26 #include "otbMRFSampler.h"
27 
28 namespace otb
29 {
51 
52 template <class TInput1, class TInput2>
53 class ITK_EXPORT MRFSamplerRandomMAP : public MRFSampler<TInput1, TInput2>
54 {
55 public:
58  typedef itk::SmartPointer<Self> Pointer;
59  typedef itk::SmartPointer<const Self> ConstPointer;
60 
61  typedef typename Superclass::InputImageNeighborhoodIterator InputImageNeighborhoodIterator;
62  typedef typename Superclass::LabelledImageNeighborhoodIterator LabelledImageNeighborhoodIterator;
63  typedef typename Superclass::LabelledImagePixelType LabelledImagePixelType;
64  typedef typename Superclass::InputImagePixelType InputImagePixelType;
65  typedef typename Superclass::EnergyFidelityType EnergyFidelityType;
66  typedef typename Superclass::EnergyRegularizationType EnergyRegularizationType;
67  typedef typename Superclass::EnergyFidelityPointer EnergyFidelityPointer;
68  typedef typename Superclass::EnergyRegularizationPointer EnergyRegularizationPointer;
69  typedef itk::Statistics::MersenneTwisterRandomVariateGenerator RandomGeneratorType;
70 
71  itkNewMacro(Self);
72 
74 
75  void SetNumberOfClasses(const unsigned int nClasses) override
76  {
77  if ((nClasses != this->m_NumberOfClasses) || (m_EnergiesInvalid == true))
78  {
79  this->m_NumberOfClasses = nClasses;
80  if (m_Energy != nullptr)
81  free(m_Energy);
82  if (m_RepartitionFunction != nullptr)
83  free(m_RepartitionFunction);
84  m_Energy = (double*)calloc(this->m_NumberOfClasses, sizeof(double));
85  m_RepartitionFunction = (double*)calloc(this->m_NumberOfClasses, sizeof(double));
86  this->Modified();
87  }
88  }
89 
90  inline int Compute(const InputImageNeighborhoodIterator& itData, const LabelledImageNeighborhoodIterator& itRegul) override
91  {
92  if (this->m_NumberOfClasses == 0)
93  {
94  itkExceptionMacro(<< "NumberOfClasses has to be greater than 0.");
95  }
96 
97  this->m_EnergyBefore = this->m_EnergyFidelity->GetValue(itData, itRegul.GetCenterPixel());
98  this->m_EnergyBefore += this->m_Lambda * this->m_EnergyRegularization->GetValue(itRegul, itRegul.GetCenterPixel());
99 
100  // Try all possible value (how to be generic ?)
101  this->m_EnergyAfter = this->m_EnergyBefore; // default values to current one
102  this->m_Value = itRegul.GetCenterPixel();
103 
104  // Compute probability for each possibility
105  double totalProba = 0.0;
106  unsigned int valueCurrent = 0;
107  for (valueCurrent = 0; valueCurrent < this->m_NumberOfClasses; ++valueCurrent)
108  {
109  this->m_EnergyCurrent = this->m_EnergyFidelity->GetValue(itData, static_cast<LabelledImagePixelType>(valueCurrent));
110  this->m_EnergyCurrent += this->m_Lambda * this->m_EnergyRegularization->GetValue(itRegul, static_cast<LabelledImagePixelType>(valueCurrent));
111 
112  m_Energy[valueCurrent] = this->m_EnergyCurrent;
113  m_RepartitionFunction[valueCurrent] = std::exp(-this->m_EnergyCurrent) + totalProba;
114  totalProba = m_RepartitionFunction[valueCurrent];
115  }
116 
117  // Pick a value according to probability
118 
119  // double select = (m_Generator->GetIntegerVariate()/(double(RAND_MAX)+1) * totalProba);
120  double select = (m_Generator->GetIntegerVariate() / (double(itk::NumericTraits<RandomGeneratorType::IntegerType>::max()) + 1) * totalProba);
121  valueCurrent = 0;
122  while ((valueCurrent < this->GetNumberOfClasses()) && (m_RepartitionFunction[valueCurrent] <= select))
123  {
124  valueCurrent++;
125  }
126 
127  if (valueCurrent == this->GetNumberOfClasses())
128  {
129  valueCurrent = this->GetNumberOfClasses() - 1;
130  }
131 
132  if (this->m_Value != static_cast<LabelledImagePixelType>(valueCurrent))
133  {
134  this->m_Value = static_cast<LabelledImagePixelType>(valueCurrent);
135  this->m_EnergyAfter = m_Energy[static_cast<unsigned int>(valueCurrent)];
136  }
137 
138  this->m_DeltaEnergy = this->m_EnergyAfter - this->m_EnergyBefore;
139 
140  return 0;
141  }
142 
144  void InitializeSeed(int seed)
145  {
146  m_Generator->SetSeed(seed);
147  }
149  {
150  m_Generator->SetSeed();
151  }
153 
154 protected:
155  // The constructor and destructor.
156  MRFSamplerRandomMAP() : m_RepartitionFunction(nullptr), m_Energy(nullptr), m_EnergiesInvalid(true)
157  {
158  m_Generator = RandomGeneratorType::GetInstance();
159  m_Generator->SetSeed();
160  }
162  {
163  if (m_Energy != nullptr)
164  free(m_Energy);
165  if (m_RepartitionFunction != nullptr)
166  free(m_RepartitionFunction);
167  }
168 
169 private:
171  double* m_Energy;
173  RandomGeneratorType::Pointer m_Generator;
174 };
175 }
176 
177 #endif
This is the base class for sampler methods used in the MRF framework.
void SetNumberOfClasses(const unsigned int nClasses) override
otb::MRFSampler< TInput1, TInput2 > Superclass
itk::SmartPointer< Self > Pointer
Superclass::InputImagePixelType InputImagePixelType
Superclass::EnergyRegularizationType EnergyRegularizationType
Superclass::EnergyFidelityType EnergyFidelityType
Superclass::InputImageNeighborhoodIterator InputImageNeighborhoodIterator
itk::SmartPointer< const Self > ConstPointer
RandomGeneratorType::Pointer m_Generator
Superclass::EnergyFidelityPointer EnergyFidelityPointer
Superclass::EnergyRegularizationPointer EnergyRegularizationPointer
itk::Statistics::MersenneTwisterRandomVariateGenerator RandomGeneratorType
Superclass::LabelledImageNeighborhoodIterator LabelledImageNeighborhoodIterator
Superclass::LabelledImagePixelType LabelledImagePixelType
int Compute(const InputImageNeighborhoodIterator &itData, const LabelledImageNeighborhoodIterator &itRegul) override
This is the base class for sampler methods used in the MRF framework.
Definition: otbMRFSampler.h:45
itk::ConstNeighborhoodIterator< TInput1 > InputImageNeighborhoodIterator
Definition: otbMRFSampler.h:52
itk::NeighborhoodIterator< TInput2 > LabelledImageNeighborhoodIterator
Definition: otbMRFSampler.h:55
TInput2::PixelType LabelledImagePixelType
Definition: otbMRFSampler.h:56
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.