1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376
|
/*
* Copyright (C) 2005-2020 Centre National d'Etudes Spatiales (CNES)
*
* This file is part of Orfeo Toolbox
*
* https://www.orfeo-toolbox.org/
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* http://www.apache.org/licenses/LICENSE-2.0
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*/
#ifndef otbMarkovRandomFieldFilter_hxx
#define otbMarkovRandomFieldFilter_hxx
#include "otbMarkovRandomFieldFilter.h"
namespace otb
{
template <class TInputImage, class TClassifiedImage>
MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::MarkovRandomFieldFilter(void)
: m_NumberOfClasses(0),
m_MaximumNumberOfIterations(50),
m_ErrorCounter(0),
m_ImageDeltaEnergy(0.0),
m_NeighborhoodRadius(1),
m_TotalNumberOfValidPixelsInOutputImage(1),
m_TotalNumberOfPixelsInInputImage(1),
m_ErrorTolerance(0.0),
m_SmoothingFactor(1.0),
m_NumberOfIterations(0),
m_Lambda(1.0),
m_ExternalClassificationSet(false),
m_StopCondition(MaximumNumberOfIterations)
{
m_Generator = RandomGeneratorType::GetInstance();
m_Generator->SetSeed();
this->SetNumberOfRequiredInputs(1);
if ((int)InputImageDimension != (int)ClassifiedImageDimension)
{
std::ostringstream msg;
msg << "Input image dimension: " << InputImageDimension << " != output image dimension: " << ClassifiedImageDimension;
throw itk::ExceptionObject(__FILE__, __LINE__, msg.str(), ITK_LOCATION);
}
m_InputImageNeighborhoodRadius.Fill(m_NeighborhoodRadius);
// m_MRFNeighborhoodWeight.resize(0);
// m_NeighborInfluence.resize(0);
// m_DummyVector.resize(0);
// this->SetMRFNeighborhoodWeight( m_DummyVector );
// this->SetDefaultMRFNeighborhoodWeight();
// srand((unsigned)time(0));
}
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::SetTrainingInput(const TrainingImageType* trainingImage)
{
m_ExternalClassificationSet = true;
// Process object is not const-correct so the const_cast is required here
this->itk::ProcessObject::SetNthInput(1, const_cast<TrainingImageType*>(trainingImage));
this->Modified();
}
template <class TInputImage, class TClassifiedImage>
const typename MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::TrainingImageType*
MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::GetTrainingInput(void)
{
if (this->GetNumberOfInputs() < 2)
{
return nullptr;
}
return static_cast<const TrainingImageType*>(this->itk::ProcessObject::GetInput(1));
}
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
{
Superclass::PrintSelf(os, indent);
os << indent << " MRF Image filter object " << std::endl;
os << indent << " Number of classes: " << m_NumberOfClasses << std::endl;
os << indent << " Maximum number of iterations: " << m_MaximumNumberOfIterations << std::endl;
os << indent << " Error tolerance for convergence: " << m_ErrorTolerance << std::endl;
os << indent << " Size of the MRF neighborhood radius:" << m_InputImageNeighborhoodRadius << std::endl;
os << indent << "StopCondition: " << m_StopCondition << std::endl;
os << indent << " Number of iterations: " << m_NumberOfIterations << std::endl;
os << indent << " Lambda: " << m_Lambda << std::endl;
} // end PrintSelf
/**
* GenerateInputRequestedRegion method.
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::GenerateInputRequestedRegion()
{
// this filter requires the all of the input images
// to be at the size of the output requested region
InputImagePointer inputPtr = const_cast<InputImageType*>(this->GetInput());
OutputImagePointer outputPtr = this->GetOutput();
inputPtr->SetRequestedRegion(outputPtr->GetRequestedRegion());
}
/**
* EnlargeOutputRequestedRegion method.
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::EnlargeOutputRequestedRegion(itk::DataObject* output)
{
// this filter requires the all of the output image to be in
// the buffer
TClassifiedImage* imgData;
imgData = dynamic_cast<TClassifiedImage*>(output);
imgData->SetRequestedRegionToLargestPossibleRegion();
}
/**
* GenerateOutputInformation method.
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::GenerateOutputInformation()
{
typename TInputImage::ConstPointer input = this->GetInput();
typename TClassifiedImage::Pointer output = this->GetOutput();
output->SetLargestPossibleRegion(input->GetLargestPossibleRegion());
}
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::GenerateData()
{
// InputImageConstPointer inputImage = this->GetInput();
// Allocate memory for the labelled images
this->Allocate();
// Branch the pipeline
this->Initialize();
// Run the Markov random field
this->ApplyMarkovRandomFieldFilter();
} // end GenerateData
/**
* Set the neighborhood radius from a single value
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::SetNeighborhoodRadius(const unsigned long radiusValue)
{
// Set up the neighbor hood
NeighborhoodRadiusType radius;
for (unsigned int i = 0; i < InputImageDimension; ++i)
{
radius[i] = radiusValue;
}
this->SetNeighborhoodRadius(radius);
} // end SetNeighborhoodRadius
/**
* Set the neighborhood radius from an array
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::SetNeighborhoodRadius(const unsigned long* radiusArray)
{
NeighborhoodRadiusType radius;
for (unsigned int i = 0; i < InputImageDimension; ++i)
{
radius[i] = radiusArray[i];
}
// Set up the neighbor hood
this->SetNeighborhoodRadius(radius);
} // end SetNeighborhoodRadius
/**
* Set the neighborhood radius from a radius object
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::SetNeighborhoodRadius(const NeighborhoodRadiusType& radius)
{
// Set up the neighbor hood
for (unsigned int i = 0; i < InputImageDimension; ++i)
{
m_InputImageNeighborhoodRadius[i] = radius[i];
m_LabelledImageNeighborhoodRadius[i] = radius[i];
}
} // end SetNeighborhoodRadius
//-------------------------------------------------------
/**
* Allocate algorithm specific resources
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::Allocate()
{
// Set the output labelled and allocate the memory
LabelledImagePointer outputPtr = this->GetOutput();
// Allocate the output buffer memory
outputPtr->SetBufferedRegion(outputPtr->GetRequestedRegion());
outputPtr->Allocate();
// Copy input data in the output buffer memory or
// initialize to random values if not set
LabelledImageRegionIterator outImageIt(outputPtr, outputPtr->GetRequestedRegion());
if (m_ExternalClassificationSet)
{
typename TrainingImageType::ConstPointer trainingImage = this->GetTrainingInput();
LabelledImageRegionConstIterator trainingImageIt(trainingImage, outputPtr->GetRequestedRegion());
while (!outImageIt.IsAtEnd())
{
LabelledImagePixelType labelvalue = static_cast<LabelledImagePixelType>(trainingImageIt.Get());
outImageIt.Set(labelvalue);
++trainingImageIt;
++outImageIt;
} // end while
}
else // set to random value
{
// srand((unsigned)time(0));
while (!outImageIt.IsAtEnd())
{
LabelledImagePixelType randomvalue = static_cast<LabelledImagePixelType>(m_Generator->GetIntegerVariate() % static_cast<int>(m_NumberOfClasses));
outImageIt.Set(randomvalue);
++outImageIt;
} // end while
}
} // Allocate
/**
* Initialize pipeline and values
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::Initialize()
{
m_ImageDeltaEnergy = 0.0;
InputImageSizeType inputImageSize = this->GetInput()->GetBufferedRegion().GetSize();
//---------------------------------------------------------------------
// Get the number of valid pixels in the output MRF image
//---------------------------------------------------------------------
m_TotalNumberOfPixelsInInputImage = 1;
m_TotalNumberOfValidPixelsInOutputImage = 1;
for (unsigned int i = 0; i < InputImageDimension; ++i)
{
m_TotalNumberOfPixelsInInputImage *= static_cast<int>(inputImageSize[i]);
m_TotalNumberOfValidPixelsInOutputImage *= (static_cast<int>(inputImageSize[i]) - 2 * m_InputImageNeighborhoodRadius[i]);
}
srand((unsigned)time(nullptr));
if (!m_EnergyRegularization)
{
itkExceptionMacro(<< "EnergyRegularization is not present");
}
if (!m_EnergyFidelity)
{
itkExceptionMacro(<< "EnergyFidelity is not present");
}
if (!m_Optimizer)
{
itkExceptionMacro(<< "Optimizer is not present");
}
if (!m_Sampler)
{
itkExceptionMacro(<< "Sampler is not present");
}
m_Sampler->SetLambda(m_Lambda);
m_Sampler->SetEnergyRegularization(m_EnergyRegularization);
m_Sampler->SetEnergyFidelity(m_EnergyFidelity);
m_Sampler->SetNumberOfClasses(m_NumberOfClasses);
}
/**
*Apply the MRF image filter
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::ApplyMarkovRandomFieldFilter()
{
// Note: error should be defined according to the number of valid pixel in the output
int maxNumPixelError = itk::Math::Round<int, double>(m_ErrorTolerance * m_TotalNumberOfPixelsInInputImage);
m_NumberOfIterations = 0;
m_ErrorCounter = m_TotalNumberOfValidPixelsInOutputImage;
while ((m_NumberOfIterations < m_MaximumNumberOfIterations) && (m_ErrorCounter >= maxNumPixelError))
{
otbMsgDevMacro(<< "Iteration No." << m_NumberOfIterations);
this->MinimizeOnce();
otbMsgDevMacro(<< "m_ErrorCounter/m_TotalNumberOfPixelsInInputImage: " << m_ErrorCounter / ((double)(m_TotalNumberOfPixelsInInputImage)));
otbMsgDevMacro(<< "m_ImageDeltaEnergy: " << m_ImageDeltaEnergy);
++m_NumberOfIterations;
}
otbMsgDevMacro(<< "m_NumberOfIterations: " << m_NumberOfIterations);
otbMsgDevMacro(<< "m_MaximumNumberOfIterations: " << m_MaximumNumberOfIterations);
otbMsgDevMacro(<< "m_ErrorCounter: " << m_ErrorCounter);
otbMsgDevMacro(<< "maxNumPixelError: " << maxNumPixelError);
// Determine stop condition
if (m_NumberOfIterations >= m_MaximumNumberOfIterations)
{
m_StopCondition = MaximumNumberOfIterations;
}
else if (m_ErrorCounter <= maxNumPixelError)
{
m_StopCondition = ErrorTolerance;
}
} // ApplyMarkovRandomFieldFilter
/**
*Apply the MRF image filter on the whole image once
*/
template <class TInputImage, class TClassifiedImage>
void MarkovRandomFieldFilter<TInputImage, TClassifiedImage>::MinimizeOnce()
{
LabelledImageNeighborhoodIterator labelledIterator(m_LabelledImageNeighborhoodRadius, this->GetOutput(), this->GetOutput()->GetLargestPossibleRegion());
InputImageNeighborhoodIterator dataIterator(m_InputImageNeighborhoodRadius, this->GetInput(), this->GetInput()->GetLargestPossibleRegion());
m_ErrorCounter = 0;
for (labelledIterator.GoToBegin(), dataIterator.GoToBegin(); !labelledIterator.IsAtEnd(); ++labelledIterator, ++dataIterator)
{
LabelledImagePixelType value;
bool changeValueBool;
m_Sampler->Compute(dataIterator, labelledIterator);
value = m_Sampler->GetValue();
changeValueBool = m_Optimizer->Compute(m_Sampler->GetDeltaEnergy());
if (changeValueBool)
{
labelledIterator.SetCenterPixel(value);
++m_ErrorCounter;
m_ImageDeltaEnergy += m_Sampler->GetDeltaEnergy();
}
}
}
} // namespace otb
#endif
|