Skip to content
Snippets Groups Projects
Commit 5762304e authored by Emmanuel Christophe's avatar Emmanuel Christophe
Browse files

ENH: Adding ImageFittingPolygonList class

parent b37be464
No related branches found
No related tags found
No related merge requests found
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#ifndef __otbImagePerturbationPolygonListFilter_h
#define __otbImagePerturbationPolygonListFilter_h
#include "otbPathListToPathListFilter.h"
#include "otbMacro.h"
#include "itkLineConstIterator.h"
namespace otb
{
/** \class ImagePerturbationPolygonListFilter
* \brief Slightly deform polygon to reach higher enery from the image
*
* <br>Limitations:</br> This filter is currently working with integer position
* for the polygon vertices. It should be optimized for continuous positions.
*
*/
template <class TPath, class TImage>
class ITK_EXPORT ImagePerturbationPolygonListFilter
: public PathListToPathListFilter<TPath>
{
public:
/** Standard typedefs */
typedef ImagePerturbationPolygonListFilter Self;
typedef PathListToPathListFilter<TPath> Superclass;
typedef itk::SmartPointer<Self> Pointer;
typedef itk::SmartPointer<const Self> ConstPointer;
/** Type macro */
itkNewMacro(Self);
/** Creation through object factory macro */
itkTypeMacro(ImagePerturbationPolygonListFilter, PathListToPathListFilter);
/** Template parameters typedefs */
typedef typename Superclass::PathType PathType;
typedef typename Superclass::PathListType PathListType;
typedef typename Superclass::PathPointerType PathPointerType;
typedef typename PathListType::Pointer PathListPointerType;
typedef typename PathListType::ConstIterator IteratorType;
typedef typename PathType::VertexType VertexType;
typedef typename PathType::VertexListType VertexListType;
typedef typename VertexListType::ConstIterator VertexListConstIteratorType;
typedef double RealType;
typedef TImage ImageType;
typedef typename ImageType::Pointer ImagePointerType;
typedef typename ImageType::ConstPointer ImageConstPointerType;
typedef itk::LineConstIterator<ImageType> LineConstIteratorType;
/**
* Set the input Likelihood image.
* \param image The Likelihood image.
*/
void SetInputImage(const ImageType * image);
/**
* Get the input Likelihood image.
* \return The input Likelihood image.
*/
const ImageType * GetInputImage(void);
/** Set/Get the search radius. */
itkSetMacro(Radius,unsigned int);
itkGetMacro(Radius,unsigned int);
/** Set/Get the number of iteration. */
itkSetMacro(NumberOfIterations,unsigned int);
itkGetMacro(NumberOfIterations,unsigned int);
protected:
/** Constructor */
ImagePerturbationPolygonListFilter();
/** Destructor */
virtual ~ImagePerturbationPolygonListFilter() {};
/** GenerateData method */
virtual void GenerateData();
/** PrintSelf method */
virtual void PrintSelf(std::ostream& os, itk::Indent indent) const;
virtual double computeValue(ImageConstPointerType image, VertexType middlePoint, VertexType previousPoint, VertexType nextPoint) const;
unsigned int m_Radius;
unsigned int m_NumberOfIterations;
private:
ImagePerturbationPolygonListFilter(const Self&); //purposely not implemented
void operator=(const Self&); //purposely not implemented
};
}// End namespace otb
#ifndef OTB_MANUAL_INSTANTIATION
#include "otbImagePerturbationPolygonListFilter.txx"
#endif
#endif
/*=========================================================================
Program: ORFEO Toolbox
Language: C++
Date: $Date$
Version: $Revision$
Copyright (c) Centre National d'Etudes Spatiales. All rights reserved.
See OTBCopyright.txt for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notices for more information.
=========================================================================*/
#ifndef __otbImagePerturbationPolygonListFilter_txx
#define __otbImagePerturbationPolygonListFilter_txx
#include "otbImagePerturbationPolygonListFilter.h"
#include "otbPolyLineImageConstIterator.h"
#include "otbMacro.h"
namespace otb
{
/**
* Constructor
*/
template <class TPath, class TImage>
ImagePerturbationPolygonListFilter<TPath, TImage>
::ImagePerturbationPolygonListFilter()
{
this->SetNumberOfRequiredInputs(2);
this->SetNumberOfInputs(2);
m_Radius=1;
m_NumberOfIterations=1;
}
template <class TPath, class TImage>
void
ImagePerturbationPolygonListFilter<TPath, TImage>
::SetInputImage(const ImageType * image)
{
this->itk::ProcessObject::SetNthInput(1,const_cast<ImageType *>(image));
}
template <class TPath, class TImage>
const typename ImagePerturbationPolygonListFilter<TPath, TImage>
::ImageType *
ImagePerturbationPolygonListFilter<TPath, TImage>
::GetInputImage(void)
{
if(this->GetNumberOfInputs()<1)
{
return 0;
}
return static_cast<const ImageType *>(this->itk::ProcessObject::GetInput(1));
}
//FIXME
//There is an issue here with integer and continous indexes
//maybe we should use the itk::LineConstIterator
template <class TPath, class TImage>
void
ImagePerturbationPolygonListFilter<TPath, TImage>
::GenerateData()
{
// I/O wiring
ImageConstPointerType inputImagePtr = this->GetInputImage();
const PathListType * inputPtr = this->GetInput();
PathListType * outputPtr = this->GetOutput();
typename ImageType::RegionType regionLargest=inputImagePtr->GetLargestPossibleRegion();
typedef itk::ImageRegionConstIteratorWithIndex<ImageType> NeighborhoodIteratorType;
typename ImageType::SizeType size;
size[0]= 2*m_Radius+1;
size[1]= 2*m_Radius+1;
typename ImageType::RegionType region;
region.SetSize(size);
typename ImageType::IndexType start;
//go through all the polygons in the list
for(IteratorType it = inputPtr->Begin(); it != inputPtr->End(); ++it)
{
PathPointerType polygon = it.Get();
if(polygon->GetVertexList()->Size()>2)
{
for (int iteration=0;iteration < m_NumberOfIterations;++iteration)
{
PathPointerType newPolygon = PathType::New();
VertexListConstIteratorType vertexIt = polygon->GetVertexList()->Begin();
//We are now going to go through all the vertex, we won't start to process
// first as we need to know the last one for that.
VertexType firstPoint = vertexIt.Value();
VertexType previousPoint = vertexIt.Value();
++vertexIt;
VertexType currentPoint = vertexIt.Value();
++vertexIt;
while (vertexIt != polygon->GetVertexList()->End())
{
VertexType nextPoint=vertexIt.Value();
/** try all the possible neighbor for the current point
* to factorize
* */
{
start[0] = static_cast<long int>(currentPoint[0]-m_Radius);
start[1] = static_cast<long int>(currentPoint[1]-m_Radius);
region.SetIndex(start);
NeighborhoodIteratorType nIt(inputImagePtr, region);
double maxValue=0.0;
VertexType maxPoint = currentPoint;
for (nIt.GoToBegin();!nIt.IsAtEnd();++nIt)
{
if(regionLargest.IsInside(nIt.GetIndex()))
{
VertexType middlePoint=static_cast<VertexType>(nIt.GetIndex());
double currentValue = computeValue(inputImagePtr, middlePoint, previousPoint, nextPoint);
if (currentValue > maxValue)
{
maxValue=currentValue;
maxPoint=middlePoint;
}
}
}
currentPoint=maxPoint;
newPolygon->AddVertex(maxPoint);
}
/** End 'to factorize' */
++vertexIt;
previousPoint=currentPoint;
currentPoint=nextPoint;
}
//We now need to process the last and the first point
VertexType nextPoint = firstPoint;
/** try all the possible neighbor for the current point
* to factorize
* */
{
start[0] = static_cast<long int>(currentPoint[0]-m_Radius);
start[1] = static_cast<long int>(currentPoint[1]-m_Radius);
region.SetIndex(start);
NeighborhoodIteratorType nIt(inputImagePtr, region);
double maxValue=0.0;
VertexType maxPoint = currentPoint;
for (nIt.GoToBegin();!nIt.IsAtEnd();++nIt)
{
if(regionLargest.IsInside(nIt.GetIndex()))
{
VertexType middlePoint=static_cast<VertexType>(nIt.GetIndex());
double currentValue = computeValue(inputImagePtr, middlePoint, previousPoint, nextPoint);
if (currentValue > maxValue)
{
maxValue=currentValue;
maxPoint=middlePoint;
}
}
}
currentPoint=maxPoint;
newPolygon->AddVertex(maxPoint);
}
/** End 'to factorize' */
previousPoint = currentPoint;
currentPoint= firstPoint;
vertexIt = newPolygon->GetVertexList()->Begin();
nextPoint=vertexIt.Value();
/** try all the possible neighbor for the current point
* to factorize
* */
{
start[0] = static_cast<long int>(currentPoint[0]-m_Radius);
start[1] = static_cast<long int>(currentPoint[1]-m_Radius);
region.SetIndex(start);
NeighborhoodIteratorType nIt(inputImagePtr, region);
double maxValue=0.0;
VertexType maxPoint = currentPoint;
for (nIt.GoToBegin();!nIt.IsAtEnd();++nIt)
{
if(regionLargest.IsInside(nIt.GetIndex()))
{
VertexType middlePoint=static_cast<VertexType>(nIt.GetIndex());
double currentValue = computeValue(inputImagePtr, middlePoint, previousPoint, nextPoint);
if (currentValue > maxValue)
{
maxValue=currentValue;
maxPoint=middlePoint;
}
}
}
currentPoint=maxPoint;
newPolygon->AddVertex(maxPoint);
}
/** End 'to factorize' */
polygon = newPolygon;//prepare the next iteration
}
}
outputPtr->PushBack(polygon);
}//going through the polygon list
}
template <class TPath, class TImage>
double
ImagePerturbationPolygonListFilter<TPath, TImage>
::computeValue(ImageConstPointerType image, VertexType middlePoint, VertexType previousPoint, VertexType nextPoint) const
{
typedef typename ImageType::IndexType IndexType;
IndexType middleIndex;
IndexType previousIndex;
IndexType nextIndex;
middleIndex[0]=static_cast<long int>(middlePoint[0]);
middleIndex[1]=static_cast<long int>(middlePoint[1]);
previousIndex[0]=static_cast<long int>(previousPoint[0]);
previousIndex[1]=static_cast<long int>(previousPoint[1]);
nextIndex[0]=static_cast<long int>(nextPoint[0]);
nextIndex[1]=static_cast<long int>(nextPoint[1]);
double currentValue = 0.0;
unsigned int count = 0;
{//compute for first segment
LineConstIteratorType itLine(image, previousIndex, middleIndex);
while (!itLine.IsAtEnd())
{
currentValue += itLine.Get();
++count;
++itLine;
}
}
{//compute for second segment
LineConstIteratorType itLine(image, nextIndex, middleIndex);
while (!itLine.IsAtEnd())
{
currentValue += itLine.Get();
++count;
++itLine;
}
}
return currentValue/count;
}
/**
* PrintSelf Method
*/
template <class TPath, class TImage>
void
ImagePerturbationPolygonListFilter<TPath, TImage>
::PrintSelf(std::ostream& os, itk::Indent indent) const
{
Superclass::PrintSelf(os, indent);
}
} // End namespace otb
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment