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

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

#include <DGtal/arithmetic/IntegerComputer.h>

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 ((concepts::CInteger< Integer >))
 
 BOOST_CONCEPT_ASSERT ((concepts::CUnsignedNumber< UnsignedInteger >))
 
 BOOST_CONCEPT_ASSERT ((concepts::CIntegralNumber< 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
 Used to store parameter a. More...
 
Integer _m_b
 Used to store parameter b. More...
 
Integer _m_a0
 Used to for successive computation in gcd. More...
 
Integer _m_a1
 Used to for successive computation in gcd. More...
 
Integer _m_q
 Used to represent a quotient. More...
 
Integer _m_r
 Used to represent a remainder. More...
 
std::vector< Integer_m_bezout [4]
 
Vector2I _m_v
 Used to represent the Bezout vector. More...
 
Vector2I _m_v0
 Used to represent vectors. More...
 
Vector2I _m_v1
 Used to represent vectors. More...
 
Integer _m_c0
 Used to represent scalar products. More...
 
Integer _m_c1
 Used to represent scalar products. More...
 
Integer _m_c2
 Used to represent scalar products. More...
 

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).
Examples
arithmetic/extended-euclid.cpp.

Definition at line 82 of file IntegerComputer.h.

Member Typedef Documentation

◆ Integer

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

Definition at line 87 of file IntegerComputer.h.

◆ IntegerParamType

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

Definition at line 88 of file IntegerComputer.h.

◆ Point2I

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

Definition at line 93 of file IntegerComputer.h.

◆ Point3I

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

Definition at line 95 of file IntegerComputer.h.

◆ Self

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

Definition at line 86 of file IntegerComputer.h.

◆ UnsignedInteger

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

Definition at line 90 of file IntegerComputer.h.

◆ UnsignedIntegerParamType

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

Definition at line 91 of file IntegerComputer.h.

◆ Vector2I

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

Definition at line 94 of file IntegerComputer.h.

◆ Vector3I

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

Definition at line 96 of file IntegerComputer.h.

Constructor & Destructor Documentation

◆ ~IntegerComputer()

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

Destructor.

◆ IntegerComputer() [1/2]

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

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).

◆ IntegerComputer() [2/2]

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

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.

Member Function Documentation

◆ abs()

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

Referenced by checkGenericPlane(), checkPlane(), checkPlaneGroupExtension(), checkWidth(), and testValidBezout().

◆ BOOST_CONCEPT_ASSERT() [1/3]

template<typename TInteger >
DGtal::IntegerComputer< TInteger >::BOOST_CONCEPT_ASSERT ( (concepts::CInteger< Integer >)  )

◆ BOOST_CONCEPT_ASSERT() [2/3]

template<typename TInteger >
DGtal::IntegerComputer< TInteger >::BOOST_CONCEPT_ASSERT ( (concepts::CIntegralNumber< UnsignedInteger >)  )

◆ BOOST_CONCEPT_ASSERT() [3/3]

template<typename TInteger >
DGtal::IntegerComputer< TInteger >::BOOST_CONCEPT_ASSERT ( (concepts::CUnsignedNumber< UnsignedInteger >)  )

◆ ceilDiv()

template<typename TInteger >
Integer DGtal::IntegerComputer< TInteger >::ceilDiv ( IntegerParamType  na,
IntegerParamType  nb 
) const

Computes the ceil value of na/nb.

Parameters
naany integer.
nbany integer.

Referenced by checkExtendWithManyPoints(), checkGenericPlane(), checkPlane(), checkPlaneGroupExtension(), checkWidth(), testCeilFloorDiv(), and unionComparisonTest().

◆ convergent()

template<typename TInteger >
Point2I DGtal::IntegerComputer< TInteger >::convergent ( const std::vector< Integer > &  quotients,
unsigned int  k 
) const

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).

Referenced by testCFrac().

◆ crossProduct()

template<typename TInteger >
Integer DGtal::IntegerComputer< TInteger >::crossProduct ( const Vector2I u,
const Vector2I v 
) const

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.

Referenced by testValidBezout().

◆ dotProduct() [1/2]

template<typename TInteger >
Integer DGtal::IntegerComputer< TInteger >::dotProduct ( const Vector2I u,
const Vector2I v 
) const

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.

Referenced by DGtal::FrechetShortcut< TIterator, TInteger >::Tools::angleVectVect(), testCoefficientIntersection(), testValidBezout(), and unionComparisonTest().

◆ dotProduct() [2/2]

template<typename TInteger >
Integer DGtal::IntegerComputer< TInteger >::dotProduct ( const Vector3I u,
const Vector3I v 
) const

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.

◆ extendedEuclid()

template<typename 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)).

The solution is chosen such that:

  • when c > 0,
    • a > 0 implies x >= 0, thus sgn(y)=-sgn(b)
    • a < 0 implies x <= 0, thus sgn(y)=-sgn(b)
  • when c < 0,
    • a > 0 implies x <= 0, thus sgn(y)=sgn(b)
    • a < 0 implies x >= 0, thus sgn(y)=sgn(b)
  • abs(x) <= abs(b * c )
  • abs(y) < abs(a * c )
Parameters
aany non-null integer.
bany non-null integer.
cany integer multiple of gcd(|a|,|b|).
Returns
a vector (x,y) solution to a x + b y = c.

Referenced by testExtendedEuclid().

◆ floorDiv()

template<typename TInteger >
Integer DGtal::IntegerComputer< TInteger >::floorDiv ( IntegerParamType  na,
IntegerParamType  nb 
) const

Computes the floor value of na/nb.

Parameters
naany integer.
nbany integer.

Referenced by testCeilFloorDiv(), testDSLSubsegment(), and unionComparisonTest().

◆ gcd()

template<typename TInteger >
Integer DGtal::IntegerComputer< TInteger >::gcd ( IntegerParamType  a,
IntegerParamType  b 
) const

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.

Referenced by testCFrac(), testDSLSubsegment(), testExtendedEuclid(), testGCD(), testInitFraction(), testPatterns(), testReducedFraction(), testSubStandardDSLQ0(), and unionComparisonTest().

◆ getCFrac() [1/2]

template<typename TInteger >
template<typename OutputIterator >
Integer DGtal::IntegerComputer< TInteger >::getCFrac ( OutputIterator  outIt,
IntegerParamType  a,
IntegerParamType  b 
) const

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.

◆ getCFrac() [2/2]

template<typename TInteger >
Integer DGtal::IntegerComputer< TInteger >::getCFrac ( std::vector< Integer > &  quotients,
IntegerParamType  a,
IntegerParamType  b 
) const

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.

Referenced by testCFrac(), and testReducedFraction().

◆ getCoefficientIntersection()

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getCoefficientIntersection ( Integer fl,
Integer ce,
const Vector2I p,
const Vector2I u,
const Vector2I N,
IntegerParamType  c 
) const

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.

Referenced by testCoefficientIntersection().

◆ getCrossProduct()

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getCrossProduct ( Integer cp,
const Vector2I u,
const Vector2I v 
) const

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.

◆ getDotProduct() [1/2]

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getDotProduct ( Integer dp,
const Vector2I u,
const Vector2I v 
) const

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.

◆ getDotProduct() [2/2]

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getDotProduct ( Integer dp,
const Vector3I u,
const Vector3I v 
) const

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.

◆ getEuclideanDiv()

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getEuclideanDiv ( Integer q,
Integer r,
IntegerParamType  a,
IntegerParamType  b 
) const

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.

◆ getFloorCeilDiv()

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getFloorCeilDiv ( Integer fl,
Integer ce,
IntegerParamType  na,
IntegerParamType  nb 
) const

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.

Referenced by testCeilFloorDiv().

◆ getGcd()

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::getGcd ( Integer g,
IntegerParamType  a,
IntegerParamType  b 
) const

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.

Referenced by testGCD().

◆ getValidBezout()

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

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.

Referenced by testValidBezout().

◆ isNegative()

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

◆ isNegativeOrZero()

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

◆ isNotZero()

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

◆ isPositive()

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

◆ isPositiveOrZero()

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

◆ isValid()

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

Checks the validity/consistency of the object.

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

◆ isZero()

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

Referenced by testCeilFloorDiv(), and testGCD().

◆ max() [1/2]

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

◆ max() [2/2]

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

◆ min() [1/2]

template<typename TInteger >
static Integer DGtal::IntegerComputer< TInteger >::min ( IntegerParamType  a,
IntegerParamType  b 
)
static
Parameters
aany integer.
bany integer.
Returns
the minimum value of a and b.

◆ min() [2/2]

template<typename TInteger >
static Integer DGtal::IntegerComputer< TInteger >::min ( IntegerParamType  a,
IntegerParamType  b,
IntegerParamType  c 
)
static
Parameters
aany integer.
bany integer.
cany integer.
Returns
the minimum value of a, b and c.

◆ operator=()

template<typename TInteger >
Self & DGtal::IntegerComputer< TInteger >::operator= ( const Self other)

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'.

◆ reduce() [1/2]

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::reduce ( Vector2I p) const

Makes p irreducible.

Parameters
pany vector in Z2.

Referenced by testValidBezout().

◆ reduce() [2/2]

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::reduce ( Vector3I p) const

Makes p irreducible.

Parameters
pany vector in Z3.

◆ selfDisplay()

template<typename TInteger >
void DGtal::IntegerComputer< TInteger >::selfDisplay ( std::ostream &  out) const

Writes/Displays the object on an output stream.

Parameters
outthe output stream where the object is written.

◆ staticGcd()

template<typename TInteger >
static Integer DGtal::IntegerComputer< TInteger >::staticGcd ( IntegerParamType  a,
IntegerParamType  b 
)
static

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.

Referenced by exhaustiveTestLatticePolytope2D().

Field Documentation

◆ _m_a

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

Used to store parameter a.

Definition at line 508 of file IntegerComputer.h.

◆ _m_a0

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

Used to for successive computation in gcd.

Definition at line 512 of file IntegerComputer.h.

◆ _m_a1

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

Used to for successive computation in gcd.

Definition at line 514 of file IntegerComputer.h.

◆ _m_b

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

Used to store parameter b.

Definition at line 510 of file IntegerComputer.h.

◆ _m_bezout

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 522 of file IntegerComputer.h.

◆ _m_c0

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

Used to represent scalar products.

Definition at line 530 of file IntegerComputer.h.

◆ _m_c1

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

Used to represent scalar products.

Definition at line 532 of file IntegerComputer.h.

◆ _m_c2

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

Used to represent scalar products.

Definition at line 534 of file IntegerComputer.h.

◆ _m_q

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

Used to represent a quotient.

Definition at line 516 of file IntegerComputer.h.

◆ _m_r

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

Used to represent a remainder.

Definition at line 518 of file IntegerComputer.h.

◆ _m_v

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

Used to represent the Bezout vector.

Definition at line 524 of file IntegerComputer.h.

◆ _m_v0

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

Used to represent vectors.

Definition at line 526 of file IntegerComputer.h.

◆ _m_v1

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

Used to represent vectors.

Definition at line 528 of file IntegerComputer.h.


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