DGtal  1.4.beta
DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder > Class Template Reference

Aim: This functor estimates normal vector for a collection of surfels using spherical accumulator based Hough voting. More...

#include <DGtal/geometry/surfaces/estimation/estimationFunctors/SphericalHoughNormalVectorEstimator.h>

Public Types

typedef TSurfel Surfel
 
typedef TEmbedder SCellEmbedder
 
typedef SCellEmbedder::RealPoint RealPoint
 
typedef RealPoint Quantity
 
typedef SimpleMatrix< double, 3, 3 > Matrix
 

Public Member Functions

 SphericalHoughNormalVectorEstimator (ConstAlias< SCellEmbedder > anEmbedder, const double h, const double minimalAspectRatio=0.001, const unsigned int nbTrials=100, const unsigned int accumulatorSize=10, const unsigned int nbAccumulators=5)
 
 SphericalHoughNormalVectorEstimator ()=delete
 
void pushSurfel (const Surfel &aSurf, const double aDistance)
 
Quantity eval ()
 
void reset ()
 

Private Member Functions

Matrix randomRotation () const
 
Quantity getNormal (const unsigned int i, const unsigned int j, const unsigned int k, double &aspect) const
 

Private Attributes

const SCellEmbeddermyEmbedder
 Alias of the geometrical embedder. More...
 
const double myH
 Grid step. More...
 
const double myAspectRatio
 Minimal aspect ratio (norm of the cross-product) to consider a given triangle. More...
 
const unsigned int myNbTrials
 Number of trials in the neignborhood. More...
 
const unsigned int mySize
 Size of the accumulator. More...
 
const unsigned int myNbAccumulators
 Number of randomly shifted spherical accumulators to consider. More...
 
std::vector< RealPointmyPoints
 vector of embedded surfels More...
 
std::vector< SphericalAccumulator< RealPoint > > myAccumulators
 Spherical Accumulators. More...
 
std::vector< MatrixmyRotations
 Random rotations. More...
 
std::vector< MatrixmyInverseRotations
 Random inverse rotations. More...
 

Detailed Description

template<typename TSurfel, typename TEmbedder>
class DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >

Aim: This functor estimates normal vector for a collection of surfels using spherical accumulator based Hough voting.

Description of template class 'SphericalHoughNormalVectorEstimator'

This functors implements [16] algorithm:

  • we first collect the surfels using the pushSurfel method;
  • we select random triples of surfels and estimate the normal vectors of the associated triangles;
  • each normal vector is added to a spherical accumulator (see SphericalAccumulator);
  • the estimated normal vector is computed from normal vectors of the bin of the accumulator with maximal vote.

To avoid aliasing artefacts of the spherical accumulator, several randomly rotated accumulators are used.

Given a random triple of surfels, a threshold on the triangle aspect ratio can be specified to discard bad aspect triangles (e.g. thin ones).

This functor is a model of concepts::CLocalEstimatorFromSurfelFunctor

Template Parameters
TSurfeltype of surfels
TEmbeddertype of functors which embed surfel to \( \mathbb{R}^3\)

Definition at line 84 of file SphericalHoughNormalVectorEstimator.h.

Member Typedef Documentation

◆ Matrix

template<typename TSurfel , typename TEmbedder >
typedef SimpleMatrix<double,3,3> DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::Matrix

Definition at line 92 of file SphericalHoughNormalVectorEstimator.h.

◆ Quantity

template<typename TSurfel , typename TEmbedder >
typedef RealPoint DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::Quantity

Definition at line 91 of file SphericalHoughNormalVectorEstimator.h.

◆ RealPoint

template<typename TSurfel , typename TEmbedder >
typedef SCellEmbedder::RealPoint DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::RealPoint

Definition at line 90 of file SphericalHoughNormalVectorEstimator.h.

◆ SCellEmbedder

template<typename TSurfel , typename TEmbedder >
typedef TEmbedder DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::SCellEmbedder

Definition at line 89 of file SphericalHoughNormalVectorEstimator.h.

◆ Surfel

template<typename TSurfel , typename TEmbedder >
typedef TSurfel DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::Surfel

Definition at line 88 of file SphericalHoughNormalVectorEstimator.h.

Constructor & Destructor Documentation

◆ SphericalHoughNormalVectorEstimator() [1/2]

template<typename TSurfel , typename TEmbedder >
DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::SphericalHoughNormalVectorEstimator ( ConstAlias< SCellEmbedder anEmbedder,
const double  h,
const double  minimalAspectRatio = 0.001,
const unsigned int  nbTrials = 100,
const unsigned int  accumulatorSize = 10,
const unsigned int  nbAccumulators = 5 
)
inline

Constructor.

Parameters
[in]anEmbedderembedder to map surfel to R^n.
[in]hgrid step
[in]minimalAspectRatiothe minimal aspect ratio of triangles to be considered in the accumulator (default=0.001).
[in]nbTrialsnumber of random triangles in the neighborhood to consider (default=100).
[in]accumulatorSizesize of the spherical accumulators (see SphericalAccumulator, default=10).
[in]nbAccumulatorsnumber of randomly rotated spherical accumulators to consider in order to avoid aliasing artefacts. (default=5)

Definition at line 104 of file SphericalHoughNormalVectorEstimator.h.

109  :
110  myEmbedder(&anEmbedder),myH(h), myAspectRatio(minimalAspectRatio),
111  myNbTrials( nbTrials), mySize(accumulatorSize) , myNbAccumulators(nbAccumulators)
112  {
113  SphericalAccumulator<RealPoint> accum(mySize);
114 
115  //We precompute the random rotations and accumulators
116  for(auto i = 0u; i < myNbAccumulators; ++i)
117  {
118  Matrix m = randomRotation();
119  myAccumulators.push_back( accum );
120  myRotations.push_back( m );
121  myInverseRotations.push_back( m.inverse() );
122  }
123  }
std::vector< SphericalAccumulator< RealPoint > > myAccumulators
Spherical Accumulators.
std::vector< Matrix > myInverseRotations
Random inverse rotations.
const SCellEmbedder * myEmbedder
Alias of the geometrical embedder.
const unsigned int myNbAccumulators
Number of randomly shifted spherical accumulators to consider.
const unsigned int myNbTrials
Number of trials in the neignborhood.
const double myAspectRatio
Minimal aspect ratio (norm of the cross-product) to consider a given triangle.

References DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myAccumulators, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myInverseRotations, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myNbAccumulators, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myRotations, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::mySize, and DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::randomRotation().

◆ SphericalHoughNormalVectorEstimator() [2/2]

template<typename TSurfel , typename TEmbedder >
DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::SphericalHoughNormalVectorEstimator ( )
delete

Disable default constructor.

Member Function Documentation

◆ eval()

template<typename TSurfel , typename TEmbedder >
Quantity DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::eval ( )
inline

Estimate normal vector using spherical accumulator voting.

Returns
the feature score

Definition at line 152 of file SphericalHoughNormalVectorEstimator.h.

153  {
154  std::default_random_engine generator;
155  std::uniform_int_distribution<int> distribution(0,
156  static_cast<int>(myPoints.size()) - 1 );
157  double aspect;
158 
159  for(auto t = 0u; t < myNbTrials ; ++t)
160  {
161  unsigned int i,j,k;
162 
163  //We pick 3 distinct point indices.
164  i = distribution(generator);
165  j = distribution(generator);
166  while ( (j = distribution(generator)) == i);
167  while (( (k = distribution(generator)) == i) || (k == j) );
168 
169  RealPoint vector = getNormal(i,j,k,aspect);
170  if ((vector.norm() > 0.00001) && (aspect > myAspectRatio))
171  {
172  //we have an admissible triangle, we push both normal vectors
173  for(auto acc = 0u; acc < myNbAccumulators; ++acc)
174  {
175  RealPoint shifted = myRotations[acc]*vector;
176  myAccumulators[acc].addDirection( shifted );
177  myAccumulators[acc].addDirection( -shifted );
178  }
179  }
180  }
181  //We return the max bin orientation summing up all accumulators vote
182  typename SphericalAccumulator<RealPoint>::Size posPhi,posTheta;
183  RealPoint vote;
184  for(auto acc = 0u; acc < myNbAccumulators; ++acc)
185  {
186  myAccumulators[acc].maxCountBin(posPhi, posTheta);
187  RealPoint dir = myInverseRotations[acc]*myAccumulators[acc].representativeDirection(posPhi, posTheta).getNormalized() ;
188 
189  //We only consider z-oriented normals (since we pushed vector and -vector)
190  if ( dir.dot(RealPoint(0,0,1)) > 0.0 )
191  vote += dir;
192  else
193  vote += -dir;
194  }
195  return vote.getNormalized();
196  }
size_t Size
Type to represent bin indexes.
std::vector< RealPoint > myPoints
vector of embedded surfels
Quantity getNormal(const unsigned int i, const unsigned int j, const unsigned int k, double &aspect) const
PointVector< 3, double > RealPoint

References DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::getNormal(), DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myAccumulators, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myAspectRatio, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myInverseRotations, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myNbAccumulators, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myNbTrials, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myPoints, and DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myRotations.

◆ getNormal()

template<typename TSurfel , typename TEmbedder >
Quantity DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::getNormal ( const unsigned int  i,
const unsigned int  j,
const unsigned int  k,
double &  aspect 
) const
inlineprivate

Computes the (unnormalized) normal vector of a triangle defined by triangle (i,j,k). The variable aspect returns the aspect ratio of the triangle.

Parameters
[in]ia first vertex index.
[in]ja second vertex index.
[in]ka third vertex index.
[out]aspectaspect ratio of the triangle.
Returns
a random rotation matrix.

Definition at line 268 of file SphericalHoughNormalVectorEstimator.h.

