DGtal 1.3.0
Loading...
Searching...
No Matches
Public Types | Private Member Functions | Private Attributes
DGtal::ShroudsRegularization< TDigitalSurfaceContainer > Class Template Reference

Aim: Implements the Shrouds Regularization algorithm of Nielson et al [90]. More...

#include <DGtal/geometry/surfaces/ShroudsRegularization.h>

Public Types

enum class  Regularization { AREA , SNAKE , SQUARED_CURVATURE }
 The enum class specifying the possible shrouds regularization. More...
 
typedef TDigitalSurfaceContainer Container
 
typedef ShroudsRegularization< ContainerSelf
 
typedef Container::KSpace KSpace
 
typedef KSpace::Space Space
 
typedef Space::RealVector RealVector
 
typedef Space::RealPoint RealPoint
 
typedef RealVector::Component Scalar
 
typedef IndexedDigitalSurface< ContainerIdxDigitalSurface
 
typedef IdxDigitalSurface::Vertex IdxVertex
 
typedef IdxDigitalSurface::Surfel IdxSurfel
 
typedef IdxVertex Vertex
 
typedef std::vector< IdxSurfelIdxSurfelRange
 
typedef std::vector< ScalarScalars
 
typedef std::vector< RealVectorRealVectors
 
typedef std::vector< RealPointRealPoints
 

Public Member Functions

Standard services (construction, initialization)
 ShroudsRegularization ()
 Default constructor. The object is not valid. More...
 
 ShroudsRegularization (CountedPtr< IdxDigitalSurface > surface)
 
void init ()
 
void setParams (double eps, double alpha=1.0, double beta=1.0)
 
std::tuple< double, double, double > getParams () const
 
Accessor services
RealPoint position (const Vertex v, const double t) const
 
RealPoint position (const Vertex v) const
 
RealPoints positions () const
 
Dimension orthDir (const Vertex v) const
 
std::pair< Vertex, Dimensionnext (const std::pair< Vertex, Dimension > &v_i) const
 
std::pair< Vertex, Dimensionprev (const std::pair< Vertex, Dimension > &v_i) const
 
Geometric services
void parameterize ()
 Computes the distances between the vertices along slices. More...
 
Scalar c1 (const std::pair< Vertex, Dimension > &v_i) const
 
std::tuple< Scalar, Scalar, Scalarc2_all (const std::pair< Vertex, Dimension > &v_i) const
 
Regularization services
std::pair< double, double > regularize (const Regularization reg=Regularization::SQUARED_CURVATURE, const double randomization=0.0, const double max_loo=0.0001, const int maxNb=100)
 
double energy (const Regularization reg=Regularization::SQUARED_CURVATURE)
 
double energySquaredCurvature ()
 
double energySnake ()
 
double energyArea ()
 
std::pair< double, double > oneStepAreaMinimization (const double randomization=0.0)
 
std::pair< double, double > oneStepSnakeMinimization (const double alpha=1.0, const double beta=1.0, const double randomization=0.0)
 
std::pair< double, double > oneStepSquaredCurvatureMinimization (const double randomization=0.0)
 
