DGtal
1.5.beta
|
Aim: This class solves Ambrosio-Tortorelli functional on a two-dimensional digital space (a 2D grid or 2D digital surface) for a piecewise smooth scalar/vector function u represented as one/several 2-form(s) and a discontinuity function v represented as a 0-form. The 2-form(s) u is a regularized approximation of an input vector data g, while v represents the set of discontinuities of u. The norm chosen for u is the \( l_2 \)-norm. More...
#include <DGtal/dec/ATSolver2D.h>
Public Member Functions | |
Standard services | |
ATSolver2D (ConstAlias< Calculus > aCalculus, int aVerbose=0) | |
ATSolver2D ()=delete | |
~ATSolver2D ()=default | |
ATSolver2D (const ATSolver2D &other)=default | |
ATSolver2D (ATSolver2D &&other)=default | |
ATSolver2D & | operator= (const ATSolver2D &other)=default |
ATSolver2D & | operator= (ATSolver2D &&other)=default |
Index | size (const int order) const |
Initialization services | |
template<typename VectorFieldInput , typename SurfelRangeConstIterator > | |
void | initInputVectorFieldU2 (const VectorFieldInput &input, SurfelRangeConstIterator itB, SurfelRangeConstIterator itE, bool normalize=false) |
template<typename ScalarFieldInput , typename SurfelRangeConstIterator > | |
void | initInputScalarFieldU2 (const ScalarFieldInput &input, SurfelRangeConstIterator itB, SurfelRangeConstIterator itE) |
template<typename ScalarVector > | |
Index | initInputVectorFieldU2 (const std::map< Surfel, ScalarVector > &input, bool normalize=false) |
template<typename Scalar > | |
Index | initInputScalarFieldU2 (const std::map< Surfel, Scalar > &input) |
void | setUp (double a, double l) |
void | setUp (double a, double l, const std::map< Surfel, Scalar > &weights) |
template<typename AlphaWeights , typename SurfelRangeConstIterator > | |
void | setUp (double a, double l, const AlphaWeights &weights, SurfelRangeConstIterator itB, SurfelRangeConstIterator itE) |
void | setEpsilon (double e) |
Optimization services | |
bool | solveOneAlternateStep () |
bool | solveForEpsilon (double eps, double n_oo_max=1e-4, unsigned int iter_max=10) |
bool | solveGammaConvergence (double eps1=2.0, double eps2=0.25, double epsr=2.0, bool compute_smallest_epsilon_map=false, double n_oo_max=1e-4, unsigned int iter_max=10) |
void | normalizeU2 () |
std::tuple< double, double, double > | diffV0 () const |
Access services | |
std::tuple< double, double, double > | checkV0 () const |
PrimalForm0 | getV0 () const |
PrimalForm1 | getV1 () const |
PrimalForm2 | getV2 () const |
PrimalForm2 | getU2 (Dimension k) const |
template<typename VectorFieldOutput , typename SurfelRangeConstIterator > | |
void | getOutputVectorFieldU2 (VectorFieldOutput &output, SurfelRangeConstIterator itB, SurfelRangeConstIterator itE) |
template<typename ScalarFieldOutput , typename SurfelRangeConstIterator > | |
void | getOutputScalarFieldU2 (ScalarFieldOutput &output, SurfelRangeConstIterator itB, SurfelRangeConstIterator itE) |
template<typename ScalarFieldOutput , typename CellRangeConstIterator > | |
void | getOutputScalarFieldV0 (ScalarFieldOutput &output, CellRangeConstIterator itB, CellRangeConstIterator itE, CellOutputPolicy policy=CellOutputPolicy::Average) |
void | updateSmallestEpsilonMap (const double threshold=.5) |
Interface services | |
void | selfDisplay (std::ostream &out) const |
bool | isValid () const |
Data Fields | |
Surfel2IndexMap | surfel2idx |
SmallestEpsilonMap | smallest_epsilon_map |
double | alpha |
double | lambda |
double | epsilon |
bool | normalize_u2 |
Indicates whether to normalize U (unit norm) at each iteration or not. More... | |
int | verbose |
Tells the verbose level. More... | |
Static Public Attributes | |
static const Dimension | dimension = KSpace::dimension |
Protected Member Functions | |
Hidden services | |
void | initOperators () |
Initializes the operators. More... | |
Protected Attributes | |
CountedConstPtrOrConstPtr< Calculus > | ptrCalculus |
A smart (or not) pointer to a calculus object. More... | |
PrimalDerivative0 | primal_D0 |
the derivative operator for primal 0-forms More... | |
PrimalDerivative1 | primal_D1 |
the derivative operator for primal 1-forms More... | |
PrimalDerivative0 | M01 |
The primal vertex to edge average operator. More... | |
PrimalDerivative1 | M12 |
The primal edge to face average operator. More... | |
PrimalAntiderivative2 | primal_AD2 |
The antiderivative of primal 2-forms. More... | |
PrimalIdentity2 | alpha_Id2 |
The alpha-weighted identity operator for primal 2-forms (stored for performance) More... | |
PrimalIdentity0 | l_1_over_4e_Id0 |
The 1/(4epsilon)-weighted identity operator for primal 0-forms (stored for performance) More... | |
std::vector< PrimalForm2 > | g2 |
The N-array of input primal 2-forms g. More... | |
std::vector< PrimalForm2 > | alpha_g2 |
The alpha-weighted N-array of input primal 2-forms g. More... | |
std::vector< PrimalForm2 > | u2 |
The N-array of regularized primal 2-forms u. More... | |
PrimalForm0 | v0 |
PrimalForm0 | former_v0 |
The primal 0-form v at the previous iteration. More... | |
PrimalForm0 | l_1_over_4e |
The primal 0-form lambda/(4epsilon) (stored for performance) More... | |
Aim: This class solves Ambrosio-Tortorelli functional on a two-dimensional digital space (a 2D grid or 2D digital surface) for a piecewise smooth scalar/vector function u represented as one/several 2-form(s) and a discontinuity function v represented as a 0-form. The 2-form(s) u is a regularized approximation of an input vector data g, while v represents the set of discontinuities of u. The norm chosen for u is the \( l_2 \)-norm.
Description of template class 'ATSolver2D'
TKSpace | any model of CCellularGridSpaceND, e.g KhalimskySpaceND |
TLinearAlgebra | any back-end for performing linear algebra, default is EigenLinearAlgebraBackend. |
Definition at line 90 of file ATSolver2D.h.
typedef DiscreteExteriorCalculus<2,dimension, LinearAlgebra> DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Calculus |
Definition at line 114 of file ATSolver2D.h.
typedef KSpace::Cell DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Cell |
Definition at line 111 of file ATSolver2D.h.
typedef HyperRectDomain<Space> DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Domain |
Definition at line 113 of file ATSolver2D.h.
typedef Calculus::Index DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Index |
Definition at line 116 of file ATSolver2D.h.
typedef TKSpace DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::KSpace |
Definition at line 95 of file ATSolver2D.h.
typedef TLinearAlgebra DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::LinearAlgebra |
Definition at line 96 of file ATSolver2D.h.
typedef EigenLinearAlgebraBackend::SolverSimplicialLDLT DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::LinearAlgebraSolver |
Definition at line 138 of file ATSolver2D.h.
typedef Calculus::PrimalAntiderivative1 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalAntiderivative1 |
Definition at line 125 of file ATSolver2D.h.
typedef Calculus::PrimalAntiderivative2 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalAntiderivative2 |
Definition at line 126 of file ATSolver2D.h.
typedef Calculus::PrimalDerivative0 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalDerivative0 |
Definition at line 123 of file ATSolver2D.h.
typedef Calculus::PrimalDerivative1 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalDerivative1 |
Definition at line 124 of file ATSolver2D.h.
typedef Calculus::PrimalForm0 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalForm0 |
Definition at line 117 of file ATSolver2D.h.
typedef Calculus::PrimalForm1 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalForm1 |
Definition at line 118 of file ATSolver2D.h.
typedef Calculus::PrimalForm2 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalForm2 |
Definition at line 119 of file ATSolver2D.h.
typedef Calculus::PrimalHodge0 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalHodge0 |
Definition at line 127 of file ATSolver2D.h.
typedef Calculus::PrimalHodge1 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalHodge1 |
Definition at line 128 of file ATSolver2D.h.
typedef Calculus::PrimalHodge2 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalHodge2 |
Definition at line 129 of file ATSolver2D.h.
typedef Calculus::PrimalIdentity0 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalIdentity0 |
Definition at line 120 of file ATSolver2D.h.
typedef Calculus::PrimalIdentity1 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalIdentity1 |
Definition at line 121 of file ATSolver2D.h.
typedef Calculus::PrimalIdentity2 DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::PrimalIdentity2 |
Definition at line 122 of file ATSolver2D.h.
typedef Space::RealVector DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::RealVector |
Definition at line 108 of file ATSolver2D.h.
typedef RealVector::Component DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Scalar |
Definition at line 109 of file ATSolver2D.h.
typedef KSpace::SCell DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::SCell |
Definition at line 110 of file ATSolver2D.h.
typedef ATSolver2D< KSpace, LinearAlgebra > DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Self |
Definition at line 97 of file ATSolver2D.h.
typedef KSpace::template SurfelMap<double>::Type DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::SmallestEpsilonMap |
Definition at line 115 of file ATSolver2D.h.
typedef DiscreteExteriorCalculusSolver<Calculus, LinearAlgebraSolver, 2, PRIMAL, 2, PRIMAL> DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::SolverU2 |
Definition at line 139 of file ATSolver2D.h.
typedef DiscreteExteriorCalculusSolver<Calculus, LinearAlgebraSolver, 0, PRIMAL, 0, PRIMAL> DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::SolverV0 |
Definition at line 140 of file ATSolver2D.h.
typedef KSpace::Space DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Space |
Definition at line 107 of file ATSolver2D.h.
typedef KSpace::Surfel DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Surfel |
Definition at line 112 of file ATSolver2D.h.
typedef KSpace::template SurfelMap<Index>::Type DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::Surfel2IndexMap |
Definition at line 130 of file ATSolver2D.h.
enum DGtal::ATSolver2D::CellOutputPolicy |
Specifies how to merge the different values of 0-form v at cell vertices when outputing the 0-form v for a range of cells (either pointels, linels, surfels).
Enumerator | |
---|---|
Average | compute average values at cell vertices |
Minimum | compute minimum value at cell vertices, |
Maximum | compute maximum value at cell vertices |
Definition at line 103 of file ATSolver2D.h.
|
inline |
Prepare an AT-solver from a valid calculus.
aCalculus | any valid calculus |
aVerbose | tells how the solver displays computing information: 0 none, 1 more, 2 even more... |
Definition at line 202 of file ATSolver2D.h.
References index(), DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators(), DGtal::PRIMAL, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::surfel2idx, DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
|
delete |
Default constructor.
|
default |
Destructor.
|
default |
Copy constructor.
other | the object to clone. |
|
default |
Move constructor.
other | the object to move. |
|
inline |
Debug method for checking if v is a scalar field between 0 and 1.
Definition at line 723 of file ATSolver2D.h.
References DGtal::Trace::info(), DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::selfDisplay().
|
inline |
Computes the norms loo, l2, l1 of (v - former_v), i.e. the evolution of discontinuity function v.
Definition at line 702 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::former_v0, and DGtal::KForm< TCalculus, order, duality >::myContainer.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveForEpsilon().
|
inline |
Given a range of surfels [itB,itE), returns in output the regularized scalar field u.
ScalarFieldOutput | the type of scalar field for output values (RandomAccess container) |
SurfelRangeConstIterator | the type of iterator for traversing a range of surfels |
[out] | output | the vector of output values (a scalar field), which should be of size length(itB,itE) |
itB | the start of the range of surfels. | |
itE | past the end of the range of surfels. |
Definition at line 798 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::surfel2idx.
Referenced by DGtal::ShortcutsGeometry< TKSpace >::getATScalarFieldApproximation().
|
inline |
Given a range of pointels, linels or 2-cells [itB,itE), returns in output the feature vector v (the average of v for linels/surfels).
ScalarFieldOutput | the type of scalar field for output values (RandomAccess container) |
CellRangeConstIterator | the type of iterator for traversing a range of cells |
[out] | output | the vector of output scalar values (a scalar field), which should be of size length(itB,itE) |
[in] | itB | the start of the range of cells. |
[in] | itE | past the end of the range of cells. |
[in] | policy | the chosen policy for outputing v values for a given cell. |
Definition at line 824 of file ATSolver2D.h.
References K, max(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus.
Referenced by DGtal::ShortcutsGeometry< TKSpace >::getATScalarFieldApproximation(), DGtal::ShortcutsGeometry< TKSpace >::getATVectorFieldApproximation(), and main().
|
inline |
Given a range of surfels [itB,itE), returns in output the regularized vector field u.
VectorFieldOutput | the type of vector field for output values (RandomAccess container) |
SurfelRangeConstIterator | the type of iterator for traversing a range of surfels |
[out] | output | the vector of output values (a scalar or vector field depending on input). |
itB | the start of the range of surfels. | |
itE | past the end of the range of surfels. |
Definition at line 772 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::size(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::surfel2idx.
Referenced by DGtal::ShortcutsGeometry< TKSpace >::getATVectorFieldApproximation(), and main().
|
inline |
k | an integer such that 0 <= k < u2.size() |
Definition at line 754 of file ATSolver2D.h.
|
inline |
Definition at line 735 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::v0.
|
inline |
Definition at line 741 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::M01, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::v0.
|
inline |
Definition at line 747 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::M01, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::M12, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::v0.
|
inline |
Given a range of surfels [itB,itE) and an input scalar field, initializes the AT 2-forms u and g. The 0-form v is itself initialized to 1 everywhere.
ScalarFieldInput | the type of scalar field for input values (RandomAccess container) |
SurfelRangeConstIterator | the type of iterator for traversing a range of surfels |
[in] | input | the input scalar field (a vector of scalar values) |
itB | the start of the range of surfels. | |
itE | past the end of the range of surfels. |
Definition at line 330 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::Trace::beginBlock(), DGtal::Trace::endBlock(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::g2, DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalize_u2, DGtal::KForm< TCalculus, order, duality >::ones(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::surfel2idx, DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
Referenced by DGtal::ShortcutsGeometry< TKSpace >::getATScalarFieldApproximation().
|
inline |
Given a map Surfel -> Scalar, initializes forms g, u and v of the AT solver. Note that there is only one 2-form u/g.
input | any map Surfel -> Scalar |
Scalar | any type representing a scalar (float, double) |
Definition at line 411 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::Trace::beginBlock(), DGtal::Trace::endBlock(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::g2, DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalize_u2, DGtal::KForm< TCalculus, order, duality >::ones(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::size(), DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
|
inline |
Given a map Surfel -> ScalarVector, initializes forms g, u and v of the AT solver. Note that there are as many 2-forms u/g as the number of dimensions of ScalarVector.
input | any map Surfel -> ScalarVector |
normalize | when 'true', the input is supposed to be a unit vector field and the solver will output a unit regularized vector field at the end of each minimization step. |
ScalarVector | any type representing a vector/array of scalars (float, double) |
Definition at line 371 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::Trace::beginBlock(), DGtal::Trace::endBlock(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::g2, DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalize_u2, DGtal::KForm< TCalculus, order, duality >::ones(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::size(), DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
|
inline |
Given a range of surfels [itB,itE) and an input vector field, initializes the AT 2-forms u and g. The 0-form v is itself initialized to 1 everywhere.
VectorFieldInput | the type of vector field for input values (RandomAccess container) |
SurfelRangeConstIterator | the type of iterator for traversing a range of surfels |
[in] | input | the input vector field (a vector of vector values) |
[in] | itB | the start of the range of surfels. |
[in] | itE | past the end of the range of surfels. |
[in] | normalize | when 'true', the input is supposed to be a unit vector field and the solver will output a unit regularized vector field at the end of each minimization step. |
Definition at line 288 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::Trace::beginBlock(), DGtal::Trace::endBlock(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::g2, DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalize_u2, DGtal::KForm< TCalculus, order, duality >::ones(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::surfel2idx, DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
Referenced by DGtal::ShortcutsGeometry< TKSpace >::getATVectorFieldApproximation(), and main().
|
inlineprotected |
Initializes the operators.
Definition at line 976 of file ATSolver2D.h.
References DGtal::Trace::beginBlock(), DGtal::Trace::endBlock(), DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::M01, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::M12, DGtal::LinearOperator< TCalculus, order_in, duality_in, order_out, duality_out >::myContainer, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::primal_AD2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::primal_D0, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::primal_D1, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ATSolver2D().
|
inline |
Checks the validity/consistency of the object.
Definition at line 962 of file ATSolver2D.h.
|
inline |
Forces the normalization of the vector u, meaning for all index i, \( \sum_{k=0}^{K-1} u[k][i]^2 = 1 \). Can be useful in some applications where you are looking for unitary vector field.
Definition at line 684 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::size().
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
default |
Move assignment operator.
other | the object to move. |
|
default |
Copy assignment operator.
other | the object to copy. |
|
inline |
Writes/Displays the object on an output stream.
out | the output stream where the object is written. |
Definition at line 949 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::checkV0().
|
inline |
Initializes the epsilon parameter of AT and precomputes the assaociated forms and operators.
e | the epsilon parameter in AT |
Definition at line 525 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::epsilon, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::l_1_over_4e, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::l_1_over_4e_Id0, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::lambda, DGtal::KForm< TCalculus, order, duality >::ones(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveForEpsilon().
|
inline |
Initializes the alpha and lambda parameters of AT.
a | the global alpha parameter |
l | the global lambda parameter |
Definition at line 444 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_Id2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::g2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::lambda, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus.
Referenced by DGtal::ShortcutsGeometry< TKSpace >::getATScalarFieldApproximation(), DGtal::ShortcutsGeometry< TKSpace >::getATVectorFieldApproximation(), and main().
|
inline |
Initializes the alpha and lambda parameters of AT, with weights on the 2-cells for the data terms.
AlphaWeights | the type of RandomAccess container for alpha weight values |
SurfelRangeConstIterator | the type of iterator for traversing a range of surfels |
[in] | a | the global alpha parameter |
[in] | l | the global lambda parameter |
[in] | weights | the vector of alpha weights for each surfel of the range [itB,itE) |
[in] | itB | the start of the range of surfels. |
[in] | itE | past the end of the range of surfels. |
Definition at line 499 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_Id2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::g2, DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::lambda, DGtal::KForm< TCalculus, order, duality >::myContainer, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::surfel2idx, DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
|
inline |
Initializes the alpha and lambda parameters of AT, with weights on the 2-cells for the data terms.
a | the global alpha parameter |
l | the global lambda parameter |
weights | the map Surfel -> Scalar that gives the weight of each 2-cell in the data terms. |
Definition at line 463 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_Id2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::g2, DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::lambda, DGtal::KForm< TCalculus, order, duality >::myContainer, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::size(), DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
|
inline |
order | the dimension of cells (0,1,2) |
Definition at line 258 of file ATSolver2D.h.
References DGtal::PRIMAL, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getOutputVectorFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputVectorFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalizeU2(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp().
|
inline |
Solves the alternate minimization of AT for a given eps. Solves for u then for v till convergence.
eps | the epsilon parameter at which AT is solved. |
n_oo_max | the alternate minimization will stop when the loo-norm of \( v^{k+1} - v^k \) is below this bound. |
iter_max | the alternate minimization will stop when the number of minimization steps exceeds iter_max. |
Definition at line 605 of file ATSolver2D.h.
References DGtal::Trace::beginBlock(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::diffV0(), DGtal::Trace::endBlock(), DGtal::Trace::info(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setEpsilon(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep(), DGtal::trace, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveGammaConvergence().
|
inline |
Solves AT by progressively decreasing epsilon from eps1 to eps2. AT is solved with solveForEpsilon at each epsilon.
eps1 | the first epsilon parameter at which AT is solved. |
eps2 | the last epsilon parameter at which AT is solved. |
epsr | the ratio (>1) used to decrease progressively epsilon. |
compute_smallest_epsilon_map | when 'true' determines for each surfel the smallest epsilon for which it is a discontinuity. |
n_oo_max | the alternate minimization will stop when the loo-norm of \( v^{k+1} - v^k \) is below this bound. |
iter_max | the alternate minimization will stop when the number of minimization steps exceeds iter_max. |
Definition at line 657 of file ATSolver2D.h.
References DGtal::Trace::beginBlock(), DGtal::Trace::endBlock(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::smallest_epsilon_map, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveForEpsilon(), DGtal::trace, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::updateSmallestEpsilonMap(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose.
Referenced by DGtal::ShortcutsGeometry< TKSpace >::getATScalarFieldApproximation(), DGtal::ShortcutsGeometry< TKSpace >::getATVectorFieldApproximation(), and main().
|
inline |
Solves one step of the alternate minimization of AT. Solves for u then for v.
Definition at line 546 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_g2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha_Id2, DGtal::Trace::beginBlock(), DGtal::DiscreteExteriorCalculusSolver< TCalculus, TLinearAlgebraSolver, order_in, duality_in, order_out, duality_out >::compute(), DGtal::dec_helper::diagonal(), DGtal::Trace::endBlock(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::epsilon, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::former_v0, DGtal::Trace::info(), DGtal::DiscreteExteriorCalculusSolver< TCalculus, TLinearAlgebraSolver, order_in, duality_in, order_out, duality_out >::isValid(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::l_1_over_4e, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::l_1_over_4e_Id0, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::lambda, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::M01, DGtal::KForm< TCalculus, order, duality >::myContainer, DGtal::DiscreteExteriorCalculusSolver< TCalculus, TLinearAlgebraSolver, order_in, duality_in, order_out, duality_out >::myLinearAlgebraSolver, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalize_u2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalizeU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::primal_AD2, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::primal_D0, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, DGtal::DiscreteExteriorCalculusSolver< TCalculus, TLinearAlgebraSolver, order_in, duality_in, order_out, duality_out >::solve(), DGtal::trace, DGtal::LinearOperator< TCalculus, order_in, duality_in, order_out, duality_out >::transpose(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::v0, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose, and DGtal::KForm< TCalculus, order, duality >::zeros().
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveForEpsilon().
|
inline |
Computes the map that stores for each surfel the smallest epsilon for which the surfel was in the discontinuity zone (more precisely, the surfel has at least two vertices that belongs to the set of discontinuity).
[in] | threshold | the threshold for discontinuity function v (below u is discontinuous, above u is continuous) |
Definition at line 903 of file ATSolver2D.h.
References DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::epsilon, K, DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ptrCalculus, and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::smallest_epsilon_map.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveGammaConvergence().
double DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::alpha |
The global coefficient alpha giving the smoothness of the reconstruction (the smaller, the smoother)
Definition at line 181 of file ATSolver2D.h.
|
protected |
The alpha-weighted N-array of input primal 2-forms g.
Definition at line 162 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputVectorFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
The alpha-weighted identity operator for primal 2-forms (stored for performance)
Definition at line 156 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
static |
Definition at line 99 of file ATSolver2D.h.
double DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::epsilon |
The global coefficient epsilon giving the width of the discontinuities (the smaller, the thinner)
Definition at line 187 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setEpsilon(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::updateSmallestEpsilonMap().
|
protected |
The primal 0-form v at the previous iteration.
Definition at line 169 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::diffV0(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
The N-array of input primal 2-forms g.
Definition at line 160 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputVectorFieldU2(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp().
|
protected |
The primal 0-form lambda/(4epsilon) (stored for performance)
Definition at line 171 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setEpsilon(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
The 1/(4epsilon)-weighted identity operator for primal 0-forms (stored for performance)
Definition at line 158 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setEpsilon(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
double DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::lambda |
The global coefficient lambda giving the length of discontinuities (the smaller, the more discontinuities)
Definition at line 184 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setEpsilon(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
The primal vertex to edge average operator.
Definition at line 150 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getV1(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getV2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
The primal edge to face average operator.
Definition at line 152 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getV2(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators().
bool DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::normalize_u2 |
Indicates whether to normalize U (unit norm) at each iteration or not.
Definition at line 189 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputVectorFieldU2(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
The antiderivative of primal 2-forms.
Definition at line 154 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
the derivative operator for primal 0-forms
Definition at line 146 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
|
protected |
the derivative operator for primal 1-forms
Definition at line 148 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators().
|
protected |
A smart (or not) pointer to a calculus object.
Definition at line 144 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ATSolver2D(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getOutputScalarFieldV0(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputVectorFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setEpsilon(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::size(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::updateSmallestEpsilonMap().
SmallestEpsilonMap DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::smallest_epsilon_map |
The map Surfel -> double telling the smallest epsilon for which the surfel was a discontinuity.
Definition at line 178 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveGammaConvergence(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::updateSmallestEpsilonMap().
Surfel2IndexMap DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::surfel2idx |
Definition at line 175 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ATSolver2D(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getOutputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getOutputVectorFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputVectorFieldU2(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp().
|
protected |
The N-array of regularized primal 2-forms u.
Definition at line 164 of file ATSolver2D.h.
|
protected |
The primal 0-form v representing the set of discontinuities S (more precisely \( 1 - \chi_S \))
Definition at line 167 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getV0(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getV1(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::getV2(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().
int DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::verbose |
Tells the verbose level.
Definition at line 191 of file ATSolver2D.h.
Referenced by DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::ATSolver2D(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::checkV0(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputScalarFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initInputVectorFieldU2(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::initOperators(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::setUp(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveForEpsilon(), DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveGammaConvergence(), and DGtal::ATSolver2D< TKSpace, TLinearAlgebra >::solveOneAlternateStep().