DGtal 2.1.0
Loading...
Searching...
No Matches
DGtal::AffineBasis< TPoint > Struct Template Reference

Aim: Utility class to determine the affine geometry of an input set of points. It provides exact results when the input is composed of lattice points, and may determine a basis, the dimension, or an orthogonal vector. More...

#include <DGtal/geometry/tools/AffineBasis.h>

Inheritance diagram for DGtal::AffineBasis< TPoint >:
[legend]

Public Types

enum struct  Type { INVALID = 0 , ECHELON_REDUCED , SHORTEST_ECHELON_REDUCED , LLL_REDUCED }
 
typedef AffineBasis< TPoint > Self
 
typedef TPoint Point
 
typedef Point::Coordinate Scalar
 
typedef std::vector< PointPoints
 
typedef AffineGeometry< PointAffine
 

Public Member Functions

standard services
 AffineBasis (const double tolerance=1e-12)
 
template<typename TInputPoint >
 AffineBasis (const std::vector< TInputPoint > &points, AffineBasis::Type type, const double delta=0.99, const double tolerance=1e-12)
 
template<typename TInputPoint >
 AffineBasis (const TInputPoint &origin, const std::vector< TInputPoint > &basis, AffineBasis::Type type, bool is_reduced=false, const double delta=0.99, const double tolerance=1e-12)
 
template<typename TInputPoint >
 AffineBasis (const TInputPoint &origin, const TInputPoint &normal, AffineBasis::Type type=Type::ECHELON_REDUCED, const double tolerance=1e-12)
 
void reduce (AffineBasis::Type type, double delta)
 
Dimension dimension () const
 
const Pointorigin () const
 
const Pointsbasis () const
 
debug and I/O services
void selfDisplay (std::ostream &out) const
 
bool isValid () const
 
std::string reductionTypeName () const
 

Data Fields

public data
Point first
 the origin of the affine basis
 
Points second
 the vector basis
 
double epsilon {1e-12}
 the accepted value below which a floating-point number is 0.
 
AffineBasis::Type _type
 the type of reduction of the basis.
 

Protected Member Functions

protected services
void normalize ()
 
void reduceAsEchelon (Type type)
 
void orderEchelonBasis ()
 Guarantees that the basis is in echelon form.
 
void reduceAsLLL (double delta, Scalar)
 
template<typename TInputPoint >
void initBasis (const std::vector< TInputPoint > &basis)
 
void sortBasis ()
 
std::size_t findIndexWithSmallestNonNullComponent (Dimension k, std::size_t i, const std::vector< Point > &basis)
 

geometry services

bool isParallel (const Self &other) const
 
std::pair< Scalar, PointrationalCoordinates (const Point &p) const
 
bool isOnAffineSpace (const Point &p) const
 
bool isParallel (const Point &w) const
 
Point recompose (Scalar d, const Point &lambda, const Point &r=Point::zero) const
 
Point recomposeVector (Scalar d, const Point &lambda, const Point &r=Point::zero) const
 
std::tuple< Scalar, Point, Pointdecompose (const Point &p) const
 
std::tuple< Scalar, Point, PointdecomposeVector (Point w) const
 
template<typename ProjectedPoint >
Scalar projectPoints (std::vector< ProjectedPoint > &result, const Points &input)
 
template<typename OtherPoint >
static void transform (OtherPoint &pp, const Point &p)
 
template<typename OtherPoint >
static void dilatedTransform (OtherPoint &pp, const Point &p, Scalar m)
 

Detailed Description

template<typename TPoint>
struct DGtal::AffineBasis< TPoint >

Aim: Utility class to determine the affine geometry of an input set of points. It provides exact results when the input is composed of lattice points, and may determine a basis, the dimension, or an orthogonal vector.

Description of template class 'AffineBasis'

Note
Useful for algorithms that requires the exact dimensionality of a set of points before processing, like QuickHull.
Warning
You can use this class to test the dimensionality of a set of points with floating point coordinates, but this approach is less robust than singular value decomposition.
Template Parameters
TPointthe type for points, which may be lattice points or points with floating-point coordinates.
#include <iostream>
#include <vector>
#include "DGtal/base/Common.h"
#include "DGtal/kernel/SpaceND.h"
#include "DGtal/geometry/tools/AffineBasis.h"
...
typedef Space::Point Point;
std::vector<Point> X = { Point{1, 0, 0}, Point{2, 1, 0}, Point{3, 2, 0}, Point{3, 1, 1}, Point{5, 2, 2}, Point{4, 2, 1} };
auto I = Affine::affineSubset( X );
auto B = Affine::affineBasis( X );
auto d = Affine::affineDimension( X );
Aim: Utility class to determine the affine geometry of an input set of points. It provides exact resu...
AffineGeometry< Point > Affine
static DGtal::int64_t affineDimension(const std::vector< TInputPoint > &X, const double tolerance=1e-12)
static std::pair< Point, Points > affineBasis(const std::vector< TInputPoint > &X, const double tolerance=1e-12)
static std::vector< Size > affineSubset(const std::vector< TInputPoint > &X, const double tolerance=1e-12)
See also
testAffineBasis.cpp
Note
Affine basis can be computed in several forms, which are more or less suited according to the targeted use:
  • ECHELON_REDUCED: the basis A is in row echelon form, useful for solving systems of the form Ax=b. A drawback is that component values may increase quickly if we ask for integers.
  • SHORTEST_ECHELON_REDUCED: the basis A is in row echelon form, but the computation first sorts the input vector from the shortest to the longest vectors before Gauss elimination. It is useful for solving systems of the form Ax=b. A drawback is also that component values may increase quickly if we ask for integers.
  • LLL_REDUCED (specific to lattice vectors): the basis A is not
  • in row echelon form but is made of almost the shortest
  • possible lattice vectors. It is implemented with the LLL
  • algorithm (see <a
  • href="https://en.wikipedia.org/wiki/Lenstra–Lenstra–Lovász_lattice_basis_reduction_algorithm">
  • Wikipedia LLL algorithm ). It is not handy for solving
  • systems but the basis vectors are in general much shorter than
  • in row echelon form.

Definition at line 111 of file AffineBasis.h.

Member Typedef Documentation

◆ Affine

template<typename TPoint >
typedef AffineGeometry< Point > DGtal::AffineBasis< TPoint >::Affine

Definition at line 117 of file AffineBasis.h.

◆ Point

template<typename TPoint >
typedef TPoint DGtal::AffineBasis< TPoint >::Point

Definition at line 114 of file AffineBasis.h.

◆ Points

template<typename TPoint >
typedef std::vector< Point > DGtal::AffineBasis< TPoint >::Points

Definition at line 116 of file AffineBasis.h.

◆ Scalar

template<typename TPoint >
typedef Point::Coordinate DGtal::AffineBasis< TPoint >::Scalar

Definition at line 115 of file AffineBasis.h.

◆ Self

template<typename TPoint >
typedef AffineBasis<TPoint> DGtal::AffineBasis< TPoint >::Self

Definition at line 113 of file AffineBasis.h.

Member Enumeration Documentation

◆ Type

template<typename TPoint >
enum struct DGtal::AffineBasis::Type
strong
Enumerator
INVALID 

invalid basis

ECHELON_REDUCED 

echelon matrix

SHORTEST_ECHELON_REDUCED 

echelon matrix starting from shortest vectors

LLL_REDUCED 

delta-LLL reduced matrix

Definition at line 119 of file AffineBasis.h.

119 {
120 INVALID = 0,
124 };
@ SHORTEST_ECHELON_REDUCED
echelon matrix starting from shortest vectors
@ LLL_REDUCED
delta-LLL reduced matrix
@ ECHELON_REDUCED
echelon matrix
@ INVALID
invalid basis

Constructor & Destructor Documentation

◆ AffineBasis() [1/4]

template<typename TPoint >
DGtal::AffineBasis< TPoint >::AffineBasis ( const double  tolerance = 1e-12)
inline

Default constructor. This create the identity basis of the space.

Parameters
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).

Definition at line 137 of file AffineBasis.h.

138 : epsilon( tolerance )
139 {
140 first = Point::zero;
141 second.resize( Point::dimension );
142 for ( auto k = 0; k < Point::dimension; k++ )
143 second[ k ] = Point::base( k );
145 }
AffineBasis::Type _type
the type of reduction of the basis.
Points second
the vector basis
double epsilon
the accepted value below which a floating-point number is 0.
Point first
the origin of the affine basis

References DGtal::AffineBasis< TPoint >::_type, DGtal::AffineBasis< TPoint >::first, DGtal::AffineBasis< TPoint >::second, and DGtal::AffineBasis< TPoint >::SHORTEST_ECHELON_REDUCED.

◆ AffineBasis() [2/4]

template<typename TPoint >
template<typename TInputPoint >
DGtal::AffineBasis< TPoint >::AffineBasis ( const std::vector< TInputPoint > &  points,
AffineBasis< TPoint >::Type  type,
const double  delta = 0.99,
const double  tolerance = 1e-12 
)
inline

Constructor from points.

Parameters
[in]pointsthe range of points belonging to the affine space.
[in]typeif Type::ECHELON_REDUCED or Type::SHORTEST_ECHELON_REDUCED, reduces the basis so that it forms a echelon matrix, otherwise computes its delta-LLL-lattice.
[in]deltathe parameter \( \delta \) of LLL-algorithm, which should be between 0.25 and 1 (value 0.99 is default in sagemath).
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).

Definition at line 164 of file AffineBasis.h.

168 : epsilon( tolerance )
169 {
170 if ( points.size() == 0 ) return;
171 first = Affine::transform( points[ 0 ] );
172 std::vector< TInputPoint > basis( points.size() - 1 );
173 for ( std::size_t i = 0; i < basis.size(); i++ )
174 basis[ i ] = ( points[ i+1 ] - first );
175 initBasis( basis );
176 reduce( type, delta );
177 }
void reduce(AffineBasis::Type type, double delta)
const Points & basis() const
void initBasis(const std::vector< TInputPoint > &basis)
static const Point & transform(const Point &w)

References DGtal::AffineBasis< TPoint >::basis(), DGtal::AffineBasis< TPoint >::first, DGtal::AffineBasis< TPoint >::initBasis(), DGtal::AffineBasis< TPoint >::reduce(), and DGtal::AffineGeometry< TPoint >::transform().

◆ AffineBasis() [3/4]

template<typename TPoint >
template<typename TInputPoint >
DGtal::AffineBasis< TPoint >::AffineBasis ( const TInputPoint &  origin,
const std::vector< TInputPoint > &  basis,
AffineBasis< TPoint >::Type  type,
bool  is_reduced = false,
const double  delta = 0.99,
const double  tolerance = 1e-12 
)
inline

Constructor from origin and basis.

Parameters
[in]originthe origin of the affine basis
[in]basisthe range of vectors forming the basis
[in]typeif Type::ECHELON_REDUCED or Type::SHORTEST_ECHELON_REDUCED, reduces the basis so that it forms a echelon matrix, otherwise computes its delta-LLL-lattice.
[in]is_reducedwhen 'true', assumes that the given basis is already reduced, otherwise it forces the reduction of the basis.
[in]deltathe parameter \( \delta \) of LLL-algorithm, which should be between 0.25 and 1 (value 0.99 is default in sagemath).
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).

Definition at line 201 of file AffineBasis.h.

207 : epsilon( tolerance )
208 {
210 initBasis( basis );
211 if ( ! is_reduced ) reduce( type, delta );
212 else _type = type;
213 }
const Point & origin() const

References DGtal::AffineBasis< TPoint >::_type, DGtal::AffineBasis< TPoint >::basis(), DGtal::AffineBasis< TPoint >::first, DGtal::AffineBasis< TPoint >::initBasis(), DGtal::AffineBasis< TPoint >::origin(), DGtal::AffineBasis< TPoint >::reduce(), and DGtal::AffineGeometry< TPoint >::transform().

◆ AffineBasis() [4/4]

template<typename TPoint >
template<typename TInputPoint >
DGtal::AffineBasis< TPoint >::AffineBasis ( const TInputPoint &  origin,
const TInputPoint &  normal,
AffineBasis< TPoint >::Type  type = Type::ECHELON_REDUCED,
const double  tolerance = 1e-12 
)
inline

Creates an affine basis going through origin and orthogonal to the lattice vector normal.

Parameters
[in]originthe origin of the affine basis
[in]normalthe lattice normal vector.
[in]typeif Type::ECHELON_REDUCED or Type::SHORTEST_ECHELON_REDUCED, then the basis will be in echelon form, otherwise it will make the vectors as short as possible, but the matrix won't be in echelon form.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).

Definition at line 230 of file AffineBasis.h.

234 : epsilon( tolerance )
235 {
237 // basis is shortened is type is LLL.
239 _type = ( type == Type::LLL_REDUCED )
241 }
static std::vector< Point > orthogonalLatticeBasis(const TInputPoint &N, bool shortened=false)

References DGtal::AffineBasis< TPoint >::_type, DGtal::AffineBasis< TPoint >::ECHELON_REDUCED, DGtal::AffineBasis< TPoint >::first, DGtal::AffineBasis< TPoint >::LLL_REDUCED, DGtal::AffineBasis< TPoint >::origin(), DGtal::AffineGeometry< TPoint >::orthogonalLatticeBasis(), DGtal::AffineBasis< TPoint >::second, and DGtal::AffineGeometry< TPoint >::transform().

Member Function Documentation

◆ basis()

template<typename TPoint >
const Points & DGtal::AffineBasis< TPoint >::basis ( ) const
inline

◆ decompose()

template<typename TPoint >
std::tuple< Scalar, Point, Point > DGtal::AffineBasis< TPoint >::decompose ( const Point p) const
inline

Decompose p as d * (p - o) = l[0]b[0] + ... l[i]b[i] + r, where r is independent from this basis B=(b[0],...,b[i]).

Parameters
[in]pany point.
Returns
(d,l,r) where d is a scalar, l/d are the rational coordinates of (p-o) in the basis, r is the remainer vector if p does not belong to the affine space.

Definition at line 401 of file AffineBasis.h.

402 {
403 return decomposeVector( p - first );
404 }
std::tuple< Scalar, Point, Point > decomposeVector(Point w) const

References DGtal::AffineBasis< TPoint >::decomposeVector(), and DGtal::AffineBasis< TPoint >::first.

Referenced by DGtal::AffineBasis< TPoint >::isOnAffineSpace(), DGtal::AffineBasis< TPoint >::projectPoints(), and DGtal::AffineBasis< TPoint >::rationalCoordinates().

◆ decomposeVector()

template<typename TPoint >
std::tuple< Scalar, Point, Point > DGtal::AffineBasis< TPoint >::decomposeVector ( Point  w) const
inline

Decompose w as d * w = l[0]b[0] + ... l[i]b[i] + r, where r is independent from B=(b[0],...,b[i]).

Parameters
[in]wany vector.
Returns
(d,l,r) where d is a scalar, l/d are the rational coordinates of w in the basis, r is the remainer vector if w does not belong to the affine space.

Definition at line 414 of file AffineBasis.h.

415 {
416 Point r;
417 Scalar alphas = 1;
418 for ( auto i = 0; i < second.size(); i++ )
419 {
420 std::pair< Scalar, Scalar > c
421 = Affine::reduceVector( w, second[ i ], i, epsilon );
422 for ( auto j = 0; j < i; j++ )
423 r[ j ] *= c.first;
424 r[ i ] = c.second;
425 alphas *= c.first;
426 }
427 return alphas >= 0
428 ? std::make_tuple( alphas, r, w )
429 : std::make_tuple( -alphas, -r, -w );
430 }
STL namespace.
static std::pair< Scalar, Scalar > reduceVector(Point &w, const Point &b, const double tolerance)

References DGtal::AffineBasis< TPoint >::epsilon, DGtal::AffineGeometry< TPoint >::reduceVector(), and DGtal::AffineBasis< TPoint >::second.

Referenced by DGtal::AffineBasis< TPoint >::decompose(), and DGtal::AffineBasis< TPoint >::isParallel().

◆ dilatedTransform()

template<typename TPoint >
template<typename OtherPoint >
static void DGtal::AffineBasis< TPoint >::dilatedTransform ( OtherPoint &  pp,
const Point p,
Scalar  m 
)
inlinestatic

Transforms the type of an input point into another one, while dilating it by a factor m.

Template Parameters
OtherPointa type of point of dimension at most Point::dimension.
Parameters
[out]ppthe output restricted point.
[in]pthe input point.
[in]mthe dilation factor.

Definition at line 492 of file AffineBasis.h.

493 {
494 BOOST_STATIC_ASSERT( OtherPoint::dimension <= Point::dimension );
495 for ( std::size_t i = 0; i < pp.dimension; ++i )
496 pp[ i ] = m * Scalar( p[ i ] );
497 }
Point::Coordinate Scalar

Referenced by DGtal::AffineBasis< TPoint >::projectPoints().

◆ dimension()

template<typename TPoint >
Dimension DGtal::AffineBasis< TPoint >::dimension ( ) const
inline
Returns
the affine dimension of the basis
Precondition
the basis is reduced.

Definition at line 271 of file AffineBasis.h.

272 {
273 return second.size();
274 }

References DGtal::AffineBasis< TPoint >::second.

Referenced by DGtal::AffineBasis< TPoint >::isParallel().

◆ findIndexWithSmallestNonNullComponent()

template<typename TPoint >
std::size_t DGtal::AffineBasis< TPoint >::findIndexWithSmallestNonNullComponent ( Dimension  k,
std::size_t  i,
const std::vector< Point > &  basis 
)
inlineprotected

Given a range of points basis, starting from rank i, find the index of the point with lowest non null k-th coefficient in absolute value, or basis.size() if every point had its k-th component null.

Parameters
[in]kthe component/coordinate of interest
[in]ithe starting index
[in]basisa range of points/vectors
Returns
starting from rank i, the index of the point with lowest non null k-th coefficient in absolute value, or basis.size() if every point had its k-th component null.

Definition at line 708 of file AffineBasis.h.

711 {
712 ASSERT( ! basis.empty() );
713 ASSERT( k < Point::dimension );
714 std::size_t index = i;
715 Scalar v = 0;
716 for ( ; index < basis.size(); index++ )
717 {
718 v = abs( basis[ index ][ k ] );
719 if ( Affine::ScalarOps::isNonZero( v, epsilon ) )
720 break;
721 }
722 for ( auto j = index + 1; j < basis.size(); j++ )
723 {
724 Scalar vj = abs( basis[ j ][ k ] );
725 if ( vj != 0 && vj < v )
726 {
727 index = j;
728 v = vj;
729 }
730 }
731 return index;
732 }
unsigned int index(DGtal::uint32_t n, unsigned int b)
Definition testBits.cpp:44

References DGtal::AffineBasis< TPoint >::basis(), DGtal::AffineBasis< TPoint >::epsilon, and index().

Referenced by DGtal::AffineBasis< TPoint >::reduceAsEchelon().

◆ initBasis()

template<typename TPoint >
template<typename TInputPoint >
void DGtal::AffineBasis< TPoint >::initBasis ( const std::vector< TInputPoint > &  basis)
inlineprotected

Removes null vectors from input set of vectors and sets the initial basis.

Template Parameters
TInputPointthe type of input points
Parameters
[in]basisa set of arbitrary vectors

Definition at line 660 of file AffineBasis.h.

661 {
662 second.reserve( basis.size() );
663 for ( auto i = 0; i < basis.size(); i++ )
664 {
665 Point b = Affine::transform( basis[ i ] );
666 if ( b != Point::zero ) second.push_back( b );
667 }
668 }

References DGtal::AffineBasis< TPoint >::basis(), DGtal::AffineBasis< TPoint >::second, and DGtal::AffineGeometry< TPoint >::transform().

Referenced by DGtal::AffineBasis< TPoint >::AffineBasis(), and DGtal::AffineBasis< TPoint >::AffineBasis().

◆ isOnAffineSpace()

template<typename TPoint >
bool DGtal::AffineBasis< TPoint >::isOnAffineSpace ( const Point p) const
inline
Parameters
[in]pany lattice point
Returns
'true' if and only if the point belongs to the affine space spanned by 'this'.

Definition at line 326 of file AffineBasis.h.

327 {
328 const auto [d, lambda, r] = decompose( p );
329 return ! Affine::ScalarOps::isNonZero( r.normInfinity(), epsilon );
330 }
std::tuple< Scalar, Point, Point > decompose(const Point &p) const

References DGtal::AffineBasis< TPoint >::decompose(), and DGtal::AffineBasis< TPoint >::epsilon.

◆ isParallel() [1/2]

template<typename TPoint >
bool DGtal::AffineBasis< TPoint >::isParallel ( const Point w) const
inline
Parameters
[in]wany lattice vector
Returns
'true' if and only if the vector w is parallel to the affine space spanned by 'this'.

Definition at line 336 of file AffineBasis.h.

337 {
338 const auto [d, lambda, r] = decomposeVector( w );
339 return ! Affine::ScalarOps::isNonZero( r.normInfinity(), epsilon );
340 }

References DGtal::AffineBasis< TPoint >::decomposeVector(), and DGtal::AffineBasis< TPoint >::epsilon.

◆ isParallel() [2/2]

template<typename TPoint >
bool DGtal::AffineBasis< TPoint >::isParallel ( const Self other) const
inline
Parameters
[in]otheran other affine basis
Returns
'true' if and only if every vector of 'this' basis is parallel to every vector of the other given basis.

Definition at line 299 of file AffineBasis.h.

300 {
301 if ( dimension() != other.dimension() ) return false;
302 if ( ( _type != Type::ECHELON_REDUCED )
304 trace.error() << "[AffineBasis::isParallel] Requires type=*_ECHELON_REDUCED\n"
305 << " type=" << reductionTypeName() << "\n" ;
306 for ( const auto& b : other.second )
307 if ( ! isParallel( b ) ) return false;
308 return true;
309 }
std::ostream & error()
Trace trace
bool isParallel(const Self &other) const
Dimension dimension() const
std::string reductionTypeName() const

References DGtal::AffineBasis< TPoint >::_type, DGtal::AffineBasis< TPoint >::dimension(), DGtal::AffineBasis< TPoint >::ECHELON_REDUCED, DGtal::Trace::error(), DGtal::AffineBasis< TPoint >::isParallel(), DGtal::AffineBasis< TPoint >::reductionTypeName(), DGtal::AffineBasis< TPoint >::second, DGtal::AffineBasis< TPoint >::SHORTEST_ECHELON_REDUCED, and DGtal::trace.

Referenced by DGtal::AffineBasis< TPoint >::isParallel().

◆ isValid()

template<typename TPoint >
bool DGtal::AffineBasis< TPoint >::isValid ( ) const
inline
Returns
'true' if and only if this object is in a consistent state. Here it means that the matrix is echelon or LLL-reduced.

Definition at line 520 of file AffineBasis.h.

521 {
522 return _type != Type::INVALID;
523 }

References DGtal::AffineBasis< TPoint >::_type, and DGtal::AffineBasis< TPoint >::INVALID.

◆ normalize()

template<typename TPoint >
void DGtal::AffineBasis< TPoint >::normalize ( )
inlineprotected

If the basis is an integer lattice, reduces the basis vectors by their gcd, otherwise normalize vectors to have 1 L2-norm.

Definition at line 554 of file AffineBasis.h.

555 {
556 for ( auto& v : second )
557 v = Affine::simplifiedVector( v );
558 }

References DGtal::AffineBasis< TPoint >::second, and DGtal::AffineGeometry< TPoint >::simplifiedVector().

Referenced by DGtal::AffineBasis< TPoint >::sortBasis().

◆ orderEchelonBasis()

template<typename TPoint >
void DGtal::AffineBasis< TPoint >::orderEchelonBasis ( )
inlineprotected

Guarantees that the basis is in echelon form.

Definition at line 596 of file AffineBasis.h.

597 {
598 auto compare = [this] ( const Point& v, const Point& w ) -> bool
599 {
600 // Note: curiously std::sort sometimes test v against itself
601 // and must return false in this case.
602 for ( auto k = 0; k < Point::dimension; ++k )
603 {
604 bool v_non_null = Affine::ScalarOps::isNonZero( v[ k ], epsilon );
605 bool w_non_null = Affine::ScalarOps::isNonZero( w[ k ], epsilon );
606 if ( v_non_null && ! w_non_null ) return true;
607 else if ( ! v_non_null && w_non_null ) return false;
608 else if ( v_non_null && w_non_null ) return false; // ==
609 }
610 return false; // 0 == 0
611 };
612 std::sort( second.begin(), second.end(), compare );
613 }
bool compare(const Range1 &pts, const Range2 &groundTruth)
Definition testFP.cpp:98

References compare(), DGtal::AffineBasis< TPoint >::epsilon, and DGtal::AffineBasis< TPoint >::second.

Referenced by DGtal::AffineBasis< TPoint >::reduceAsEchelon().

◆ origin()

template<typename TPoint >
const Point & DGtal::AffineBasis< TPoint >::origin ( ) const
inline
Returns
the origin of this affine basis.

Definition at line 277 of file AffineBasis.h.

278 {
279 return first;
280 }

References DGtal::AffineBasis< TPoint >::first.

Referenced by DGtal::AffineBasis< TPoint >::AffineBasis(), and DGtal::AffineBasis< TPoint >::AffineBasis().

◆ projectPoints()

template<typename TPoint >
template<typename ProjectedPoint >
Scalar DGtal::AffineBasis< TPoint >::projectPoints ( std::vector< ProjectedPoint > &  result,
const Points input 
)
inline

Projects the range of points input onto the affine basis and outputs it in result. A consistent choice for ProjectedPoint is to match the affine dimension.

Template Parameters
ProjectedPointa type of point that matches the affine dimension.
Parameters
[out]resultthe range of projected points.
[in]inputthe range of input points
Returns
the maximum dilation of projected rational coordinates that was used to create integer coordinates.

Definition at line 444 of file AffineBasis.h.

446 {
447 Scalar lcm = 1;
448 std::vector< Scalar > denoms ( input.size() );
449 std::vector< Point > lambdas( input.size() );
450 Point r;
451 // get points rational coordinates.
452 for ( std::size_t i = 0; i < input.size(); i++ )
453 std::tie( denoms[ i ], lambdas[ i ], r ) = decompose( input[ i ] );
454 // compute ppcm
455 for ( std::size_t i = 0; i < denoms.size(); i++ )
456 {
457 const Scalar d = denoms[ i ];
458 if ( d == 1 ) continue;
459 lcm = Affine::ScalarOps::lcmPositive( lcm, d );
460 }
461 // project points
462 result.resize( input.size() );
463 for ( std::size_t i = 0; i < input.size(); i++ )
464 dilatedTransform( result[ i ], lambdas[ i ], lcm / denoms[ i ] );
465 return lcm;
466 }
static void dilatedTransform(OtherPoint &pp, const Point &p, Scalar m)

References DGtal::AffineBasis< TPoint >::decompose(), and DGtal::AffineBasis< TPoint >::dilatedTransform().

◆ rationalCoordinates()

template<typename TPoint >
std::pair< Scalar, Point > DGtal::AffineBasis< TPoint >::rationalCoordinates ( const Point p) const
inline
Parameters
[in]pany lattice point
Returns
its rational coordinates L/d as a pair (d, L) in this affine basis, where L is a lattice vector and d is the common denominator for the components of L.

Definition at line 316 of file AffineBasis.h.

317 {
318 const auto [d, lambda, remainder] = decompose( p );
319 return std::make_pair( d, lambda );
320 }

References DGtal::AffineBasis< TPoint >::decompose().

◆ recompose()

template<typename TPoint >
Point DGtal::AffineBasis< TPoint >::recompose ( Scalar  d,
const Point lambda,
const Point r = Point::zero 
) const
inline

Given (d, lambda, r), recompose the point p such that d * (p - o) = l[0]b[0] + ... l[i]b[i] + r. Given only the rational coordinates (d,lambda), the point is assumed to lie on this affine space (i.e. r == 0 ).

Parameters
[in]dthe common denominator that defines the rational coordinates with lambda.
[in]lambdathe numerators that defines the rational coordinates with d.
[in]rthe remainder vector (r is independent from this basis B=(b[0],...,b[i])), which represents the displacement from p to this affine space.
Returns
the point p such that d * (p - o) = B lambda + r, for B the vectors of this basis (arranged as rows in matrix) and o its origin.

Definition at line 360 of file AffineBasis.h.

362 {
363 return first + recomposeVector( d, lambda, r );
364 }
Point recomposeVector(Scalar d, const Point &lambda, const Point &r=Point::zero) const

References DGtal::AffineBasis< TPoint >::first, and DGtal::AffineBasis< TPoint >::recomposeVector().

◆ recomposeVector()

template<typename TPoint >
Point DGtal::AffineBasis< TPoint >::recomposeVector ( Scalar  d,
const Point lambda,
const Point r = Point::zero 
) const
inline

Given (d, lambda, r), recompose the vector w such that d * w = l[0]b[0] + ... l[i]b[i] + r. Given only the rational coordinates (d,lambda), the vector is assumed to be parallel to this affine space (i.e. r == 0 ).

Parameters
[in]dthe common denominator that defines the rational coordinates with lambda.
[in]lambdathe numerators that defines the rational coordinates with d.
[in]rthe remainder vector (r is independent from this basis B=(b[0],...,b[i])), which represents the displacement vector to this vector space.
Returns
the point p such that d * w = B lambda + r, for B the vectors of this basis (arranged as rows in matrix) and o its origin.

Definition at line 384 of file AffineBasis.h.

386 {
387 Point w = r;
388 for ( std::size_t i = 0; i < second.size(); i++ )
389 w += lambda[ i ] * second[ i ];
390 return w / d;
391 }

References DGtal::AffineBasis< TPoint >::second.

Referenced by DGtal::AffineBasis< TPoint >::recompose().

◆ reduce()

template<typename TPoint >
void DGtal::AffineBasis< TPoint >::reduce ( AffineBasis< TPoint >::Type  type,
double  delta 
)
inline

Reduces the basis into a set of a linearly independent vectors, and in the desired reduced form.

Parameters
[in]typethe desired type of matrix reduction.
[in]deltathe parameter \( \delta \) of LLL-algorithm, which should be between 0.25 and 1 (value 0.99 is default in sagemath).

Definition at line 252 of file AffineBasis.h.

253 {
254 if ( type == Type::SHORTEST_ECHELON_REDUCED )
255 sortBasis();
257 reduceAsEchelon( type );
258 else if ( type == Type::LLL_REDUCED )
259 {
260 std::vector< Point > X( second.size()+1 );
261 X[ 0 ] = first;
262 for ( auto i = 0; i < second.size(); i++ )
263 X[ i+1 ] = second[ i ] + first;
265 reduceAsLLL( delta, (Scalar) 0 );
266 }
267 }
void reduceAsEchelon(Type type)
void reduceAsLLL(double delta, Scalar)

References DGtal::AffineGeometry< TPoint >::affineBasis(), DGtal::AffineBasis< TPoint >::ECHELON_REDUCED, DGtal::AffineBasis< TPoint >::epsilon, DGtal::AffineBasis< TPoint >::first, DGtal::AffineBasis< TPoint >::LLL_REDUCED, DGtal::AffineBasis< TPoint >::reduceAsEchelon(), DGtal::AffineBasis< TPoint >::reduceAsLLL(), DGtal::AffineBasis< TPoint >::second, DGtal::AffineBasis< TPoint >::SHORTEST_ECHELON_REDUCED, and DGtal::AffineBasis< TPoint >::sortBasis().

Referenced by DGtal::AffineBasis< TPoint >::AffineBasis(), and DGtal::AffineBasis< TPoint >::AffineBasis().

◆ reduceAsEchelon()

template<typename TPoint >
void DGtal::AffineBasis< TPoint >::reduceAsEchelon ( Type  type)
inlineprotected

Reduces the basis so that each basis vector is normalized, removes linearly dependent vectors, and builds a echelon matrix.

Definition at line 562 of file AffineBasis.h.

563 {
564 std::vector< bool > is_independent( second.size(), false );
565 std::vector< std::vector< Scalar > > U( second.size() );
566 Dimension k = 0; // the current column to put in echelon form.
567 for ( std::size_t i = 0; i < second.size(); i++ )
568 {
569 std::size_t row = findIndexWithSmallestNonNullComponent( k, i, second );
570 if ( row != i && row != second.size() )
571 std::swap( second[ i ], second[ row ] );
572 Point& w = second[ i ];
573 // check if this vector is independent from the previous ones
574 is_independent[ i ] = true;
575 for ( std::size_t j = 0; j < i; j++ )
576 if ( is_independent[ j ] )
578 if ( ! Affine::ScalarOps::isNonZero( w.normInfinity(), epsilon ) )
579 is_independent[ i ] = false; // not independent, forget it
580 else
581 { // independent, make sure it is reduced.
583 k++;
584 }
585 }
586 Points new_basis;
587 for ( std::size_t i = 0; i < second.size(); i++ )
588 if ( is_independent[ i ] )
589 new_basis.push_back( second[ i ] );
590 std::swap( second, new_basis );
592 _type = type;
593 }
DGtal::uint32_t Dimension
Definition Common.h:119
std::vector< Point > Points
void orderEchelonBasis()
Guarantees that the basis is in echelon form.
std::size_t findIndexWithSmallestNonNullComponent(Dimension k, std::size_t i, const std::vector< Point > &basis)
static Point simplifiedVector(Point v)

References DGtal::AffineBasis< TPoint >::_type, DGtal::AffineBasis< TPoint >::epsilon, DGtal::AffineBasis< TPoint >::findIndexWithSmallestNonNullComponent(), DGtal::AffineBasis< TPoint >::orderEchelonBasis(), DGtal::AffineGeometry< TPoint >::reduceVector(), DGtal::AffineBasis< TPoint >::second, and DGtal::AffineGeometry< TPoint >::simplifiedVector().

Referenced by DGtal::AffineBasis< TPoint >::reduce().

◆ reduceAsLLL()

template<typename TPoint >
void DGtal::AffineBasis< TPoint >::reduceAsLLL ( double  delta,
Scalar   
)
inlineprotected

Reduces the basis so that each basis vector is normalized, then computes its delta-LLL-reduction lattice, and removes linearly dependent vectors.

Parameters
[in]deltathe parameter \( \delta \) of LLL-algorithm, which should be between 0.25 and 1 (value 0.99 is default in sagemath).

keep only independent vectors in basis

Definition at line 622 of file AffineBasis.h.

623 {
624 if constexpr( std::is_floating_point< Scalar >::value == true )
625 {
626 trace.error() << "[AffineBasis::reduceAsLLL]"
627 << " It has no meaning to use LLL algorithm on matrix with double coefficients\n";
629 return;
630 }
631 for ( auto& v : second )
632 v = Affine::simplifiedVector( v );
633 std::vector< std::vector< Scalar > > B( second.size() );
634 for ( auto i = 0; i < second.size(); i++ )
635 {
636 B[ i ] = std::vector< Scalar >( Point::dimension );
637 for ( auto j = 0; j < Point::dimension; j++ )
638 B[ i ][ j ] = second[ i ][ j ];
639 }
642 second.clear();
643 for ( std::size_t i = 0; i < B.size(); i++ )
644 {
645 Point b;
646 for ( auto j = 0; j < Point::dimension; j++ )
647 b[ j ] = B[ i ][ j ];
648 if ( b != Point::zero )
649 second.push_back( Affine::simplifiedVector( b ) );
650 }
652 }
void reduceBasisWithLLL(std::vector< std::vector< TComponent > > &B, TDouble delta=0.75)

References DGtal::AffineBasis< TPoint >::_type, DGtal::Trace::error(), DGtal::AffineBasis< TPoint >::INVALID, DGtal::AffineBasis< TPoint >::LLL_REDUCED, DGtal::functions::reduceBasisWithLLL(), DGtal::AffineBasis< TPoint >::second, DGtal::AffineGeometry< TPoint >::simplifiedVector(), and DGtal::trace.

Referenced by DGtal::AffineBasis< TPoint >::reduce().

◆ reductionTypeName()

template<typename TPoint >
std::string DGtal::AffineBasis< TPoint >::reductionTypeName ( ) const
inline
Returns
the type of matrix as a string

Definition at line 526 of file AffineBasis.h.

527 {
528 if ( _type == Type::INVALID ) return "INVALID";
529 else if ( _type == Type::ECHELON_REDUCED ) return "ECHELON_REDUCED";
530 else if ( _type == Type::SHORTEST_ECHELON_REDUCED ) return "SHORTEST_ECHELON_REDUCED";
531 else if ( _type == Type::LLL_REDUCED ) return "LLL_REDUCED";
532 else return "";
533 }

References DGtal::AffineBasis< TPoint >::_type, DGtal::AffineBasis< TPoint >::ECHELON_REDUCED, DGtal::AffineBasis< TPoint >::INVALID, DGtal::AffineBasis< TPoint >::LLL_REDUCED, and DGtal::AffineBasis< TPoint >::SHORTEST_ECHELON_REDUCED.

Referenced by DGtal::AffineBasis< TPoint >::isParallel(), and DGtal::AffineBasis< TPoint >::selfDisplay().

◆ selfDisplay()

template<typename TPoint >
void DGtal::AffineBasis< TPoint >::selfDisplay ( std::ostream &  out) const
inline

Displays this object on the output stream.

Parameters
[in,out]outany output stream.

Definition at line 509 of file AffineBasis.h.

510 {
511 out << "[ AffineBasis o=" << first
512 << " type=" << reductionTypeName()
513 << " B=";
514 for ( auto b : second ) std::cout << "\n " << b;
515 out << " ]";
516 }

References DGtal::AffineBasis< TPoint >::first, DGtal::AffineBasis< TPoint >::reductionTypeName(), and DGtal::AffineBasis< TPoint >::second.

Referenced by DGtal::operator<<().

◆ sortBasis()

template<typename TPoint >
void DGtal::AffineBasis< TPoint >::sortBasis ( )
inlineprotected

Simplifies vectors, removes duplicates and puts smallest candidate basis vectors before longest.

Definition at line 672 of file AffineBasis.h.

673 {
674 // Reduces all vectors
675 normalize();
676 // Purge duplicates
677 std::sort( second.begin(), second.end() );
678 second.erase( std::unique( second.begin(), second.end() ), second.end() );
679 // Sort according to size of components.
680 auto compare = []( const Point& u, const Point& v ) -> bool
681 {
682 const auto n1_u = u.norm1();
683 const auto n1_v = v.norm1();
684 if ( n1_u < n1_v ) return true;
685 else if ( n1_v < n1_u ) return false;
686 const auto noo_u = u.normInfinity();
687 const auto noo_v = v.normInfinity();
688 if ( noo_u < noo_v ) return true;
689 else if ( noo_v < noo_u ) return false;
690 return u < v;
691 };
692 std::sort( second.begin(), second.end(), compare );
693 }

References compare(), DGtal::AffineBasis< TPoint >::normalize(), and DGtal::AffineBasis< TPoint >::second.

Referenced by DGtal::AffineBasis< TPoint >::reduce().

◆ transform()

template<typename TPoint >
template<typename OtherPoint >
static void DGtal::AffineBasis< TPoint >::transform ( OtherPoint &  pp,
const Point p 
)
inlinestatic

Transforms the type of an input point into another one.

Template Parameters
OtherPointa type of point of dimension at most Point::dimension.
Parameters
[out]ppthe output restricted point.
[in]pthe input point.

Definition at line 475 of file AffineBasis.h.

476 {
477 BOOST_STATIC_ASSERT( OtherPoint::dimension <= Point::dimension );
478 typedef typename OtherPoint::Coordinate Scalar;
479 for ( std::size_t i = 0; i < pp.dimension; ++i )
480 pp[ i ] = Scalar( p[ i ] );
481 }

Field Documentation

◆ _type

◆ epsilon

◆ first

◆ second


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