void enforceBounds ()
 Forces t to stay in ]0,1[. More...
 

Protected Member Functions

Internal methods
void precomputeTopology ()
 

Private Member Functions

 BOOST_CONCEPT_ASSERT ((concepts::CDigitalSurfaceContainer< TDigitalSurfaceContainer >))
 

Private Attributes

CountedPtr< IdxDigitalSurfacemyPtrIdxSurface
 the indexed digital surface (internal surface representation). More...
 
const KSpacemyPtrK
 A const pointer to the cellular space in which lives the digital surface. More...
 
Scalar myEpsilon
 
Scalar myAlpha
 The alpha parameter for Snake first order regularization (~ area) More...
 
Scalar myBeta
 The beta parameter for Snake second order regularization (~ curvature) More...
 
Vertex myInvalid
 the index of the invalid vertex. More...
 
Scalars myT
 
RealPoints myInsV
 
RealPoints myOutV
 
std::vector< DimensionmyOrthDir
 the direction axis of each dual edge. More...
 
std::vector< VertexmyNext [3]
 for each vertex, its successor on the slice of given axis direction. More...
 
std::vector< VertexmyPrev [3]
 for each vertex, its predessor on the slice of given axis direction. More...
 
Scalars myNextD [3]
 
Scalars myPrevD [3]
 

Detailed Description

template<typename TDigitalSurfaceContainer>
class DGtal::ShroudsRegularization< TDigitalSurfaceContainer >

Aim: Implements the Shrouds Regularization algorithm of Nielson et al [90].

Description of template class 'ShroudsRegularization'

Represents a decomposition of a closed digital surface into three stacks of slices, one per dimension. Gives methods to optimize the positions of vertices and smooth the resulting surface, generally according to curvature-related regularizers.

See also
Digital surface regularization using the Shrouds algorithm
Note
This method is limited to closed digital surfaces.
Template Parameters
TDigitalSurfaceContainerany digital surface container (a model concepts::CDigitalSurfaceContainer), for instance a SetOfSurfels.
See also
testShroudsRegularization.cpp

Definition at line 73 of file ShroudsRegularization.h.

Member Typedef Documentation

◆ Container

template<typename TDigitalSurfaceContainer >
typedef TDigitalSurfaceContainer DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::Container

Definition at line 79 of file ShroudsRegularization.h.

◆ IdxDigitalSurface

template<typename TDigitalSurfaceContainer >
typedef IndexedDigitalSurface< Container > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::IdxDigitalSurface

Definition at line 86 of file ShroudsRegularization.h.

◆ IdxSurfel

template<typename TDigitalSurfaceContainer >
typedef IdxDigitalSurface::Surfel DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::IdxSurfel

Definition at line 88 of file ShroudsRegularization.h.

◆ IdxSurfelRange

template<typename TDigitalSurfaceContainer >
typedef std::vector< IdxSurfel > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::IdxSurfelRange

Definition at line 90 of file ShroudsRegularization.h.

◆ IdxVertex

template<typename TDigitalSurfaceContainer >
typedef IdxDigitalSurface::Vertex DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::IdxVertex

Definition at line 87 of file ShroudsRegularization.h.

◆ KSpace

template<typename TDigitalSurfaceContainer >
typedef Container::KSpace DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::KSpace

Definition at line 81 of file ShroudsRegularization.h.

◆ RealPoint

template<typename TDigitalSurfaceContainer >
typedef Space::RealPoint DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::RealPoint

Definition at line 84 of file ShroudsRegularization.h.

◆ RealPoints

template<typename TDigitalSurfaceContainer >
typedef std::vector< RealPoint > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::RealPoints

Definition at line 93 of file ShroudsRegularization.h.

◆ RealVector

template<typename TDigitalSurfaceContainer >
typedef Space::RealVector DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::RealVector

Definition at line 83 of file ShroudsRegularization.h.

◆ RealVectors

template<typename TDigitalSurfaceContainer >
typedef std::vector< RealVector > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::RealVectors

Definition at line 92 of file ShroudsRegularization.h.

◆ Scalar

template<typename TDigitalSurfaceContainer >
typedef RealVector::Component DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::Scalar

Definition at line 85 of file ShroudsRegularization.h.

◆ Scalars

template<typename TDigitalSurfaceContainer >
typedef std::vector< Scalar > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::Scalars

Definition at line 91 of file ShroudsRegularization.h.

◆ Self

template<typename TDigitalSurfaceContainer >
typedef ShroudsRegularization< Container > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::Self

Definition at line 80 of file ShroudsRegularization.h.

◆ Space

template<typename TDigitalSurfaceContainer >
typedef KSpace::Space DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::Space

Definition at line 82 of file ShroudsRegularization.h.

◆ Vertex

template<typename TDigitalSurfaceContainer >
typedef IdxVertex DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::Vertex

Definition at line 89 of file ShroudsRegularization.h.

Member Enumeration Documentation

◆ Regularization

template<typename TDigitalSurfaceContainer >
enum class DGtal::ShroudsRegularization::Regularization
strong

The enum class specifying the possible shrouds regularization.

Enumerator
AREA 
SNAKE 
SQUARED_CURVATURE 

Definition at line 96 of file ShroudsRegularization.h.

Constructor & Destructor Documentation

◆ ShroudsRegularization() [1/2]

template<typename TDigitalSurfaceContainer >
DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::ShroudsRegularization ( )
inline

Default constructor. The object is not valid.

Definition at line 104 of file ShroudsRegularization.h.

105 : myPtrIdxSurface( nullptr ), myPtrK( nullptr )
106 {}
CountedPtr< IdxDigitalSurface > myPtrIdxSurface
the indexed digital surface (internal surface representation).
const KSpace * myPtrK
A const pointer to the cellular space in which lives the digital surface.

◆ ShroudsRegularization() [2/2]

template<typename TDigitalSurfaceContainer >
DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::ShroudsRegularization ( CountedPtr< IdxDigitalSurface surface)
inline

Constructor from (closed) surface. Also calls methods precomputeTopology and init.

Parameters
surfacea counted pointer on an indexed digital surface.
Note
Complexity is linear in the number of surfels of surface.

Definition at line 113 of file ShroudsRegularization.h.

114 : myPtrIdxSurface( surface ),
115 myPtrK( &surface->container().space() ),
116 myEpsilon( 0.0001 ), myAlpha( 1.0 ), myBeta( 1.0 )
117 {
119 init();
120 }
Scalar myAlpha
The alpha parameter for Snake first order regularization (~ area)
Scalar myBeta
The beta parameter for Snake second order regularization (~ curvature)

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::init(), and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::precomputeTopology().

Member Function Documentation

◆ BOOST_CONCEPT_ASSERT()

template<typename TDigitalSurfaceContainer >
DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::BOOST_CONCEPT_ASSERT ( (concepts::CDigitalSurfaceContainer< TDigitalSurfaceContainer >)  )
private

◆ c1()

template<typename TDigitalSurfaceContainer >
Scalar DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::c1 ( const std::pair< Vertex, Dimension > &  v_i) const
inline
Parameters
v_ia pair (vertex,tangent direction)
Returns
the coefficient for centered first-order finite difference.
Note
We have y'i ~= c1 * (y{i+1} - y_{i-1} )

Definition at line 261 of file ShroudsRegularization.h.

262 {
263 const Scalar din = myNextD[ v_i.second ][ v_i.first ];
264 const Scalar dip = myPrevD[ v_i.second ][ v_i.first ];
265 return 1.0 / (din + dip);
266 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNextD, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrevD.

◆ c2_all()

template<typename TDigitalSurfaceContainer >
std::tuple< Scalar, Scalar, Scalar > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::c2_all ( const std::pair< Vertex, Dimension > &  v_i) const
inline
Parameters
v_ia pair (vertex,tangent direction)
Returns
the coefficients for centered second-order finite difference.
Note
We have y''i ~= c2<0> * y{i+1} - c2<1> * y_{i} + c2<2> * y_{i-1}

Definition at line 271 of file ShroudsRegularization.h.

272 {
273 const Scalar din = myNextD[ v_i.second ][ v_i.first ];
274 const Scalar dip = myPrevD[ v_i.second ][ v_i.first ];
275 return std::make_tuple( 2.0 / ( din * ( din + dip ) ),
276 2.0 / ( din * dip ),
277 2.0 / ( dip * ( din + dip ) ) );
278 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNextD, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrevD.

◆ energy()

template<typename TDigitalSurfaceContainer >
double DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energy ( const Regularization  reg = Regularization::SQUARED_CURVATURE)
inline
Parameters
regyour choice of regularization in AREA (not nice in 3D), SNAKE (not so bad), SQUARED_CURVATURE (works best)
Returns
the current energy of the shrouds, according to the chosen regularization process.

Definition at line 335 of file ShroudsRegularization.h.

336 {
338 return energySquaredCurvature();
339 else if ( reg == Regularization::SNAKE )
340 return energySnake();
341 else
342 return energyArea();
343 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energyArea(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energySnake(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energySquaredCurvature(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::SNAKE, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::SQUARED_CURVATURE.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::regularize(), and TEST_CASE().

◆ energyArea()

template<typename TDigitalSurfaceContainer >
double DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energyArea ( )
Returns
the current energy associated to the area regularization process.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energy().

◆ energySnake()

template<typename TDigitalSurfaceContainer >
double DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energySnake ( )
Returns
the current energy associated to the snake regularization process.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energy().

◆ energySquaredCurvature()

template<typename TDigitalSurfaceContainer >
double DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energySquaredCurvature ( )
Returns
the current energy associated to the squared curvature regularization process.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energy().

◆ enforceBounds()

template<typename TDigitalSurfaceContainer >
void DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::enforceBounds ( )

Forces t to stay in ]0,1[.

◆ getParams()

template<typename TDigitalSurfaceContainer >
std::tuple< double, double, double > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::getParams ( ) const
inline

Retrieves the parameters that affect the output regularized shape.

Returns
a tuple (eps, alpha, beta) where eps is the bounds for varying the positions of vertices in ]0,1[, alpha is the parameter for Snake first order regularization (~ area), beta is the parameter for Snake second order regularization (~ curvature).

Definition at line 167 of file ShroudsRegularization.h.

168 {
169 return std::make_tuple( myEpsilon, myAlpha, myBeta );
170 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myAlpha, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myBeta, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myEpsilon.

◆ init()

template<typename TDigitalSurfaceContainer >
void DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::init ( )
inline

Prepares the shroud for optimization. Must be called before any call to regularize, oneStepAreaMinimization, oneStepSnakeMinimization, oneStepSquaredCurvatureMinimization.

Note
Complexity is linear in the number of surfels of surface.

Definition at line 128 of file ShroudsRegularization.h.

129 {
130 const auto embedder = CanonicSCellEmbedder<KSpace>( *myPtrK );
131 const auto nbV = myPtrIdxSurface->nbVertices();
132 myT.resize( nbV );
133 myInsV.resize( nbV );
134 myOutV.resize( nbV );
135 for ( Dimension i = 0; i < 3; ++i )
136 {
137 myPrevD[ i ].resize( nbV );
138 myNextD[ i ].resize( nbV );
139 }
140 for ( Vertex v = 0; v < myT.size(); ++v )
141 {
142 const auto s = myPtrIdxSurface->surfel( v );
143 const auto k = myPtrK->sOrthDir( s );
144 myInsV[ v ] = embedder( myPtrK->sDirectIncident( s, k ) );
145 myOutV[ v ] = embedder( myPtrK->sIndirectIncident( s, k ) );
146 }
147 }
DGtal::uint32_t Dimension
Definition: Common.h:137
TriMesh::Vertex Vertex

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myInsV, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNextD, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOutV, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrevD, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPtrIdxSurface, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPtrK, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myT.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::ShroudsRegularization(), and TEST_CASE().

◆ next()

template<typename TDigitalSurfaceContainer >
std::pair< Vertex, Dimension > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::next ( const std::pair< Vertex, Dimension > &  v_i) const
inline

Useful to navigate tangentially along a slice.

Parameters
v_ia pair (vertex,tangent direction)
Returns
the next vertex and associated tangent direction along the slice.

Definition at line 216 of file ShroudsRegularization.h.

217 {
218 const Vertex vn = myNext[ v_i.second ][ v_i.first ];
219 const Dimension in = myOrthDir[ vn ] == v_i.second
220 ? myOrthDir[ v_i.first ] : v_i.second;
221 return std::make_pair( vn, in );
222 }
std::vector< Vertex > myNext[3]
for each vertex, its successor on the slice of given axis direction.
std::vector< Dimension > myOrthDir
the direction axis of each dual edge.

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNext, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOrthDir.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::precomputeTopology().

◆ oneStepAreaMinimization()

template<typename TDigitalSurfaceContainer >
std::pair< double, double > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::oneStepAreaMinimization ( const double  randomization = 0.0)

Smooths the shape according to the minimization of area.

Parameters
randomizationif greater than 0.0 add some perturbation to the solution. May be used for quitting a local minima.
Returns
the pair of \( l_\infty \) and \( l_2 \) norms of vertex displacements.
Note
Not very nice.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::regularize().

◆ oneStepSnakeMinimization()

template<typename TDigitalSurfaceContainer >
std::pair< double, double > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::oneStepSnakeMinimization ( const double  alpha = 1.0,
const double  beta = 1.0,
const double  randomization = 0.0 
)

Smooths the shape according to the minimization of elastic + thin plate energy (like snakes).

\[ E^{snk}(C) = \int_C \alpha (x'(s)^2 + y'(s)^2) + \beta (x''(s)^2 + y''(s)^2) ds, \]

for \( C=(x(s),y(s)) \) and boundary constraints.

Parameters
alphaparameter for first order regularization (~ area)
betaparameter for second order regularization (~ curvature)
randomizationif greater than 0.0 add some perturbation to the solution. May be used for quitting a local minima.
Returns
the pair of \( l_\infty \) and \( l_2 \) norms of vertex displacements.

Euler-Lagrange equations leads to necessary conditions: for all s, \( x''''(s) = 0 \) and \( y''''(s) = 0 \).

A vertex \( v(s) \) has a position \( X \). Next and previous vertices are denoted by \( Xn, Xnn, Xp, Xpp \). If i the tangent direction and k the orthogonal direction at \( X \), then \( x(s)=X[i], y(s)=X[k] \).

At a vertex \( v(s) \), the position is parameterized by t. Varying t only modifies the vertical position \( X[k] \), i.e. y. So \( y(t) = (1-t)*I[v][k] + t*O[v][k] \) (for \( I[v] \) and \( O[v] \) inside and outside voxel positions).

\( x' ~= ( Xn[i] - Xp[i] ) / (d(X,Xn)+d(X,Xp) \) does not depend on t. \( y' ~= ( Xn[k] - Xp[k] ) / (d(X,Xn)+d(X,Xp) \) does not depend on t. \( x'' ~= ( c_0 * Xn[i] - c_1 * X[i] + c_2 * Xp[i] ) \) does not depend on t.

for \( c_0, c_1, c_2 \) constants depending on \( d(X,Xn) \) and \( d(X,Xp) \). \( y'' ~= ( c_0 * Xn[k] - c_1 * X[k] + c_2 * Xp[k] ) \) solely depend on t.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::regularize().

◆ oneStepSquaredCurvatureMinimization()

template<typename TDigitalSurfaceContainer >
std::pair< double, double > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::oneStepSquaredCurvatureMinimization ( const double  randomization = 0.0)

Smooths the shape according to the minimization of squared curvature.

\[ E^{\kappa^2}(C) = \int_C (x'(s) y''(s) + x''(s) y'(s))^2 / (x'(s)^2 + y'(s)^2)^3 ds, \]

for \( C=(x(s),y(s)) \) and boundary constraints.

Parameters
randomizationif greater than 0.0 add some perturbation to the solution. May be used for quitting a local minima.
Returns
the pair of \( l_\infty \) and \( l_2 \) norms of vertex displacements.

Euler-Lagrange equations leads to rather complex necessary conditions. At each position (we omit parameter s for making things more readable):

\[ 24*x'^3*x''^3*y' + x'''*( - 13*x'^4*x''*y' - 14*x'^2*x''*y'^3 - x''*y'^5 ) + y'''*( 8*x'^5*x'' + 4*x'^3*x''*y'^2 - 4*x'*x''*y'^4 ) + x''''*( x'^5*y' + 2*x'^3*y'^3 + x'*y'^5 ) + y''*( 12*x'^4*y'*y''' + 12*x'^2*y'^3*y''' - 24*x'^4*x''^2 + 5*x'^5*x''' + 51*x'^2*x''^2*y'^2 - 2*x'^3*x'''*y'^2 + 3*x''^2*y'^4 - 7*x'*x'''*y'^4 ) + y''^2*( -54*x'^3*x''*y' + 18*x'*x''*y'^3 ) + y''^3*( 3*x'^4 - 21*x'^2*y'^2 ) + y''''*( - x'^6 - 2*x'^4*y'^2 - x'^2*y'^4 ) == 0 \]

A vertex \( v(s) \) has a position \( X \). Next and previous vertices are denoted by \( Xn, Xnn, Xp, Xpp \). If i the tangent direction and k the orthogonal direction at \( X \), then \( x(s)=X[i], y(s)=X[k] \).

At a vertex \( v(s) \), the position is parameterized by t. Varying t only modifies the vertical position \( X[k] \), i.e. y. So \( y(t) = (1-t)*I[v][k] + t*O[v][k] \) (for \( I[v] \) and \( O[v] \) inside and outside voxel positions).

The Euler-Lagrange are factorized as above in order to isolate at best the parameter t of vertex \( v(s) \).

  • The first four lines do not depend on t (except \( y''' \) but it is neglected and the old value is used).
  • \( y''^2 \) and \( y''^3 \) are linarized as \( y'' * y''_{old} \) and \( y'' * y''_{old}^2 \)
  • \( y'''' ~= ( c_0 * Yn - c_1 * y'' + c_2 * Yp ) \)
  • all \( y'' \) are sumed up, while the rest gives a constant.

    \( x' ~= ( Xn[i] - Xp[i] ) / (d(X,Xn)+d(X,Xp) \) does not depend on t. \( y' ~= ( Xn[k] - Xp[k] ) / (d(X,Xn)+d(X,Xp) \) does not depend on t. \( x'' ~= ( c_0 * Xn[i] - c_1 * X[i] + c_2 * Xp[i] ) \) does not depend on t.

for \( c_0, c_1, c_2 \) constants depending on \( d(X,Xn) \) and \( d(X,Xp) \). \( y'' ~= ( c_0 * Xn[k] - c_1 * X[k] + c_2 * Xp[k] ) \) solely depend on t. where \( X[k] = y(t) = (1-t)*I[v][k] + t*O[v][k]. \) And so on.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::regularize().

◆ orthDir()

template<typename TDigitalSurfaceContainer >
Dimension DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::orthDir ( const Vertex  v) const
inline

Useful to find the two tangent directions to a given vertex (i.e. (orthDir(v)+1)%3 and (orthDir(v)+2)%3).

Returns
the orthogonal direction to vertex v.

Definition at line 208 of file ShroudsRegularization.h.

209 {
210 return myOrthDir[ v ];
211 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOrthDir.

◆ parameterize()

template<typename TDigitalSurfaceContainer >
void DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::parameterize ( )
inline

Computes the distances between the vertices along slices.

Definition at line 243 of file ShroudsRegularization.h.

244 {
245 Scalar d = 0.0;
246 // Vertex n = 0;
247 for ( Dimension i = 0; i < 3; ++i )
248 for ( Vertex v = 0; v < myT.size(); ++v )
249 {
250 if ( myNext[ i ][ v ] == myInvalid ) continue; // not a valid slice
251 myNextD[ i ][ v ] = ( position( myNext[ i ][ v ] ) - position( v ) ).norm();
252 myPrevD[ i ][ v ] = ( position( myPrev[ i ][ v ] ) - position( v ) ).norm();
253 d += myNextD[ i ][ v ];
254 // n += 1;
255 }
256 }
RealPoint position(const Vertex v, const double t) const
Vertex myInvalid
the index of the invalid vertex.
std::vector< Vertex > myPrev[3]
for each vertex, its predessor on the slice of given axis direction.

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myInvalid, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNext, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNextD, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrev, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrevD, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myT, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::position().

◆ position() [1/2]

template<typename TDigitalSurfaceContainer >
RealPoint DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::position ( const Vertex  v) const
inline
Parameters
vany valid vertex.
Returns
its position.

Definition at line 189 of file ShroudsRegularization.h.

190 {
191 const auto t = myT[ v ];
192 return (1-t) * myInsV[ v ] + t * myOutV[ v ];
193 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myInsV, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOutV, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myT.

◆ position() [2/2]

template<typename TDigitalSurfaceContainer >
RealPoint DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::position ( const Vertex  v,
const double  t 
) const
inline
Parameters
vany valid vertex.
tsome adjustement parameter
Returns
the position of vertex v for this parameter t.

Definition at line 182 of file ShroudsRegularization.h.

183 {
184 return (1-t) * myInsV[ v ] + t * myOutV[ v ];
185 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myInsV, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOutV.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::parameterize(), and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::positions().

◆ positions()

template<typename TDigitalSurfaceContainer >
RealPoints DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::positions ( ) const
inline
Returns
the vector of vertex positions

Definition at line 196 of file ShroudsRegularization.h.

197 {
198 RealPoints result( myT.size() );
199 for ( Vertex v = 0; v < myT.size(); ++v )
200 result[ v ] = position( v );
201 return result;
202 }
std::vector< RealPoint > RealPoints

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myT, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::position().

Referenced by TEST_CASE().

◆ precomputeTopology()

template<typename TDigitalSurfaceContainer >
void DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::precomputeTopology ( )
inlineprotected

This method precomputes the neighbors of each vertex along each crossing curves. Must be called at shroud initialization.

Definition at line 476 of file ShroudsRegularization.h.

477 {
478 typedef typename Container::Tracker Tracker;
479 typedef DigitalSurface2DSlice< Tracker > Slice;
480 myT = Scalars( myPtrIdxSurface->nbVertices(), 0.5 );
481 myInvalid = myT.size();
482 myOrthDir.resize( myT.size() );
483 for ( Dimension i = 0; i < 3; ++i )
484 {
485 myNext[ i ] = std::vector<Vertex>( myT.size(), myInvalid );
486 myPrev[ i ] = std::vector<Vertex>( myT.size(), myInvalid );
487 }
488 // for each vertex, extracts its two slices.
489 for ( Vertex v = 0; v < myT.size(); ++v )
490 {
491 auto surf = myPtrIdxSurface->surfel( v );
492 Dimension k = myPtrK->sOrthDir( surf );
493 myOrthDir[ v ] = k;
494 for ( Dimension i = 0; i < 3; ++i )
495 {
496 if ( k == i ) continue; // not a valid slice
497 if ( myNext[ i ][ v ] != myInvalid ) continue; // already computed
498 Tracker* tracker = myPtrIdxSurface->container().newTracker( surf );
499 Slice slice( tracker, i );
500 if ( ! slice.isClosed() ) {
501 trace.error() << "[ShroudsRegularization::precomputeTopology]"
502 << " Shrouds works solely on closed surfaces."
503 << std::endl;
504 return;
505 }
506 auto start = slice.cstart();
507 auto next = start;
508 auto prev = next++;
509 Vertex vp = v;
510 do
511 {
512 auto sp = *prev;
513 auto sn = *next;
514 Dimension in = myPtrK->sOrthDir( sn ) == k ? i : k;
515 Dimension ip = myPtrK->sOrthDir( sp ) == k ? i : k;
516 Vertex vn = myPtrIdxSurface->getVertex( sn );
517 myNext[ ip ][ vp ] = vn;
518 myPrev[ in ][ vn ] = vp;
519 prev = next++;
520 vp = vn;
521 }
522 while ( prev != start );
523 delete tracker;
524 }
525 }
526 }
std::pair< Vertex, Dimension > prev(const std::pair< Vertex, Dimension > &v_i) const
std::pair< Vertex, Dimension > next(const std::pair< Vertex, Dimension > &v_i) const
std::ostream & error()
Trace trace
Definition: Common.h:154

References DGtal::Trace::error(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myInvalid, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNext, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOrthDir, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrev, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPtrIdxSurface, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPtrK, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myT, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::next(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::prev(), and DGtal::trace.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::ShroudsRegularization().

◆ prev()

template<typename TDigitalSurfaceContainer >
std::pair< Vertex, Dimension > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::prev ( const std::pair< Vertex, Dimension > &  v_i) const
inline

Useful to navigate tangentially along a slice.

Parameters
v_ia pair (vertex,tangent direction)
Returns
the previous vertex and associated tangent direction along the slice.

Definition at line 227 of file ShroudsRegularization.h.

228 {
229 const Vertex vp = myPrev[ v_i.second ][ v_i.first ];
230 const Dimension ip = myOrthDir[ vp ] == v_i.second
231 ? myOrthDir[ v_i.first ] : v_i.second;
232 return std::make_pair( vp, ip );
233 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOrthDir, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrev.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::precomputeTopology().

◆ regularize()

template<typename TDigitalSurfaceContainer >
std::pair< double, double > DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::regularize ( const Regularization  reg = Regularization::SQUARED_CURVATURE,
const double  randomization = 0.0,
const double  max_loo = 0.0001,
const int  maxNb = 100 
)
inline

Generic method for shrouds regularization.

Parameters
regyour choice of regularization in AREA (not nice in 3D), SNAKE (not so bad), SQUARED_CURVATURE (works best)
randomizationif greater than 0.0 add some perturbation to the solution. May be used for quitting a local minima.
max_loothe maximum \( l_\infty \)-norm of one step of regularization before stopping the regularization.
maxNbthe maximum number of optimization steps.
See also
oneStepAreaMinimization
oneStepSnakeMinimization
oneStepSquaredCurvatureMinimization

Definition at line 304 of file ShroudsRegularization.h.

308 {
309 double loo = 0.0;
310 double l2 = 0.0;
311 int nb = 0;
312 double r = 0.5;
313 do {
314 std::tie( loo, l2 ) =
317 : reg == Regularization::SNAKE
320 if ( nb % 50 == 0 )
321 trace.info() << "[Shrouds iteration " << nb
322 << " E=" << energy( reg )
323 << "] dx <= " << loo << " l2=" << l2 << std::endl;
324 r *= 0.9;
325 nb += 1;
326 } while ( loo > max_loo && nb < maxNb );
327 return std::make_pair( loo, l2 );
328 }
Regularization
The enum class specifying the possible shrouds regularization.
std::pair< double, double > oneStepAreaMinimization(const double randomization=0.0)
double energy(const Regularization reg=Regularization::SQUARED_CURVATURE)
std::pair< double, double > oneStepSquaredCurvatureMinimization(const double randomization=0.0)
std::pair< double, double > oneStepSnakeMinimization(const double alpha=1.0, const double beta=1.0, const double randomization=0.0)
std::ostream & info()

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::energy(), DGtal::Trace::info(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myAlpha, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myBeta, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::oneStepAreaMinimization(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::oneStepSnakeMinimization(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::oneStepSquaredCurvatureMinimization(), DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::SNAKE, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::SQUARED_CURVATURE, and DGtal::trace.

Referenced by TEST_CASE().

◆ setParams()

template<typename TDigitalSurfaceContainer >
void DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::setParams ( double  eps,
double  alpha = 1.0,
double  beta = 1.0 
)
inline

Sets some parameters that affect the output regularized shape.

Parameters
epsthe bounds for varying the positions of vertices in ]0,1[
alphaparameter for Snake first order regularization (~ area)
betaparameter for Snake second order regularization (~ curvature)

Definition at line 154 of file ShroudsRegularization.h.

155 {
156 myEpsilon = eps;
157 myAlpha = alpha;
158 myBeta = beta;
159 }

References DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myAlpha, DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myBeta, and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myEpsilon.

Field Documentation

◆ myAlpha

template<typename TDigitalSurfaceContainer >
Scalar DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myAlpha
private

◆ myBeta

template<typename TDigitalSurfaceContainer >
Scalar DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myBeta
private

◆ myEpsilon

template<typename TDigitalSurfaceContainer >
Scalar DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myEpsilon
private

The limiting bounds for the displacement of vertices along their unit dual edge: \( \lbrack \epsilon, 1-\epsilon \rbrack \)

Definition at line 539 of file ShroudsRegularization.h.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::getParams(), and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::setParams().

◆ myInsV

template<typename TDigitalSurfaceContainer >
RealPoints DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myInsV
private

the vector of dual points lying inside (inside extremity of each dual edge).

Definition at line 552 of file ShroudsRegularization.h.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::init(), and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::position().

◆ myInvalid

template<typename TDigitalSurfaceContainer >
Vertex DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myInvalid
private

◆ myNext

template<typename TDigitalSurfaceContainer >
std::vector<Vertex> DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNext[3]
private

◆ myNextD

template<typename TDigitalSurfaceContainer >
Scalars DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myNextD[3]
private

◆ myOrthDir

template<typename TDigitalSurfaceContainer >
std::vector<Dimension> DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOrthDir
private

◆ myOutV

template<typename TDigitalSurfaceContainer >
RealPoints DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myOutV
private

the vector of dual points lying outside (outside extremity of each dual edge).

Definition at line 555 of file ShroudsRegularization.h.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::init(), and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::position().

◆ myPrev

template<typename TDigitalSurfaceContainer >
std::vector<Vertex> DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrev[3]
private

◆ myPrevD

template<typename TDigitalSurfaceContainer >
Scalars DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPrevD[3]
private

◆ myPtrIdxSurface

template<typename TDigitalSurfaceContainer >
CountedPtr<IdxDigitalSurface> DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPtrIdxSurface
private

the indexed digital surface (internal surface representation).

Definition at line 534 of file ShroudsRegularization.h.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::init(), and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::precomputeTopology().

◆ myPtrK

template<typename TDigitalSurfaceContainer >
const KSpace* DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myPtrK
private

A const pointer to the cellular space in which lives the digital surface.

Definition at line 536 of file ShroudsRegularization.h.

Referenced by DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::init(), and DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::precomputeTopology().

◆ myT

template<typename TDigitalSurfaceContainer >
Scalars DGtal::ShroudsRegularization< TDigitalSurfaceContainer >::myT
private

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