Skip to content
Snippets Groups Projects
Commit d2a06dc0 authored by Jonathan Guinet's avatar Jonathan Guinet
Browse files

ENH Optionnal attributes parameters methods addition for ConnectedComponentsegmentation class.

parent 1191c3c2
No related branches found
No related tags found
No related merge requests found
......@@ -81,13 +81,18 @@ public:
/** Get the compute perimeter flag */
bool GetComputePerimeter() const;
/** Set the compute perimeter flag */
void SetComputeFlusser(bool flag);
/** Get the compute perimeter flag */
bool GetComputeFlusser() const;
/** Set the polygonalisation flag */
void SetComputePolygon(bool flag);
/** Get the polygonalisation flag */
bool GetComputePolygon() const;
/** Set the compute feret diameter flag */
void SetComputeFeretDiameter(bool flag);
......@@ -140,6 +145,9 @@ private:
/** Do we compute the perimeter ? */
bool m_ComputePerimeter;
/** Do we compute flusser moments ? */
bool m_ComputeFlusser;
/** Do we polygonise ? */
bool m_ComputePolygon;
......@@ -232,13 +240,22 @@ public:
itkBooleanMacro(ComputePerimeter);
/**
* Set/Get whether the polygonalisation process should be computed or not. The defaut value
* Set/Get whether the polygonalisation process should be computed or not. The default value
* is true, to assure backward compatibility.
*/
void SetComputePolygon(bool flag);
bool GetComputePolygon() const;
itkBooleanMacro(ComputePolygon);
/**
* Set/Get whether the Flussrer moments should be computed or not. The default value
* is true, to assure backward compatibility.
*/
void SetComputeFlusser(bool flag);
bool GetComputeFlusser() const;
itkBooleanMacro(ComputeFlusser);
/** Set/get the ReducedAttributesSet flag */
void SetReducedAttributeSet(bool flag);
bool GetReducedAttributeSet() const;
......
......@@ -40,6 +40,7 @@ template <class TLabelObject, class TLabelImage>
ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
::ShapeAttributesLabelObjectFunctor() : m_ComputeFeretDiameter(false),
m_ComputePerimeter(false),
m_ComputeFlusser(true),
m_ComputePolygon(true),
m_ReducedAttributeSet(true),
m_PerimeterCalculator(NULL),
......@@ -93,7 +94,7 @@ ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
return m_ComputePerimeter;
}
/** Set the compute perimeter flag */
/** Set the compute polygon flag */
template <class TLabelObject, class TLabelImage>
void
ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
......@@ -102,7 +103,7 @@ ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
m_ComputePolygon = flag;
}
/** Get the compute perimeter flag */
/** Get the compute polygon flag */
template <class TLabelObject, class TLabelImage>
bool
ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
......@@ -111,6 +112,23 @@ ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
return m_ComputePolygon;
}
/** Set the compute Flusser flag */
template <class TLabelObject, class TLabelImage>
void
ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
::SetComputeFlusser(bool flag)
{
m_ComputeFlusser = flag;
}
/** Get the compute Flusser flag */
template <class TLabelObject, class TLabelImage>
bool
ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
::GetComputeFlusser() const
{
return m_ComputeFlusser;
}
/** Set the compute feret diameter flag */
template <class TLabelObject, class TLabelImage>
......@@ -463,59 +481,63 @@ ShapeAttributesLabelObjectFunctor<TLabelObject, TLabelImage>
ellipsoidSize[i] = 2.0 * equivalentRadius * vcl_sqrt(principalMoments[i] / edet);
}
// Flusser moments (only make sense when ImageDimension == 2)
if (LabelObjectType::ImageDimension == 2)
{
// complex centered moments
std::complex<double> c11, c20, c12, c21, c22, c30, c31, c40;
c11 = c20 = c12 = c21 = c22 = c30 = c31 = c40 = std::complex<double>(0, 0);
for (lit = lineContainer.begin(); lit != lineContainer.end(); lit++)
{
const typename LabelObjectType::IndexType& idx = lit->GetIndex();
unsigned long length = lit->GetLength();
//
long endIdx0 = idx[0] + length;
for (typename LabelObjectType::IndexType iidx = idx; iidx[0] < endIdx0; iidx[0]++)
if (m_ComputeFlusser)
{
// Flusser moments (only make sense when ImageDimension == 2)
if (LabelObjectType::ImageDimension == 2)
{
// complex centered moments
std::complex<double> c11, c20, c12, c21, c22, c30, c31, c40;
c11 = c20 = c12 = c21 = c22 = c30 = c31 = c40 = std::complex<double>(0, 0);
for (lit = lineContainer.begin(); lit != lineContainer.end(); lit++)
{
typename TLabelImage::PointType cPP;
m_LabelImage->TransformIndexToPhysicalPoint(iidx, cPP);
cPP -= physicalCentroid.GetVectorFromOrigin();
std::complex<double> xpiy(cPP[0], cPP[1]); // x + i y
std::complex<double> xmiy(cPP[0], -cPP[1]); // x - i y
c11 += xpiy * xmiy;
c20 += xpiy * xpiy;
c12 += xpiy * xmiy * xmiy;
c21 += xpiy * xpiy * xmiy;
c30 += xpiy * xpiy * xpiy;
c22 += xpiy * xpiy * xmiy * xmiy;
c31 += xpiy * xpiy * xpiy * xmiy;
c40 += xpiy * xpiy * xpiy * xpiy;
const typename LabelObjectType::IndexType& idx = lit->GetIndex();
unsigned long length = lit->GetLength();
//
long endIdx0 = idx[0] + length;
for (typename LabelObjectType::IndexType iidx = idx; iidx[0] < endIdx0; iidx[0]++)
{
typename TLabelImage::PointType cPP;
m_LabelImage->TransformIndexToPhysicalPoint(iidx, cPP);
cPP -= physicalCentroid.GetVectorFromOrigin();
std::complex<double> xpiy(cPP[0], cPP[1]); // x + i y
std::complex<double> xmiy(cPP[0], -cPP[1]); // x - i y
c11 += xpiy * xmiy;
c20 += xpiy * xpiy;
c12 += xpiy * xmiy * xmiy;
c21 += xpiy * xpiy * xmiy;
c30 += xpiy * xpiy * xpiy;
c22 += xpiy * xpiy * xmiy * xmiy;
c31 += xpiy * xpiy * xpiy * xmiy;
c40 += xpiy * xpiy * xpiy * xpiy;
}
}
}
// normalize
c11 /= physicalSize * physicalSize;
c20 /= physicalSize * physicalSize;
c12 /= vcl_pow(physicalSize, 5.0 / 2);
c21 /= vcl_pow(physicalSize, 5.0 / 2);
c30 /= vcl_pow(physicalSize, 5.0 / 2);
c22 /= vcl_pow(physicalSize, 3);
c31 /= vcl_pow(physicalSize, 3);
c40 /= vcl_pow(physicalSize, 3);
lo->SetAttribute("SHAPE::Flusser01", c11.real());
lo->SetAttribute("SHAPE::Flusser02", (c21 * c12).real());
lo->SetAttribute("SHAPE::Flusser03", (c20 * vcl_pow(c12, 2)).real());
lo->SetAttribute("SHAPE::Flusser04", (c20 * vcl_pow(c12, 2)).imag());
lo->SetAttribute("SHAPE::Flusser05", (c30 * vcl_pow(c12, 3)).real());
lo->SetAttribute("SHAPE::Flusser06", (c30 * vcl_pow(c12, 3)).imag());
lo->SetAttribute("SHAPE::Flusser07", c22.real());
lo->SetAttribute("SHAPE::Flusser08", (c31 * vcl_pow(c12, 2)).real());
lo->SetAttribute("SHAPE::Flusser09", (c31 * vcl_pow(c12, 2)).imag());
lo->SetAttribute("SHAPE::Flusser10", (c40 * vcl_pow(c12, 4)).real());
lo->SetAttribute("SHAPE::Flusser11", (c40 * vcl_pow(c12, 4)).imag());
// normalize
c11 /= physicalSize * physicalSize;
c20 /= physicalSize * physicalSize;
c12 /= vcl_pow(physicalSize, 5.0 / 2);
c21 /= vcl_pow(physicalSize, 5.0 / 2);
c30 /= vcl_pow(physicalSize, 5.0 / 2);
c22 /= vcl_pow(physicalSize, 3);
c31 /= vcl_pow(physicalSize, 3);
c40 /= vcl_pow(physicalSize, 3);
lo->SetAttribute("SHAPE::Flusser01", c11.real());
lo->SetAttribute("SHAPE::Flusser02", (c21 * c12).real());
lo->SetAttribute("SHAPE::Flusser03", (c20 * vcl_pow(c12, 2)).real());
lo->SetAttribute("SHAPE::Flusser04", (c20 * vcl_pow(c12, 2)).imag());
lo->SetAttribute("SHAPE::Flusser05", (c30 * vcl_pow(c12, 3)).real());
lo->SetAttribute("SHAPE::Flusser06", (c30 * vcl_pow(c12, 3)).imag());
lo->SetAttribute("SHAPE::Flusser07", c22.real());
lo->SetAttribute("SHAPE::Flusser08", (c31 * vcl_pow(c12, 2)).real());
lo->SetAttribute("SHAPE::Flusser09", (c31 * vcl_pow(c12, 2)).imag());
lo->SetAttribute("SHAPE::Flusser10", (c40 * vcl_pow(c12, 4)).real());
lo->SetAttribute("SHAPE::Flusser11", (c40 * vcl_pow(c12, 4)).imag());
}
}
// Set the attributes
......@@ -790,6 +812,27 @@ ShapeAttributesLabelMapFilter<TImage, TLabelImage>
return this->GetFunctor().GetComputePolygon();
}
template<class TImage, class TLabelImage>
void
ShapeAttributesLabelMapFilter<TImage, TLabelImage>
::SetComputeFlusser(bool flag)
{
if (this->GetFunctor().GetComputeFlusser() != flag)
{
this->GetFunctor().SetComputeFlusser(flag);
this->Modified();
}
}
template<class TImage, class TLabelImage>
bool
ShapeAttributesLabelMapFilter<TImage, TLabelImage>
::GetComputeFlusser() const
{
return this->GetFunctor().GetComputeFlusser();
}
template<class TImage, class TLabelImage>
void
......
......@@ -113,6 +113,7 @@ public:
typedef typename RelabelComponentFilterType::ObjectSizeType ObjectSizeType;
/* Set the mathematical expression used for the mask */
itkSetStringMacro(MaskExpression);
......@@ -137,6 +138,43 @@ public:
/* Get the mathematical expression for filtering labelobjects */
itkGetStringMacro(OBIAExpression);
/* Set shape reduced set attributes flag for object attributes computing */
itkSetMacro(ShapeReducedSetOfAttributes, bool);
/* Get shape reduced set attributes flag for object attributes computing */
itkGetMacro(ShapeReducedSetOfAttributes,bool);
/* Set stat reduced set attributes flag for object attributes computing */
itkSetMacro(StatsReducedSetOfAttributes,bool);
/* Get stat reduced set attributes flag for object attributes computing */
itkGetMacro(StatsReducedSetOfAttributes,bool);
/* Set compute polygon flag for object attributes computing */
itkSetMacro(ComputePolygon,bool);
/* Get compute polygon flag for object attributes computing */
itkGetMacro(ComputePolygon,bool);
/* Set compute Flusser flag for object attributes computing */
itkSetMacro(ComputeFlusser,bool);
/* Get compute Flusser flag for object attributes computing */
itkGetMacro(ComputeFlusser,bool);
/* Set compute perimeter flag for object attributes computing */
itkSetMacro(ComputePerimeter,bool);
/* Get compute perimeter flag for object attributes computing */
itkGetMacro(ComputePerimeter,bool);
/* Set compute feret diameter flag for object attributes computing */
itkSetMacro(ComputeFeretDiameter,bool);
/* Get compute FeretdDiameter flag for object attributes computing */
itkGetMacro(ComputeFeretDiameter,bool);
protected:
PersistentConnectedComponentSegmentationOBIAToVectorDataFilter();
......@@ -150,6 +188,14 @@ private:
std::string m_ConnectedComponentExpression;
std::string m_OBIAExpression;
// attributes
bool m_ShapeReducedSetOfAttributes;
bool m_StatsReducedSetOfAttributes;
bool m_ComputeFlusser;
bool m_ComputePolygon;
bool m_ComputeFeretDiameter;
bool m_ComputePerimeter;
virtual VectorDataPointerType ProcessTile();
};
......
......@@ -27,7 +27,8 @@ namespace otb {
template<class TVImage, class TLabelImage, class TMaskImage, class TOutputVectorData>
PersistentConnectedComponentSegmentationOBIAToVectorDataFilter<TVImage, TLabelImage, TMaskImage, TOutputVectorData>
::PersistentConnectedComponentSegmentationOBIAToVectorDataFilter()
: m_MinimumObjectSize(2)
: m_MinimumObjectSize(2),m_ShapeReducedSetOfAttributes(false),m_StatsReducedSetOfAttributes(false),
m_ComputeFlusser(false),m_ComputePolygon(false),m_ComputeFeretDiameter(false),m_ComputePerimeter(false)
{
}
......@@ -98,16 +99,17 @@ PersistentConnectedComponentSegmentationOBIAToVectorDataFilter<TVImage, TLabelIm
// shape attributes computation
typename ShapeLabelMapFilterType::Pointer shapeLabelMapFilter = ShapeLabelMapFilterType::New();
shapeLabelMapFilter->SetInput(labelImageToLabelMap->GetOutput());
shapeLabelMapFilter->SetReducedAttributeSet(true);
shapeLabelMapFilter->SetComputePolygon(false);
shapeLabelMapFilter->SetComputePerimeter(false);
shapeLabelMapFilter->SetComputeFeretDiameter(false);
shapeLabelMapFilter->SetReducedAttributeSet(m_ShapeReducedSetOfAttributes);
shapeLabelMapFilter->SetComputePolygon(m_ComputePolygon);
shapeLabelMapFilter->SetComputePerimeter(m_ComputePerimeter);
shapeLabelMapFilter->SetComputeFeretDiameter(m_ComputeFeretDiameter);
shapeLabelMapFilter->SetComputeFlusser(m_ComputeFlusser);
// band stat attributes computation
typename RadiometricLabelMapFilterType::Pointer radiometricLabelMapFilter = RadiometricLabelMapFilterType::New();
radiometricLabelMapFilter->SetInput(shapeLabelMapFilter->GetOutput());
radiometricLabelMapFilter->SetFeatureImage(extract->GetOutput());
radiometricLabelMapFilter->SetReducedAttributeSet(true);
radiometricLabelMapFilter->SetReducedAttributeSet(m_StatsReducedSetOfAttributes);
// OBIA Filtering using shape and radiometric object characteristics
typename LabelObjectOpeningFilterType::Pointer opening = LabelObjectOpeningFilterType::New();
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment