DGtal  0.6.devel
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | Static Public Member Functions | Private Attributes
DGtal::IntegerComputer< TInteger > Class Template Reference

#include <IntegerComputer.h>

Inheritance diagram for DGtal::IntegerComputer< TInteger >:
Inheritance graph
[legend]
Collaboration diagram for DGtal::IntegerComputer< TInteger >:
Collaboration graph
[legend]

Public Types

typedef IntegerComputer< TInteger > Self
typedef NumberTraits< TInteger >
::SignedVersion 
Integer
typedef NumberTraits< Integer >
::ParamType 
IntegerParamType
typedef NumberTraits< TInteger >
::UnsignedVersion 
UnsignedInteger
typedef NumberTraits
< UnsignedInteger >::ParamType 
UnsignedIntegerParamType
typedef SpaceND< 2, Integer >
::Point 
Point2I
typedef SpaceND< 2, Integer >
::Vector 
Vector2I
typedef SpaceND< 3, Integer >
::Point 
Point3I
typedef SpaceND< 3, Integer >
::Vector 
Vector3I

Public Member Functions

 BOOST_CONCEPT_ASSERT ((CInteger< Integer >))
 BOOST_CONCEPT_ASSERT ((CUnsignedInteger< UnsignedInteger >))
 ~IntegerComputer ()
 IntegerComputer ()
 IntegerComputer (const Self &other)
Selfoperator= (const Self &other)
void getEuclideanDiv (Integer &q, Integer &r, IntegerParamType a, IntegerParamType b) const
Integer floorDiv (IntegerParamType na, IntegerParamType nb) const
Integer ceilDiv (IntegerParamType na, IntegerParamType nb) const
void getFloorCeilDiv (Integer &fl, Integer &ce, IntegerParamType na, IntegerParamType nb) const
Integer gcd (IntegerParamType a, IntegerParamType b) const
void getGcd (Integer &g, IntegerParamType a, IntegerParamType b) const
Integer getCFrac (std::vector< Integer > &quotients, IntegerParamType a, IntegerParamType b) const
template<typename OutputIterator >
Integer getCFrac (OutputIterator outIt, IntegerParamType a, IntegerParamType b) const
Point2I convergent (const std::vector< Integer > &quotients, unsigned int k) const
void reduce (Vector2I &p) const
Integer crossProduct (const Vector2I &u, const Vector2I &v) const
void getCrossProduct (Integer &cp, const Vector2I &u, const Vector2I &v) const
Integer dotProduct (const Vector2I &u, const Vector2I &v) const
void getDotProduct (Integer &dp, const Vector2I &u, const Vector2I &v) const
Vector2I extendedEuclid (IntegerParamType a, IntegerParamType b, IntegerParamType c) const
void getCoefficientIntersection (Integer &fl, Integer &ce, const Vector2I &p, const Vector2I &u, const Vector2I &N, IntegerParamType c) const
void getValidBezout (Vector2I &v, const Point2I &A, const Vector2I &u, const Vector2I &N, IntegerParamType c, const Vector2I &N2, IntegerParamType c2, bool compute_v=true) const
void reduce (Vector3I &p) const
Integer dotProduct (const Vector3I &u, const Vector3I &v) const
void getDotProduct (Integer &dp, const Vector3I &u, const Vector3I &v) const
void selfDisplay (std::ostream &out) const
bool isValid () const

Static Public Member Functions

static bool isZero (IntegerParamType a)
static bool isNotZero (IntegerParamType a)
static bool isPositive (IntegerParamType a)
static bool isNegative (IntegerParamType a)
static bool isPositiveOrZero (IntegerParamType a)
static bool isNegativeOrZero (IntegerParamType a)
static Integer abs (IntegerParamType a)
static Integer max (IntegerParamType a, IntegerParamType b)
static Integer max (IntegerParamType a, IntegerParamType b, IntegerParamType c)
static Integer min (IntegerParamType a, IntegerParamType b)
static Integer min (IntegerParamType a, IntegerParamType b, IntegerParamType c)
static Integer staticGcd (IntegerParamType a, IntegerParamType b)

Private Attributes

Integer _m_a
Integer _m_b
Integer _m_a0
Integer _m_a1
Integer _m_q
Integer _m_r
std::vector< Integer_m_bezout [4]
Vector2I _m_v
Vector2I _m_v0
Vector2I _m_v1
Integer _m_c0
Integer _m_c1
Integer _m_c2

Detailed Description

template<typename TInteger>
class DGtal::IntegerComputer< TInteger >

Aim: This class gathers several types and methods to make computation with integers.

Description of template class 'IntegerComputer'

This class is especially useful with using big integers (like GMP), since a substantial part of the execution time cames from the allocation/desallocation of integers. The idea is that the user instantiate once this object and computes gcd, bezout, continued fractions with it.

To be thread-safe, each thread must instantiate an IntegerComputer.

It is a model of boost::CopyConstructible, boost::DefaultConstructible, boost::Assignable. All its member data are mutable.

It is a backport of ImaGene.

Template Parameters:
TIntegerany model of integer (CInteger), like int, long int, int64_t, BigInteger (when GMP is installed).

Definition at line 82 of file IntegerComputer.h.


Member Typedef Documentation

template<typename TInteger>
typedef NumberTraits<TInteger>::SignedVersion DGtal::IntegerComputer< TInteger >::Integer

Definition at line 87 of file IntegerComputer.h.

template<typename TInteger>
typedef NumberTraits<Integer>::ParamType DGtal::IntegerComputer< TInteger >::IntegerParamType

Definition at line 88 of file IntegerComputer.h.

template<typename TInteger>
typedef SpaceND<2,Integer>::Point DGtal::IntegerComputer< TInteger >::Point2I

Definition at line 93 of file IntegerComputer.h.

template<typename TInteger>
typedef SpaceND<3,Integer>::Point DGtal::IntegerComputer< TInteger >::Point3I

Definition at line 95 of file IntegerComputer.h.

template<typename TInteger>
typedef IntegerComputer<TInteger> DGtal::IntegerComputer< TInteger >::Self

Definition at line 86 of file IntegerComputer.h.

template<typename TInteger>
typedef NumberTraits<TInteger>::UnsignedVersion DGtal::IntegerComputer< TInteger >::UnsignedInteger

Definition at line 90 of file IntegerComputer.h.

template<typename TInteger>
typedef NumberTraits<UnsignedInteger>::ParamType DGtal::IntegerComputer< TInteger >::UnsignedIntegerParamType

Definition at line 91 of file IntegerComputer.h.

template<typename TInteger>
typedef SpaceND<2,Integer>::Vector DGtal::IntegerComputer< TInteger >::Vector2I

Definition at line 94 of file IntegerComputer.h.

template<typename TInteger>
typedef SpaceND<3,Integer>::Vector DGtal::IntegerComputer< TInteger >::Vector3I

Definition at line 96 of file IntegerComputer.h.


Constructor & Destructor Documentation

template<typename TInteger >
DGtal::IntegerComputer< TInteger >::~IntegerComputer ( )
inline

Destructor.

Definition at line 44 of file IntegerComputer.ih.

{}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::IntegerComputer ( )
inline

Constructor. Each thread must have its own instance for all computations. Such object stores several local variables to limit the number of memory allocations.

Does nothing (member data are allocated, but their values are not used except during the execution of methods).

Definition at line 49 of file IntegerComputer.ih.

{}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::IntegerComputer ( const Self other)
inline

Copy constructor.

Does nothing (member data are allocated, but their values are not used except during the execution of methods).

Parameters:
otherthe object to clone.

Definition at line 54 of file IntegerComputer.ih.

{}

Member Function Documentation

template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::abs ( IntegerParamType  a)
inlinestatic
Parameters:
aany integer.
Returns:
its absolute value.

Definition at line 125 of file IntegerComputer.ih.

{
if ( isPositiveOrZero( a ) )
return a;
else
return -a;
}
template<typename TInteger>
DGtal::IntegerComputer< TInteger >::BOOST_CONCEPT_ASSERT ( (CInteger< Integer >)  )
template<typename TInteger>
DGtal::IntegerComputer< TInteger >::BOOST_CONCEPT_ASSERT ( (CUnsignedInteger< UnsignedInteger >)  )
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::ceilDiv ( IntegerParamType  na,
IntegerParamType  nb 
) const
inline

Computes the ceil value of na/nb.

Parameters:
naany integer.
nbany integer.

Definition at line 213 of file IntegerComputer.ih.

{
_m_a = na;
_m_b = nb;
if ( isNegative( _m_b ) )
{
}
// if ( isNegative( _m_a ) && isNegative( _m_b ) )
// {
// _m_a=-_m_a;
// _m_b=-_m_b;
// }
// else if ( isNegative( _m_b ) )
// {
// _m_a=-_m_a;
// _m_b=-_m_b;
// }
if ( isNegative( _m_a ) || isZero( _m_a % _m_b ) )
return _m_a/_m_b;
else
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Point2I DGtal::IntegerComputer< TInteger >::convergent ( const std::vector< Integer > &  quotients,
unsigned int  k 
) const
inline

Returns the fraction corresponding to the given quotients, more precisely its k-th principal convergent. When k >= quotients.size() - 1, it is the inverse of the function getCFrac.

Parameters:
quotientsthe sequence of partial quotients.
kthe desired partial convergent.
Returns:
the corresponding fraction p_k / q_k as a point (p_k, q_k).

Definition at line 365 of file IntegerComputer.ih.

{
_m_bezout[ 0 ].clear(); // p
_m_bezout[ 1 ].clear(); // q
if ( k >= quotients.size() )
k = (quotients.size() - 1);
for ( unsigned int i = 0; i <= k; ++i )
{
_m_bezout[ 0 ].push_back( quotients[ i ] * _m_bezout[ 0 ][ i + 1 ]
+ _m_bezout[ 0 ][ i ] );
_m_bezout[ 1 ].push_back( quotients[ i ] * _m_bezout[ 1 ][ i + 1 ]
+ _m_bezout[ 1 ][ i ] );
}
_m_v[ 0 ] = _m_bezout[ 0 ].back();
_m_v[ 1 ] = _m_bezout[ 1 ].back();
return _m_v;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::crossProduct ( const Vector2I u,
const Vector2I v 
) const
inline

Computes and returns the cross product of u and v.

Parameters:
uany vector in Z2.
vany vector in Z2.
Returns:
the cross product of u and v.

Definition at line 405 of file IntegerComputer.ih.

Referenced by DGtal::FrechetShortcut< TIterator, TInteger >::Backpath::updateIntervals().

{
_m_a0 = u[ 0 ] * v[ 1 ];
_m_a1 = u[ 1 ] * v[ 0 ];
return _m_a0 - _m_a1;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::dotProduct ( const Vector2I u,
const Vector2I v 
) const
inline

Computes and returns the dot product of u and v.

Parameters:
uany vector in Z2.
vany vector in Z2.
Returns:
the dot product of u and v.

Definition at line 427 of file IntegerComputer.ih.

Referenced by DGtal::FrechetShortcut< TIterator, TInteger >::Tools::angleVectVect(), DGtal::FrechetShortcut< TIterator, TInteger >::Backpath::updateIntervals(), and DGtal::FrechetShortcut< TIterator, TInteger >::Backpath::updateOcculters().

{
_m_a0 = u[ 0 ] * v[ 0 ];
_m_a1 = u[ 1 ] * v[ 1 ];
return _m_a0 + _m_a1;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::dotProduct ( const Vector3I u,
const Vector3I v 
) const
inline

Computes and returns the dot product of u and v.

Parameters:
uany vector in Z2.
vany vector in Z2.
Returns:
the dot product of u and v.

Definition at line 554 of file IntegerComputer.ih.

{
_m_a0 = u[ 0 ] * v[ 0 ];
_m_a0 += u[ 1 ] * v[ 1 ];
_m_a0 += u[ 2 ] * v[ 2 ];
return _m_a0;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Vector2I DGtal::IntegerComputer< TInteger >::extendedEuclid ( IntegerParamType  a,
IntegerParamType  b,
IntegerParamType  c 
) const

Returns a solution of the Diophantine equation: a x + b y = c. The integer c should be a multiple of the gcd of a and b. Uses the extended Euclid algorithm to do it, whose complexity is bounded by max(log(a),log(b)).

Parameters:
aany integer.
bany integer.
cany integer (see above).
Returns:
the vector (x,y) solution to a x + b y = c.

Definition at line 448 of file IntegerComputer.ih.

References DGtal::abs().

{
if( isZero( a ) ) return Point2I( NumberTraits<Integer>::ZERO, b * c );
if( isZero( b ) ) return Point2I( a * c, NumberTraits<Integer>::ZERO );
for ( unsigned int i = 0; i < 4; ++i )
_m_bezout[ i ].clear();
_m_bezout[ 0 ].push_back( abs( a ) );
_m_bezout[ 0 ].push_back( abs( b ) );
unsigned int k = 0; // index of the iteration during the computation.
while( isNotZero( _m_bezout[ 0 ][ k+1 ] ) )
{
_m_bezout[ 1 ].push_back( _m_bezout[ 0 ][ k ]
/ _m_bezout[ 0 ][ k+1 ] );
_m_bezout[ 0 ].push_back( _m_bezout[ 0 ][ k ]
% _m_bezout[ 0 ][ k+1 ] );
_m_bezout[ 2 ].push_back( _m_bezout[ 2 ][ k ]
- _m_bezout[ 1 ][ k+1 ]
* _m_bezout[ 2 ][ k+1 ] );
_m_bezout[ 3 ].push_back( _m_bezout[ 3 ][ k ]
- _m_bezout[ 1 ][ k+1 ]
*_m_bezout[ 3 ][ k+1 ] );
++k;
}
_m_v[ 0 ] = _m_bezout[ 2 ][ k ];
_m_v[ 1 ] = _m_bezout[ 3 ][ k ];
if ( isNegative( a ) ) _m_v[ 0 ] = -_m_v[ 0 ];
if ( isNegative( b ) ) _m_v[ 1 ] = -_m_v[ 1 ];
//_m_v *= c;
return _m_v;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::floorDiv ( IntegerParamType  na,
IntegerParamType  nb 
) const
inline

Computes the floor value of na/nb.

Parameters:
naany integer.
nbany integer.

Definition at line 184 of file IntegerComputer.ih.

Referenced by DGtal::ClosedIntegerHalfPlane< TSpace >::ClosedIntegerHalfPlane().

{
_m_a = na;
_m_b = nb;
if ( isNegative( _m_b ) )
{
}
// if ( isNegative( _m_a ) && isNegative( _m_b ) )
// {
// _m_a=-_m_a;
// _m_b=-_m_b;
// }
// else if ( isNegative( _m_b ) )
// {
// _m_a=-_m_a;
// _m_b=-_m_b;
// }
if ( isPositive( _m_a ) || isZero( _m_a % _m_b ) )
return _m_a/_m_b;
else
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::gcd ( IntegerParamType  a,
IntegerParamType  b 
) const
inline

Returns the greatest common divisor of a and b (a and b may be either positive or negative).

Parameters:
aany integer.
bany integer.
Returns:
the gcd of a and b.

Definition at line 284 of file IntegerComputer.ih.

References DGtal::abs().

Referenced by DGtal::ClosedIntegerHalfPlane< TSpace >::ClosedIntegerHalfPlane(), DGtal::SternBrocot< TInteger, TQuotient >::fraction(), and DGtal::Pattern< TFraction >::Pattern().

{
_m_a = abs( a );
_m_b = abs( b );
_m_a0 = max( _m_a, _m_b );
_m_a1 = min( _m_a, _m_b );
while ( isNotZero( _m_a1 ) )
{
}
return _m_a0;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::getCFrac ( std::vector< Integer > &  quotients,
IntegerParamType  a,
IntegerParamType  b 
) const
inline

Computes and push_backs the simple continued fraction of a / b. (positive) and b (positive). For instance, 5/13=[0;2,1,1,2], which is exactly what is pushed at the back of quotients.

Parameters:
quotients(modifies) adds to the back of the vector the quotients of the continued fraction of a/b.
aany positive integer.
bany positive integer.
Returns:
the gcd of a and b.

Definition at line 323 of file IntegerComputer.ih.

{
ASSERT( isPositiveOrZero( a ) && isPositiveOrZero( b ) );
_m_a0 = a;
_m_a1 = b;
while ( isNotZero( _m_a1 ) )
{
quotients.push_back( _m_q );
}
return _m_a0;
}
template<typename TInteger >
template<typename OutputIterator >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::getCFrac ( OutputIterator  outIt,
IntegerParamType  a,
IntegerParamType  b 
) const
inline

Computes and outputs the quotients of the simple continued fraction of a / b. (positive) and b (positive). For instance, 5/13=[0;2,1,1,2], which is exactly what is outputed with the OutputIterator outIt.

Template Parameters:
OutputIteratora model of boost::OutputIterator
Parameters:
outItan instance of output iterator that is used to write the successive quotients of the continued fraction of a/b.
aany positive integer.
bany positive integer.
Returns:
the gcd of a and b.

Definition at line 344 of file IntegerComputer.ih.

{
BOOST_CONCEPT_ASSERT(( boost::OutputIterator< OutputIterator, Integer > ));
ASSERT( isPositiveOrZero( a ) && isPositiveOrZero( b ) );
_m_a0 = a;
_m_a1 = b;
while ( isNotZero( _m_a1 ) )
{
*outIt++ = _m_q;
}
return _m_a0;
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getCoefficientIntersection ( Integer fl,
Integer ce,
const Vector2I p,
const Vector2I u,
const Vector2I N,
IntegerParamType  c 
) const
inline

Computes the floor (fl) and the ceiling (ce) value of the real number k such that p + k u lies on the supporting line of the linear constraint N.p <= c.

Otherwise said: (u.N) fl <= c - p.N < (u.N) ce

Parameters:
flthe greatest integer such that (u.N) fl <= c - p.N
cethe smallest integer such that c - p.N < (u.N) ce
pany vector in Z2
uany vector in Z2 in the same quadrant as N.
Nany vector in Z2 in the same quadrant as u.
cany integer.

Definition at line 494 of file IntegerComputer.ih.

{
getDotProduct( _m_c0, p, N );
getDotProduct( _m_c1, u, N );
_m_c2 = c - _m_c0;
fl = floorDiv( _m_c2, _m_c1 );
ce = ceilDiv( _m_c2, _m_c1 );
// getFloorCeilDiv( fl, ce, _m_c2, _m_c1 );
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getCrossProduct ( Integer cp,
const Vector2I u,
const Vector2I v 
) const
inline

Computes the cross product of u and v.

Parameters:
cp(returns) the cross product of u and v.
uany vector in Z2.
vany vector in Z2.

Definition at line 416 of file IntegerComputer.ih.

{
cp = u[ 0 ] * v[ 1 ];
cp -= u[ 1 ] * v[ 0 ];
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getDotProduct ( Integer dp,
const Vector2I u,
const Vector2I v 
) const
inline

Computes and returns the dot product of u and v.

Parameters:
dp(returns) the dot product of u and v.
uany vector in Z2.
vany vector in Z2.

Definition at line 438 of file IntegerComputer.ih.

Referenced by DGtal::ClosedIntegerHalfPlane< TSpace >::ClosedIntegerHalfPlane().

{
dp = u[ 0 ] * v[ 0 ];
dp += u[ 1 ] * v[ 1 ];
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getDotProduct ( Integer dp,
const Vector3I u,
const Vector3I v 
) const
inline

Computes and returns the dot product of u and v.

Parameters:
dp(returns) the dot product of u and v.
uany vector in Z3.
vany vector in Z3.

Definition at line 566 of file IntegerComputer.ih.

{
dp = u[ 0 ] * v[ 0 ];
dp += u[ 1 ] * v[ 1 ];
dp += u[ 2 ] * v[ 2 ];
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getEuclideanDiv ( Integer q,
Integer r,
IntegerParamType  a,
IntegerParamType  b 
) const
inline

Computes the euclidean division of a/b, returning quotient and remainder. May be faster than computing separately quotient and remainder, depending on the integral type in use.

Todo:
Specialize it for GMP as mpz_fdiv_qr (q.get_mpz_t(), r.get_mpz_t(), a.get_mpz_t(), b.get_mpz_t());
Parameters:
q(returns) the quotient of a/b.
r(returns) the remainder of a/b.
aany integer.
bany integer.

Definition at line 173 of file IntegerComputer.ih.

Referenced by DGtal::LightSternBrocot< TInteger, TQuotient, TMap >::Fraction::Fraction(), and DGtal::LighterSternBrocot< TInteger, TQuotient, TMap >::Fraction::Fraction().

{
q = a / b;
r = a % b;
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getFloorCeilDiv ( Integer fl,
Integer ce,
IntegerParamType  na,
IntegerParamType  nb 
) const
inline

Computes the floor and ceil value of na/nb.

Parameters:
fl(returns) the floor value of na/nb.
ce(returns) the ceil value of na/nb.
naany integer.
nbany integer.

Definition at line 242 of file IntegerComputer.ih.

{
_m_a = na;
_m_b = nb;
if ( isNegative( _m_b ) )
{
}
fl = ce = _m_a/_m_b;
if ( isNotZero( _m_a % _m_b ) )
{
if ( isNegativeOrZero( _m_a ) ) --fl;
if ( isPositiveOrZero( _m_a ) ) ++ce;
}
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getGcd ( Integer g,
IntegerParamType  a,
IntegerParamType  b 
) const
inline

Returns the greatest common divisor of a and b (a and b may be either positive or negative).

Parameters:
g(returns) the gcd of a and b.
aany integer.
bany integer.

Definition at line 303 of file IntegerComputer.ih.

References DGtal::abs().

{
// std::cerr << "gcd(" << a << ", " << b << ")=";
_m_a = abs( a );
_m_b = abs( b );
_m_a0 = max( _m_a, _m_b );
_m_a1 = min( _m_a, _m_b );
while ( isNotZero( _m_a1 ) )
{
}
g = _m_a0;
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getValidBezout ( Vector2I v,
const Point2I A,
const Vector2I u,
const Vector2I N,
IntegerParamType  c,
const Vector2I N2,
IntegerParamType  c2,
bool  compute_v = true 
) const
inline

Compute the valid bezout vector v of u such that A+v satifies the constraints C2 and such that A+v+u doesn't satify the constraint C2.

(A+v).N2 <= c2, (A+v+u).N2 > c2.

The constraint (N,c) is used like this: if the Bezout is such that (A+v).N > c, then v <- -v. Thus v is oriented toward the constraint.

If ! compute_v, v is used as is. The constraint N.(A+v) <= c should be satisfied.

Parameters:
v(modifies) a Bezout vector for u, with the constraints above.
Aany point in Z2.
uany vector in Z2.
Nany vector in Z2, defining the first constraint.
cthe integer for the first constraint.
N2any vector in Z2, defining the second constraint.
c2the integer for the second constraint.
compute_vtells if v should be recomputed (true) or is already given (false), default to true.

Definition at line 512 of file IntegerComputer.ih.

{
if ( compute_v )
{
v = extendedEuclid( -u[ 1 ], u[ 0 ], NumberTraits<Integer>::ONE );
_m_v0 = A + v;
if ( _m_c0 > c )
{
v[ 0 ] = -v[ 0 ];
v[ 1 ] = -v[ 1 ];
_m_v0 = A + v;
}
}
else _m_v0 = A + v;
_m_v0, u, N2, c2 );
_m_v1 = u * _m_a0; // floor value
v += _m_v1;
ASSERT( N2.dot( A + v ) <= c2 );
ASSERT( N2.dot( A + v + u ) > c2 );
}
template<typename TInteger >
bool DGtal::IntegerComputer< TInteger >::isNegative ( IntegerParamType  a)
inlinestatic
Parameters:
aany integer.
Returns:
true if a < 0.

Definition at line 98 of file IntegerComputer.ih.

{
return a < NumberTraits<Integer>::ZERO;
}
template<typename TInteger >
bool DGtal::IntegerComputer< TInteger >::isNegativeOrZero ( IntegerParamType  a)
inlinestatic
Parameters:
aany integer.
Returns:
true if a <= 0.

Definition at line 116 of file IntegerComputer.ih.

{
return a <= NumberTraits<Integer>::ZERO;
}
template<typename TInteger >
bool DGtal::IntegerComputer< TInteger >::isNotZero ( IntegerParamType  a)
inlinestatic
Parameters:
aany integer.
Returns:
true if a != 0.

Definition at line 80 of file IntegerComputer.ih.

{
}
template<typename TInteger >
bool DGtal::IntegerComputer< TInteger >::isPositive ( IntegerParamType  a)
inlinestatic
Parameters:
aany integer.
Returns:
true if a > 0.

Definition at line 89 of file IntegerComputer.ih.

template<typename TInteger >
bool DGtal::IntegerComputer< TInteger >::isPositiveOrZero ( IntegerParamType  a)
inlinestatic
Parameters:
aany integer.
Returns:
true if a >= 0.

Definition at line 107 of file IntegerComputer.ih.

{
}
template<typename TInteger >
bool DGtal::IntegerComputer< TInteger >::isValid ( ) const
inline

Checks the validity/consistency of the object.

Returns:
'true' if the object is valid, 'false' otherwise.

Definition at line 597 of file IntegerComputer.ih.

{
return true;
}
template<typename TInteger >
bool DGtal::IntegerComputer< TInteger >::isZero ( IntegerParamType  a)
inlinestatic
Parameters:
aany integer.
Returns:
true if a == 0.

Definition at line 71 of file IntegerComputer.ih.

{
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::max ( IntegerParamType  a,
IntegerParamType  b 
)
inlinestatic
Parameters:
aany integer.
bany integer.
Returns:
the maximum value of a and b.

Definition at line 137 of file IntegerComputer.ih.

{
return ( a >= b ) ? a : b;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::max ( IntegerParamType  a,
IntegerParamType  b,
IntegerParamType  c 
)
inlinestatic
Parameters:
aany integer.
bany integer.
cany integer.
Returns:
the maximum value of a, b and c.

Definition at line 146 of file IntegerComputer.ih.

{
return max( max(a,b), c );
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::min ( IntegerParamType  a,
IntegerParamType  b 
)
inlinestatic
Parameters:
aany integer.
bany integer.
Returns:
the minimum value of a and b.

Definition at line 155 of file IntegerComputer.ih.

{
return ( a <= b ) ? a : b;
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::min ( IntegerParamType  a,
IntegerParamType  b,
IntegerParamType  c 
)
inlinestatic
Parameters:
aany integer.
bany integer.
cany integer.
Returns:
the minimum value of a, b and c.

Definition at line 164 of file IntegerComputer.ih.

{
return min( min(a,b), c );
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger > & DGtal::IntegerComputer< TInteger >::operator= ( const Self other)
inline

Assignment.

Does nothing (member data are allocated, but their values are not used except during the execution of methods).

Parameters:
otherthe object to copy.
Returns:
a reference on 'this'.

Definition at line 60 of file IntegerComputer.ih.

{
return *this;
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::reduce ( Vector2I p) const
inline

Makes p irreducible.

Parameters:
pany vector in Z2.

Definition at line 394 of file IntegerComputer.ih.

References DGtal::gcd().

{
_m_a = gcd( p[ 0 ], p[ 1 ] );
p /= _m_a;
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::reduce ( Vector3I p) const
inline

Makes p irreducible.

Parameters:
pany vector in Z3.

Definition at line 543 of file IntegerComputer.ih.

References DGtal::gcd().

{
_m_a = gcd( p[ 0 ], p[ 1 ] );
_m_r = gcd( _m_a, p[ 2 ] );
p /= _m_r;
}
template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::selfDisplay ( std::ostream &  out) const
inline

Writes/Displays the object on an output stream.

Parameters:
outthe output stream where the object is written.

Definition at line 585 of file IntegerComputer.ih.

{
out << "[IntegerComputer]";
}
template<typename TInteger >
DGtal::IntegerComputer< TInteger >::Integer DGtal::IntegerComputer< TInteger >::staticGcd ( IntegerParamType  a,
IntegerParamType  b 
)
inlinestatic

Returns the greatest common divisor of a and b (a and b may be either positive or negative).

Parameters:
aany integer.
bany integer.
Returns:
the gcd of a and b.

Definition at line 264 of file IntegerComputer.ih.

References DGtal::abs().

{
Integer _m_a = abs( a );
Integer _m_b = abs( b );
Integer _m_a0 = max( _m_a, _m_b );
Integer _m_a1 = min( _m_a, _m_b );
while ( isNotZero( _m_a1 ) )
{
_m_r = _m_a0 % _m_a1;
_m_a0 = _m_a1;
_m_a1 = _m_r;
}
return _m_a0;
}

Field Documentation

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_a
mutableprivate

Used to store parameter a.

Definition at line 497 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_a0
mutableprivate

Used to for successive computation in gcd.

Definition at line 501 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_a1
mutableprivate

Used to for successive computation in gcd.

Definition at line 503 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_b
mutableprivate

Used to store parameter b.

Definition at line 499 of file IntegerComputer.h.

template<typename TInteger>
std::vector<Integer> DGtal::IntegerComputer< TInteger >::_m_bezout[4]
mutableprivate

Used to represent the variables during an extended Euclid computation, [0] are remainders, [1] are quotients, [2] are successive a, [3] are successive b.

Definition at line 511 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_c0
mutableprivate

Used to represent scalar products.

Definition at line 519 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_c1
mutableprivate

Used to represent scalar products.

Definition at line 521 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_c2
mutableprivate

Used to represent scalar products.

Definition at line 523 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_q
mutableprivate

Used to represent a quotient.

Definition at line 505 of file IntegerComputer.h.

template<typename TInteger>
Integer DGtal::IntegerComputer< TInteger >::_m_r
mutableprivate

Used to represent a remainder.

Definition at line 507 of file IntegerComputer.h.

template<typename TInteger>
Vector2I DGtal::IntegerComputer< TInteger >::_m_v
mutableprivate

Used to represent the Bezout vector.

Definition at line 513 of file IntegerComputer.h.

template<typename TInteger>
Vector2I DGtal::IntegerComputer< TInteger >::_m_v0
mutableprivate

Used to represent vectors.

Definition at line 515 of file IntegerComputer.h.

template<typename TInteger>
Vector2I DGtal::IntegerComputer< TInteger >::_m_v1
mutableprivate

Used to represent vectors.

Definition at line 517 of file IntegerComputer.h.


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