272  {
273  ASSERT( i < myPoints.size());
274  ASSERT( j < myPoints.size());
275  ASSERT( k < myPoints.size());
276 
277  const RealPoint v = myPoints[i] - myPoints[j];
278  const RealPoint u = myPoints[i] - myPoints[k];
279  const RealPoint w = myPoints[j] - myPoints[k];
280 
281  //aspect ratio
282  const double a = u.norm() , b = v.norm();
283  const double c = w.norm();
284  const double s = (a+b+c)/2.0;
285  double denom = (8.0*(s-a)*(s-b)*(s-c));
286  if ( std::abs( denom ) <= std::abs( denom ) * 1e-6 )
287  aspect = 0.;
288  else
289  aspect = a*b*c / denom;
290 
291  return v.crossProduct(u);
292  }

References DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myPoints.

Referenced by DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::eval().

◆ pushSurfel()

template<typename TSurfel , typename TEmbedder >
void DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::pushSurfel ( const Surfel aSurf,
const double  aDistance 
)
inline

Add the geometrical embedding of a surfel to the point list and update the normal spherical hough voting.

Parameters
aSurfa surfel to add
aDistancedistance of aSurf to the neighborhood boundary (NOT USED HERE)

Definition at line 138 of file SphericalHoughNormalVectorEstimator.h.

140  {
141  BOOST_VERIFY(aDistance == aDistance);
142  RealPoint p = myH * ( myEmbedder->operator()(aSurf) );
143  myPoints.push_back(p);
144  }

References DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myEmbedder, DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myH, and DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myPoints.

◆ randomRotation()

template<typename TSurfel , typename TEmbedder >
Matrix DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::randomRotation ( ) const
inlineprivate
Returns
a random rotation matrix.

Definition at line 215 of file SphericalHoughNormalVectorEstimator.h.

216  {
217  const double theta = (rand()+0.f)/RAND_MAX * 2* M_PI;
218  const double phi = (rand()+0.f)/RAND_MAX * 2* M_PI;
219  const double psi = (rand()+0.f)/RAND_MAX * 2* M_PI;
220  Matrix Rt;
221  Rt.setComponent(0,0,1);
222  Rt.setComponent(1,0,0);
223  Rt.setComponent(2,0,0);
224  Rt.setComponent(0,1,0);
225  Rt.setComponent(1,1,cos(theta));
226  Rt.setComponent(2,1,-sin(theta));
227  Rt.setComponent(0,2,0);
228  Rt.setComponent(1,2,sin(theta));
229  Rt.setComponent(2,2,cos(theta));
230 
231  Matrix Rph;
232  Rph.setComponent(0,0,cos(phi));
233  Rph.setComponent(1,0,0);
234  Rph.setComponent(2,0,sin(phi));
235  Rph.setComponent(0,1,0);
236  Rph.setComponent(1,1,1);
237  Rph.setComponent(2,1,0);
238  Rph.setComponent(0,2,-sin(phi));
239  Rph.setComponent(1,2,0);
240  Rph.setComponent(2,2,cos(phi));
241 
242  Matrix Rps;
243  Rps.setComponent(0,0,cos(psi));
244  Rps.setComponent(1,0,-sin(psi));
245  Rps.setComponent(2,0,0);
246  Rps.setComponent(0,1,sin(psi));
247  Rps.setComponent(1,1,cos(psi));
248  Rps.setComponent(2,1,0);
249  Rps.setComponent(0,2,0);
250  Rps.setComponent(1,2,0);
251  Rps.setComponent(2,2,1);
252 
253  return Rt*Rph*Rps;
254  }
void setComponent(const DGtal::Dimension i, const DGtal::Dimension j, const Component &aValue)

References DGtal::SimpleMatrix< TComponent, TM, TN >::setComponent().

Referenced by DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::SphericalHoughNormalVectorEstimator().

◆ reset()

template<typename TSurfel , typename TEmbedder >
void DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::reset ( )
inline

Field Documentation

◆ myAccumulators

◆ myAspectRatio

template<typename TSurfel , typename TEmbedder >
const double DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myAspectRatio
private

Minimal aspect ratio (norm of the cross-product) to consider a given triangle.

Definition at line 301 of file SphericalHoughNormalVectorEstimator.h.

Referenced by DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::eval().

◆ myEmbedder

template<typename TSurfel , typename TEmbedder >
const SCellEmbedder* DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myEmbedder
private

◆ myH

template<typename TSurfel , typename TEmbedder >
const double DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myH
private

◆ myInverseRotations

◆ myNbAccumulators

◆ myNbTrials

template<typename TSurfel , typename TEmbedder >
const unsigned int DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::myNbTrials
private

Number of trials in the neignborhood.

Definition at line 304 of file SphericalHoughNormalVectorEstimator.h.

Referenced by DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::eval().

◆ myPoints

◆ myRotations

◆ mySize

template<typename TSurfel , typename TEmbedder >
const unsigned int DGtal::functors::SphericalHoughNormalVectorEstimator< TSurfel, TEmbedder >::mySize
private

The documentation for this class was generated from the following file: