DGtal 2.1.0
Loading...
Searching...
No Matches
DGtal::AffineGeometry< 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/AffineGeometry.h>

Public Types

typedef TPoint Point
 
typedef Point::Coordinate Scalar
 
typedef Point::Container Container
 
typedef std::size_t Size
 
typedef std::vector< PointPoints
 type for range of points.
 
typedef std::vector< SizeSizes
 type for range of sizes.
 
typedef DGtal::detail::AffineGeometryPointOperations< dimension, Scalar, ContainerPointOps
 
typedef DGtal::detail::AffineGeometryScalarOperations< ScalarScalarOps
 

Static Public Member Functions

static affine services
template<typename TInputPoint >
static DGtal::int64_t affineDimension (const std::vector< TInputPoint > &X, const double tolerance=1e-12)
 
template<typename TInputPoint >
static std::vector< SizeaffineSubset (const std::vector< TInputPoint > &X, const double tolerance=1e-12)
 
template<typename TInputPoint , typename TIndexRange >
static std::vector< SizeaffineSubset (const std::vector< TInputPoint > &X, const TIndexRange &I, const double tolerance=1e-12)
 
template<typename TInputPoint >
static std::pair< Point, PointsaffineBasis (const std::vector< TInputPoint > &X, const double tolerance=1e-12)
 
template<typename TInputPoint , typename TIndexRange >
static std::pair< TInputPoint, PointsaffineBasis (const std::vector< TInputPoint > &X, const TIndexRange &I, const double tolerance=1e-12)
 
static Point independentVector (const Points &basis, const double tolerance=1e-12)
 
template<typename TOtherPoint >
static TOtherPoint independentVector (const Points &basis, const double tolerance=1e-12)
 
static void completeBasis (Points &basis, bool normal_vector, const double tolerance=1e-12)
 
template<typename TInputPoint >
static std::vector< PointorthogonalLatticeBasis (const TInputPoint &N, bool shortened=false)
 
static specific services

(Used internally)

static Point reductionOnBasis (const Point &v, const Points &basis, const double tolerance)
 
static bool addIfIndependent (Points &basis, const Point &v, const double tolerance)
 
static std::pair< Scalar, ScalarreduceVector (Point &w, const Point &b, const double tolerance)
 
static std::pair< Scalar, ScalarreduceVector (Point &w, const Point &b, Dimension start, const double tolerance)
 
static const Pointtransform (const Point &w)
 
template<typename TInputPoint >
static Point transform (const TInputPoint &w)
 
static Point orthogonalVector (const Points &basis)
 
template<typename TInternalNumber >
static Point orthogonalVector (const Points &basis)
 
static Point simplifiedVector (Point v)
 

Static Public Attributes

static const Size dimension = Point::dimension
 

Detailed Description

template<typename TPoint>
struct DGtal::AffineGeometry< 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 'AffineGeometry'

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/AffineGeometry.h"
...
typedef Space::Point Point;
typedef AffineGeometry< Point > Affine;
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...
See also
testAffineGeometry.cpp

Definition at line 453 of file AffineGeometry.h.

Member Typedef Documentation

◆ Container

template<typename TPoint >
typedef Point::Container DGtal::AffineGeometry< TPoint >::Container

Definition at line 457 of file AffineGeometry.h.

◆ Point

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

Definition at line 455 of file AffineGeometry.h.

◆ PointOps

template<typename TPoint >
typedef DGtal::detail::AffineGeometryPointOperations< dimension, Scalar, Container > DGtal::AffineGeometry< TPoint >::PointOps

Definition at line 463 of file AffineGeometry.h.

◆ Points

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

type for range of points.

Definition at line 459 of file AffineGeometry.h.

◆ Scalar

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

Definition at line 456 of file AffineGeometry.h.

◆ ScalarOps

template<typename TPoint >
typedef DGtal::detail::AffineGeometryScalarOperations< Scalar > DGtal::AffineGeometry< TPoint >::ScalarOps

Definition at line 464 of file AffineGeometry.h.

◆ Size

template<typename TPoint >
typedef std::size_t DGtal::AffineGeometry< TPoint >::Size

Definition at line 458 of file AffineGeometry.h.

◆ Sizes

template<typename TPoint >
typedef std::vector< Size > DGtal::AffineGeometry< TPoint >::Sizes

type for range of sizes.

Definition at line 460 of file AffineGeometry.h.

Member Function Documentation

◆ addIfIndependent()

template<typename TPoint >
static bool DGtal::AffineGeometry< TPoint >::addIfIndependent ( Points basis,
const Point v,
const double  tolerance 
)
inlinestatic

Checks if the vector v can be decomposed as a linear combination of the vectors of the basis. If yes, returns 'false', otherwise returns 'true' and adds the obtained reduced vector to the basis after normalizing it.

Parameters
[in,out]basisa range of vectors forming a (partial or not) basis of the space.
[in]vany vector.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
'true' iff the vector v was not a linear combination of the vectors of the basis and the basis is then enriched by the direction indicated by v, otherwise returns 'false'.
Note
If the points have integer coordinates, the normalization of the new basis vector is its reduction by its gcd, otherwise the vector has unit L2-norm.

Definition at line 834 of file AffineGeometry.h.

837 {
838 Point w = reductionOnBasis( v, basis, tolerance );
839 Scalar x;
841 if ( ScalarOps::isNonZero( x, tolerance ) )
842 {
843 // Useful to reduce the norm of vectors for lattice vectors
844 // and necessary for real vectors so that `tolerance` keeps
845 // the same meaning.
847 basis.push_back( w );
848 return true;
849 }
850 return false;
851 }
void getSquaredNormL2(TOutput &n, const std::vector< T > &a)
static Point reductionOnBasis(const Point &v, const Points &basis, const double tolerance)
static void normalizeVector(Point &w, TInteger)
static bool isNonZero(Integer x, double)

References DGtal::functions::getSquaredNormL2(), DGtal::detail::AffineGeometryScalarOperations< TScalar >::isNonZero(), DGtal::detail::AffineGeometryPointOperations< dim, TEuclideanRing, TContainer >::normalizeVector(), and DGtal::AffineGeometry< TPoint >::reductionOnBasis().

Referenced by DGtal::AffineGeometry< TPoint >::affineBasis(), DGtal::AffineGeometry< TPoint >::affineBasis(), DGtal::AffineGeometry< TPoint >::affineSubset(), and DGtal::AffineGeometry< TPoint >::affineSubset().

◆ affineBasis() [1/2]

template<typename TPoint >
template<typename TInputPoint >
static std::pair< Point, Points > DGtal::AffineGeometry< TPoint >::affineBasis ( const std::vector< TInputPoint > &  X,
const double  tolerance = 1e-12 
)
inlinestatic

Given a range of points X, returns a point and a range of vectors forming an affine basis containing X.

Template Parameters
TInputPointthe type of input points (may be less precise).
Parameters
[in]Xthe range of input points (may be lattice points or not).
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
a point and a range of vectors forming an affine basis containing X.
Note
Complexity is O( m n^2 ), where m=Cardinal(X) and n=dimension.

Definition at line 589 of file AffineGeometry.h.

590 {
591 Points basis; //< direction vectors
592 Size m = X.size();
593 // Process trivial cases.
594 if ( m == 0 ) return std::make_pair( Point::zero, basis );
595 const Point o = transform( X[ 0 ] );
596 if ( m == 1 ) return std::make_pair( o, basis );
597 // Process general case.
598 basis.reserve( dimension );
599 for ( Size i = 1; i < m; i++ )
600 {
601 Point v = transform( X[ i ] ) - o;
602 if ( addIfIndependent( basis, v, tolerance ) )
603 if ( basis.size() >= dimension ) break;
604 }
605 return std::make_pair( o, basis );
606 }
static const Point & transform(const Point &w)
static const Size dimension
static bool addIfIndependent(Points &basis, const Point &v, const double tolerance)
std::vector< Point > Points
type for range of points.
HalfEdgeDataStructure::Size Size

References DGtal::AffineGeometry< TPoint >::addIfIndependent(), DGtal::AffineGeometry< TPoint >::dimension, and DGtal::AffineGeometry< TPoint >::transform().

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

◆ affineBasis() [2/2]

template<typename TPoint >
template<typename TInputPoint , typename TIndexRange >
static std::pair< TInputPoint, Points > DGtal::AffineGeometry< TPoint >::affineBasis ( const std::vector< TInputPoint > &  X,
const TIndexRange &  I,
const double  tolerance = 1e-12 
)
inlinestatic

Given a range of points X and a range of indices I describing its subset of interest, returns a point and a range of vectors forming an affine basis containing the subset X[I].

Template Parameters
TInputPointthe type of input points (may be less precise).
Parameters
[in]Xthe range of input points (may be lattice points or not).
[in]Ithe range of indices within X that specifies the subset of interest.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
a point and a range of vectors forming an affine basis containing X[I].
Note
Complexity is O( m n^2 ), where m=Cardinal(I) and n=dimension.

Definition at line 627 of file AffineGeometry.h.

630 {
631 Points basis; //< direction vectors
632 Size m = I.size();
633 // Process trivial cases.
634 if ( m == 0 ) return std::make_pair( TInputPoint::zero, basis );
635 if ( m == 1 ) return std::make_pair( X[ I[ 0 ] ], basis );
636 // Process general case.
637 const Point o = transform( X[ I[ 0 ] ] );
638 basis.reserve( dimension );
639 for ( Size i = 1; i < m; i++ )
640 {
641 Point v = transform( X[ I[ i ] ] ) - o;
642 if ( addIfIndependent( basis, v, tolerance ) )
643 if ( basis.size() > dimension ) break;
644 }
645 return std::make_pair( X[ I[ 0 ] ], basis );
646 }

References DGtal::AffineGeometry< TPoint >::addIfIndependent(), DGtal::AffineGeometry< TPoint >::dimension, and DGtal::AffineGeometry< TPoint >::transform().

◆ affineDimension()

template<typename TPoint >
template<typename TInputPoint >
static DGtal::int64_t DGtal::AffineGeometry< TPoint >::affineDimension ( const std::vector< TInputPoint > &  X,
const double  tolerance = 1e-12 
)
inlinestatic

Given a range of points X, returns the affine dimension of its spanned affine subspace.

Template Parameters
TInputPointthe type of input points (may be less precise).
Parameters
[in]Xthe range of input points (may be lattice points or not).
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
the affine dimension of X (-1 is empty set, 0 is a point, 1 is a line, etc)

Definition at line 486 of file AffineGeometry.h.

487 {
488 return DGtal::int64_t( affineSubset( X, tolerance ).size() ) - 1;
489 }
std::int64_t int64_t
signed 94-bit integer.
Definition BasicTypes.h:73
static std::vector< Size > affineSubset(const std::vector< TInputPoint > &X, const double tolerance=1e-12)

References DGtal::AffineGeometry< TPoint >::affineSubset().

Referenced by DGtal::functions::computeAffineDimension(), and SCENARIO().

◆ affineSubset() [1/2]

template<typename TPoint >
template<typename TInputPoint >
static std::vector< Size > DGtal::AffineGeometry< TPoint >::affineSubset ( const std::vector< TInputPoint > &  X,
const double  tolerance = 1e-12 
)
inlinestatic

Given a range of points X, returns a subset of these points that form an affine basis of X. Equivalently it is a simplex whose affine space spans all the points of X.

Template Parameters
TInputPointthe type of input points (may be less precise).
Parameters
[in]Xthe range of input points (may be lattice points or not).
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
a subset of these points as a range of indices.
Note
Complexity is \(O( m n^2 )\), where m=Cardinal(X) and n=dimension.

Definition at line 508 of file AffineGeometry.h.

509 {
510 Size m = X.size();
511 // Process trivial cases.
512 if ( m == 0 ) return { };
513 if ( m == 1 ) return { Size( 0 ) };
514 // Process general case.
515 Points basis; //< direction vectors
516 Sizes chosen; //< selected points
517 chosen.push_back( 0 ); //< reference point (first one, as it may be any one)
518 for ( Size i = 1; i < m; i++ )
519 {
520 Point v = transform( X[ i ] - X[ 0 ] );
521 if ( addIfIndependent( basis, v, tolerance ) )
522 chosen.push_back( i );
523 if ( chosen.size() > dimension ) break;
524 }
525 return chosen;
526 }
std::vector< Size > Sizes
type for range of sizes.

References DGtal::AffineGeometry< TPoint >::addIfIndependent(), DGtal::AffineGeometry< TPoint >::dimension, and DGtal::AffineGeometry< TPoint >::transform().

Referenced by DGtal::AffineGeometry< TPoint >::affineDimension(), DGtal::functions::computeAffineSubset(), DGtal::functions::computeAffineSubset(), DGtal::functions::getOrthogonalVector(), and DGtal::QuickHull< TKernel >::pickInitialSimplex().

◆ affineSubset() [2/2]

template<typename TPoint >
template<typename TInputPoint , typename TIndexRange >
static std::vector< Size > DGtal::AffineGeometry< TPoint >::affineSubset ( const std::vector< TInputPoint > &  X,
const TIndexRange &  I,
const double  tolerance = 1e-12 
)
inlinestatic

Given a range of points X and a subset of it given by indices I, returns a subset of these points that form an affine basis of X[I]. Equivalently it is a simplex whose affine space spans all the points of X[I].

Template Parameters
TInputPointthe type of input points (may be less precise).
TIndexRangeany type of range of indices that specifies the subset of points of interest.
Parameters
[in]Xthe range of input points (may be lattice points or not).
[in]Ithe range of indices specifying the subset of interest.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
a subset of these points as a range of indices.
Note
Complexity is \(O( m n^2 )\), where m=Cardinal(X) and n=dimension.

Definition at line 551 of file AffineGeometry.h.

554 {
555 Size m = I.size();
556 // Process trivial cases.
557 if ( m == 0 ) return { };
558 if ( m == 1 ) return { I[ 0 ] };
559 // Process general case.
560 Points basis; //< direction vectors
561 Sizes chosen; //< selected points
562 chosen.push_back( I[ 0 ] ); //< reference point (first one, as it may be any one)
563 for ( Size i = 1; i < m; i++ )
564 {
565 Point v = transform( X[ I[ i ] ] - X[ I[ 0 ] ] );
566 if ( addIfIndependent( basis, v, tolerance ) )
567 chosen.push_back( I[ i ] );
568 if ( chosen.size() > dimension ) break;
569 }
570 return chosen;
571 }

References DGtal::AffineGeometry< TPoint >::addIfIndependent(), DGtal::AffineGeometry< TPoint >::dimension, and DGtal::AffineGeometry< TPoint >::transform().

◆ completeBasis()

template<typename TPoint >
static void DGtal::AffineGeometry< TPoint >::completeBasis ( Points basis,
bool  normal_vector,
const double  tolerance = 1e-12 
)
inlinestatic

Complete the vectors basis with independent vectors so as to form a basis of the space. When normal_vector is true, the last added vector is guaranteed to be orthogonal to all the previous vectors.

Note
In 3D, given two independent vectors as input, then the added vector is the cross product of these two vectors. In nD, it is thus a generalization of the cross product.
Parameters
[in,out]basisa range of independent vectors of size less than dimension.
[in]normal_vectorwhen 'true', the last vector of the basis is orthogonal to all the other vectors of the basis, otherwise it is just an independent canonic vector.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).

Definition at line 737 of file AffineGeometry.h.

740 {
741 if ( basis.size() >= dimension ) return;
742 while ( ( basis.size() + 1 ) < dimension )
743 { // add an independent vector
744 const Point u = independentVector( basis, tolerance );
745 basis.push_back( u );
746 }
747 if ( normal_vector )
748 basis.push_back( orthogonalVector( basis ) );
749 else
750 basis.push_back( independentVector( basis, tolerance ) );
751 }
static Point independentVector(const Points &basis, const double tolerance=1e-12)
static Point orthogonalVector(const Points &basis)

References DGtal::AffineGeometry< TPoint >::dimension, DGtal::AffineGeometry< TPoint >::independentVector(), and DGtal::AffineGeometry< TPoint >::orthogonalVector().

Referenced by DGtal::functions::getCompleteBasis().

◆ independentVector() [1/2]

template<typename TPoint >
static Point DGtal::AffineGeometry< TPoint >::independentVector ( const Points basis,
const double  tolerance = 1e-12 
)
inlinestatic

Given a partial basis of vectors, returns a new vector that is independent.

Parameters
[in]basisa range of independent vectors that defines a partial basis of the space.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
a canonic unit vector independent of all vectors of basis, or the null vector if the basis was not partial.

Definition at line 660 of file AffineGeometry.h.

661 {
662 // If basis has already d independent vectors, then there is no
663 // other independent vector.
664 if ( basis.size() >= dimension ) return Point::zero;
665 // At least one trivial canonic vector should be independant.
666 Dimension k = 0;
667 for ( ; k < dimension; k++ )
668 {
669 Point e_k = Point::base( k );
670 Point w = reductionOnBasis( e_k, basis, tolerance );
671 Scalar sql2;
673 if ( ScalarOps::isNonZero( sql2, tolerance ) )
674 return e_k;
675 }
676 trace.error() << "[AffineGeometry::independentVector]"
677 << " Unable to find independent vector." << std::endl;
678 return Point::zero;
679 }
std::ostream & error()
DGtal::uint32_t Dimension
Definition Common.h:119
Trace trace

References DGtal::AffineGeometry< TPoint >::dimension, DGtal::Trace::error(), DGtal::functions::getSquaredNormL2(), DGtal::detail::AffineGeometryScalarOperations< TScalar >::isNonZero(), DGtal::AffineGeometry< TPoint >::reductionOnBasis(), and DGtal::trace.

Referenced by DGtal::AffineGeometry< TPoint >::completeBasis(), and DGtal::functions::computeIndependentVector().

◆ independentVector() [2/2]

template<typename TPoint >
template<typename TOtherPoint >
static TOtherPoint DGtal::AffineGeometry< TPoint >::independentVector ( const Points basis,
const double  tolerance = 1e-12 
)
inlinestatic

Given a partial basis of vectors, returns a new vector that is independent.

Template Parameters
TOtherPointany type of point
Parameters
[in]basisa range of independent vectors that defines a partial basis of the space.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
a canonic unit vector independent of all vectors of basis, or the null vector if the basis was not partial.

Definition at line 696 of file AffineGeometry.h.

697 {
698 // If basis has already d independent vectors, then there is no
699 // other independent vector.
700 if ( basis.size() >= dimension ) return TOtherPoint::zero;
701 // At least one trivial canonic vector should be independant.
702 Dimension k = 0;
703 for ( ; k < dimension; k++ )
704 {
705 Point e_k = Point::base( k );
706 Point w = reductionOnBasis( e_k, basis, tolerance );
707 Scalar sql2;
709 if ( ScalarOps::isNonZero( sql2, tolerance ) )
710 return TOtherPoint::base( k );
711 }
712 trace.error() << "[AffineGeometry::independentVector]"
713 << " Unable to find independent vector." << std::endl;
714 return TOtherPoint::zero;
715 }

References DGtal::AffineGeometry< TPoint >::dimension, DGtal::Trace::error(), DGtal::functions::getSquaredNormL2(), DGtal::detail::AffineGeometryScalarOperations< TScalar >::isNonZero(), DGtal::AffineGeometry< TPoint >::reductionOnBasis(), and DGtal::trace.

◆ orthogonalLatticeBasis()

template<typename TPoint >
template<typename TInputPoint >
static std::vector< Point > DGtal::AffineGeometry< TPoint >::orthogonalLatticeBasis ( const TInputPoint &  N,
bool  shortened = false 
)
inlinestatic

Given a primitive lattice vector N, returns a possible basis for its orthogonal d-1 dimensional lattice.

Template Parameters
TOtherPointany type of lattice point
Parameters
[in]Na non null primitive lattice vector
[in]shortenedwhen 'false', the basis is in row echelon form, when 'true', the basis is shortened (with pairwise shortenings, may not be optimal, except in 3D).
Returns
the d-1 vectors forming a basis of the orthogonal lattice to N.

Definition at line 768 of file AffineGeometry.h.

769 {
770 std::vector< Scalar > Np( N.begin(), N.end() );
771 Np.resize( Point::dimension );
773 if ( shortened )
775 std::vector< Point > R( B.size() );
776 for ( std::size_t i = 0; i < B.size(); i++ )
777 for ( std::size_t j = 0; j < Point::dimension; j++ )
778 R[ i ][ j ] = B[ i ][ j ];
779 return R;
780 }
std::size_t shortenBasis(std::vector< std::vector< TComponent > > &B)
std::vector< std::vector< TComponent > > computeOrthogonalLattice(std::vector< TComponent > N)

References DGtal::functions::computeOrthogonalLattice(), DGtal::R, and DGtal::functions::shortenBasis().

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

◆ orthogonalVector() [1/2]

template<typename TPoint >
static Point DGtal::AffineGeometry< TPoint >::orthogonalVector ( const Points basis)
inlinestatic

Given d-1 independent vectors in dD, returns a vector that is orthogonal to each of them.

Note
In 3D, given two independent vectors as input, then the added vector is the (reduced) cross product of these two vectors. In nD, it is thus a kind of generalization of the cross product.
Parameters
[in]basisa range of independent vectors of size dimension-1.
Returns
a vector of coefficients (represented with the given number type), or the null vector if the basis is not d-1-dimensional.

Definition at line 962 of file AffineGeometry.h.

963 {
964 Point w;
965 const std::size_t n = dimension;
966 if ( ( basis.size() + 1 ) != dimension ) return w;
967 const std::size_t m = basis.size();
968 SimpleMatrix< Scalar, dimension-1, dimension> A;
969 for ( std::size_t i = 0; i < m; ++i )
970 for ( std::size_t j = 0; j < n; ++j )
971 A( i, j ) = basis[ i ][ j ];
972 for ( std::size_t col = 0; col < n; ++col)
973 { // construct sub-matrix removing column col
974 SimpleMatrix< Scalar, dimension-1, dimension-1> M;
975 for ( std::size_t i = 0; i < m; ++i )
976 {
977 std::size_t c = 0;
978 for ( std::size_t j = 0; j < n; ++j)
979 if ( j != col ) M( i, c++ ) = A( i, j );
980 }
982 if ( (col+dimension) % 2 == 0 ) w[ col ] = -w[ col ];
983 }
984 return simplifiedVector( w );
985 }
void getDeterminantBareiss(TInternalNumber &result, const SimpleMatrix< TComponent, TN, TN > &matrix)
Point::Coordinate Scalar
static Point simplifiedVector(Point v)

References DGtal::AffineGeometry< TPoint >::dimension, DGtal::functions::getDeterminantBareiss(), and DGtal::AffineGeometry< TPoint >::simplifiedVector().

Referenced by DGtal::AffineGeometry< TPoint >::completeBasis(), and DGtal::functions::computeOrthogonalVectorToBasis().

◆ orthogonalVector() [2/2]

template<typename TPoint >
template<typename TInternalNumber >
static Point DGtal::AffineGeometry< TPoint >::orthogonalVector ( const Points basis)
inlinestatic

Given d-1 independent vectors in dD, returns a vector that is orthogonal to each of them.

Note
In 3D, given two independent vectors as input, then the added vector is the (reduced) cross product of these two vectors. In nD, it is thus a kind of generalization of the cross product.
Parameters
[in]basisa range of independent vectors of size dimension-1.
Returns
a vector of coefficients (represented with the given number type), or the null vector if the basis is not d-1-dimensional.

Definition at line 1002 of file AffineGeometry.h.

1003 {
1004 Point w;
1005 typedef PointVector< dimension, TInternalNumber > InternalPoint;
1006 InternalPoint W;
1007 const std::size_t n = dimension;
1008 if ( ( basis.size() + 1 ) != dimension ) return w;
1009 const std::size_t m = basis.size();
1010 SimpleMatrix< Scalar, dimension-1, dimension> A;
1011 for ( std::size_t i = 0; i < m; ++i )
1012 for ( std::size_t j = 0; j < n; ++j )
1013 A( i, j ) = basis[ i ][ j ];
1014 for ( std::size_t col = 0; col < n; ++col)
1015 { // construct sub-matrix removing column col
1016 SimpleMatrix< Scalar, dimension-1, dimension-1> M;
1017 for ( std::size_t i = 0; i < m; ++i )
1018 {
1019 std::size_t c = 0;
1020 for ( std::size_t j = 0; j < n; ++j)
1021 if ( j != col ) M( i, c++ ) = A( i, j );
1022 }
1023 functions::getDeterminantBareiss( W[ col ], M );
1024 if ( (col+dimension) % 2 == 0 ) W[ col ] = -W[ col ];
1025 }
1027 return transform( W );
1028 }

References DGtal::AffineGeometry< TPoint >::dimension, DGtal::functions::getDeterminantBareiss(), DGtal::AffineGeometry< TPoint >::simplifiedVector(), DGtal::PointVector< dim, TEuclideanRing, TContainer >::size(), and DGtal::AffineGeometry< TPoint >::transform().

◆ reduceVector() [1/2]

template<typename TPoint >
static std::pair< Scalar, Scalar > DGtal::AffineGeometry< TPoint >::reduceVector ( Point w,
const Point b,
const double  tolerance 
)
inlinestatic

Reduces vector w by the vector b

Parameters
[in,out]wany vector
[in]ba non-null vector of the current basis.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
the coefficients (alpha,beta) for reduction such that alpha * w - beta * b is the returned reduced vector, alpha >= 0.
Note
This reduction is an elementary Gauss elimination.

Definition at line 867 of file AffineGeometry.h.

868 {
869 Scalar mul_w = 0;
870 Scalar mul_b = 0;
871 Size n = w.size();
872 // Find index of first non null pivot in b.
873 Size lead = n;
874 for ( Size j = 0; j < n; j++)
875 if ( ScalarOps::isNonZero( b[j], tolerance ) ) { lead = j; break; }
876 if ( lead == n ) return std::make_pair( mul_w, mul_b ); // b is null vector
877
878 std::tie( mul_w, mul_b ) = ScalarOps::getMultipliers( w[ lead ], b[ lead ] );
879
880 for (Size j = 0; j < n; j++)
881 w[j] = mul_w * w[j] - mul_b * b[j];
882 return std::make_pair( mul_w, mul_b );
883 }
static std::pair< Integer, Integer > getMultipliers(Integer a, Integer b)

References DGtal::detail::AffineGeometryScalarOperations< TScalar >::getMultipliers(), and DGtal::detail::AffineGeometryScalarOperations< TScalar >::isNonZero().

Referenced by DGtal::AffineBasis< TPoint >::decomposeVector(), DGtal::AffineBasis< TPoint >::reduceAsEchelon(), and DGtal::AffineGeometry< TPoint >::reductionOnBasis().

◆ reduceVector() [2/2]

template<typename TPoint >
static std::pair< Scalar, Scalar > DGtal::AffineGeometry< TPoint >::reduceVector ( Point w,
const Point b,
Dimension  start,
const double  tolerance 
)
inlinestatic

Reduces vector w by the vector b

Parameters
[in,out]wany vector
[in]ba non-null vector of the current basis.
[in]startthe starting component index.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
the coefficients (alpha,beta) for reduction such that alpha * w - beta * b is the returned reduced vector, alpha >= 0.
Note
This reduction is an elementary Gauss elimination.

Definition at line 900 of file AffineGeometry.h.

902 {
903 Scalar mul_w = 0;
904 Scalar mul_b = 0;
905 Size n = w.size();
906 // Find index of first non null pivot in b.
907 Size lead = n;
908 for ( Size j = start; j < n; j++)
909 if ( ScalarOps::isNonZero( b[j], tolerance ) ) { lead = j; break; }
910 if ( lead == n ) return std::make_pair( mul_w, mul_b ); // b is null vector
911
912 std::tie( mul_w, mul_b ) = ScalarOps::getMultipliers( w[ lead ], b[ lead ] );
913
914 for (Size j = 0; j < n; j++)
915 w[j] = mul_w * w[j] - mul_b * b[j];
916 return std::make_pair( mul_w, mul_b );
917 }

References DGtal::detail::AffineGeometryScalarOperations< TScalar >::getMultipliers(), and DGtal::detail::AffineGeometryScalarOperations< TScalar >::isNonZero().

◆ reductionOnBasis()

template<typename TPoint >
static Point DGtal::AffineGeometry< TPoint >::reductionOnBasis ( const Point v,
const Points basis,
const double  tolerance 
)
inlinestatic

Reduces the vector v on the (partial or not) basis of vectors basis.

Parameters
[in]vany vector
[in]basisa range of vectors forming a (partial or not) basis of the space.
[in]tolerancethe accepted oo-norm below which the vector is null (used only for points with float/double coordinates).
Returns
the part of v that cannot be expressed as a linear combination of vectors of basis, and hence a null vector if v is a linear combination of the vectors of the basis.

Definition at line 804 of file AffineGeometry.h.

807 {
808 Point w( v );
809 for ( const auto& b : basis ) reduceVector( w, b, tolerance );
810 return w;
811 }
static std::pair< Scalar, Scalar > reduceVector(Point &w, const Point &b, const double tolerance)

References DGtal::AffineGeometry< TPoint >::reduceVector().

Referenced by DGtal::AffineGeometry< TPoint >::addIfIndependent(), and DGtal::AffineGeometry< TPoint >::independentVector().

◆ simplifiedVector()

template<typename TPoint >
static Point DGtal::AffineGeometry< TPoint >::simplifiedVector ( Point  v)
inlinestatic

Given a vector, returns the aligned vector with its component simplified by the gcd of all components (when the components are integers) or the aligned vector with unit L2-norm (when the components are floating-point numbers).

Parameters
[in]vany vector.
Returns
a simplified vector aligned with v.

Definition at line 1039 of file AffineGeometry.h.

1040 {
1041 Scalar x;
1044 return v;
1045 }

References DGtal::functions::getSquaredNormL2(), and DGtal::detail::AffineGeometryPointOperations< dim, TEuclideanRing, TContainer >::normalizeVector().

Referenced by DGtal::functions::computeSimplifiedVector(), DGtal::AffineBasis< TPoint >::normalize(), DGtal::AffineGeometry< TPoint >::orthogonalVector(), DGtal::AffineBasis< TPoint >::reduceAsEchelon(), and DGtal::AffineBasis< TPoint >::reduceAsLLL().

◆ transform() [1/2]

template<typename TPoint >
static const Point & DGtal::AffineGeometry< TPoint >::transform ( const Point w)
inlinestatic

Transform a vector or point into the representation chosen for AffineGeometry class. This overloading takes care of the trivial case, where there is no need to change the representation.

Parameters
[in]wany vector or point
Returns
the same vector or point (does nothing).

Definition at line 929 of file AffineGeometry.h.

930 {
931 return w;
932 }

Referenced by DGtal::AffineBasis< TPoint >::AffineBasis(), DGtal::AffineGeometry< TPoint >::affineBasis(), DGtal::AffineGeometry< TPoint >::affineBasis(), DGtal::AffineBasis< TPoint >::AffineBasis(), DGtal::AffineBasis< TPoint >::AffineBasis(), DGtal::AffineGeometry< TPoint >::affineSubset(), DGtal::AffineGeometry< TPoint >::affineSubset(), DGtal::AffineBasis< TPoint >::initBasis(), and DGtal::AffineGeometry< TPoint >::orthogonalVector().

◆ transform() [2/2]

template<typename TPoint >
template<typename TInputPoint >
static Point DGtal::AffineGeometry< TPoint >::transform ( const TInputPoint &  w)
inlinestatic

Transform a vector or point into the representation chosen for AffineGeometry class. This overloading takes care of the generic case, where a transformation of each component is required.

Parameters
[in]wany vector or point
Returns
the same vector or point, but in the represention chosen for this class.

Definition at line 943 of file AffineGeometry.h.

944 {
945 return PointOps::cast( w );
946 }
static Point cast(const TOtherPoint &other)
Specialized version to cast points to other point type.

References DGtal::detail::AffineGeometryPointOperations< dim, TEuclideanRing, TContainer >::cast().

Field Documentation

◆ dimension


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