Lattice Tester Guide  1.0-9
Software Package For Testing The Uniformity Of Integral Lattices In The Real Space
LatticeTester Namespace Reference

Lattice namespace. More...

Namespaces

 CoordinateSets
 A namespace containing different implementation of sets of coordinates.
 
 Random
 This class generates random numbers (in fact pseudo-random numbers).
 

Classes

class  BasisConstruction
 This class implements general methods to perform a lattice basis construction from a set of vectors, as well as general methods to obtain the dual lattice basis depending on the current lattice basis. More...
 
class  Config
 This class contains a configuration that can be passed to a LatticeAnalysis object to perform a computation. More...
 
class  Coordinates
 This is basically a std::set<std::size_t>. More...
 
class  IntFactor
 The objects of this class are the "prime" factors in the decomposition of a positive integer. More...
 
class  IntLattice
 This class is a skeleton for the implementation of different types of lattices. More...
 
class  IntLatticeBasis
 This class represents a lattice and its basis and offers tools to do basic manipulations on lattice bases. More...
 
class  Lacunary
 This class represents a set of indices with lacunary values. More...
 
class  LatticeAnalysis
 Objects of this class can perform various tests on lattices. More...
 
class  NormaBestBound
 This class implements upper bounds on the lenght of the shortest nonzero vector in a lattice. More...
 
class  NormaBestLat
 This class implements upper bounds on the lenght of the shortest nonzero vector in a lattice. More...
 
class  NormaLaminated
 This class implements upper bounds on the lenght of the shortest nonzero vector in a lattice. More...
 
class  Normalizer
 Classes which inherit from this base class are used in implementing bounds on the length of the shortest nonzero vector in a lattice. More...
 
class  NormaMinkL1
 This class implements theoretical bounds on the length of the shortest nonzero vector in a lattice, based on the densest sphere packing in space. More...
 
class  NormaMinkowski
 This class implements Minkowski’s theoretical LOWER bound on the length of the shortest non-zero vector in a lattice. More...
 
class  NormaPalpha
 This class implements theoretical bounds on the values of \(P_{\alpha}\) for a lattice (see class Palpha). More...
 
class  NormaRogers
 This class implements upper bounds on the lenght of the shortest nonzero vector in a lattice. More...
 
class  OrderDependentWeights
 Order-dependent weights. More...
 
class  ParamReader
 Utility class that can be used to read different kind of data from a file. More...
 
class  PODWeights
 Product and order-dependent (POD) weights. More...
 
class  ProductWeights
 Product weights. More...
 
class  ProjectionDependentWeights
 Projection-dependent weights. More...
 
class  Rank1Lattice
 This class implements a general rank 1 lattice. More...
 
class  Reducer
 This class implements (or wraps from NTL) all the functions that are needed to reduce a basis. More...
 
class  Types
 Sets standard typedef’s for the types that can be used in LatticeTester. More...
 
class  UniformWeights
 This class is used to implement the same weight for all projections. More...
 
class  Weights
 Abstract class representing Weights for figures of merit. More...
 
class  Writer
 This is an abstract class that represents an interface to Writer classes. More...
 
class  WriterRes
 This class is a simple implementation of the Writer abstract class to write in plain text format on the stream. More...
 

Enumerations

enum  NormType
 Indicates which norm is used to measure the length of vectors. More...
 
enum  OutputType
 Indicates in which form and where the results will be sent. More...
 
enum  ProblemType
 An enum listing the problems that LatticeTester can solve. More...
 
enum  PrecisionType
 Indicates in which precision the NTL algorithms will be perfoms : DOUBLE – double QUADRUPLE – quad_float (quasi quadruple precision) this is useful when roundoff errors can cause problems EXPONENT – xdouble (extended exponent doubles) this is useful when numbers get too big ARBITRARY – RR (arbitrary precision floating point) this is useful for large precision and magnitudes Generally speaking, the choice DOUBLE will be the fastest, but may be prone to roundoff errors and/or overflow. More...
 
enum  PrimeType
 Indicates whether an integer is prime, probably prime, composite or its status is unknown (or don’t care). More...
 
enum  CriterionType
 Gives the merit criterion for ranking generators or lattices. More...
 
enum  NormaType
 Indicates which normalization is used to compute \(S_t\) in the spectral test, for each dimension \(t\). More...
 
enum  CalcType
 Indicates which type of calculation is considered for the \(P_{\alpha}\) test. More...
 
enum  PreReductionType
 A list of all the possible lattice reduction implemented in LatticeTester. More...
 

Functions

std::int64_t lFactorial (int t)
 Calculates \(t!\), the factorial of \(t\). More...
 
double Digamma (double x)
 Returns the value of the logarithmic derivative of the Gamma function \(\psi(x) = \Gamma'(x) / \Gamma(x)\). More...
 
double BernoulliPoly (int n, double x)
 Evaluates the Bernoulli polynomial \(B_n(x)\) of degree \(n\) at \(x\). More...
 
double Harmonic (std::int64_t n)
 Computes the \(n\)-th harmonic number \(H_n = \sum_{j=1}^n 1/j\). More...
 
double Harmonic2 (std::int64_t n)
 Computes the sum

\[ \sideset{}{'}\sum_{-n/2<j\le n/2}\; \frac 1{|j|}, \]

where the symbol \(\sum^\prime\) means that the term with \(j=0\) is excluded from the sum. More...

 
double FourierC1 (double x, std::int64_t n)
 Computes and returns the value of the series (see [10])

\[ S(x, n) = \sum_{j=1}^{n} \frac{\cos(2\pi j x)}{j}. \]

Restrictions: \(n>0\) and \(0 \le x \le 1\). More...

 
double FourierE1 (double x, std::int64_t n)
 Computes and returns the value of the series

\[ G(x, n) = \sideset{}{'}\sum_{-n/2<h\le n/2}\; \frac{e^{2\pi i h x}}{|h|}, \]

where the symbol \(\sum^\prime\) means that the term with \(h=0\) is excluded from the sum, and assuming that the imaginary part of \(G(x, n)\) vanishes. More...

 
template<typename T >
void swap9 (T &x, T &y)
 Takes references to two variables of a generic type and swaps their content. More...
 
toString functions

Useful functions for printing the enum constants in this module.

Returns the enum constants in this module as strings.

std::string toStringNorm (NormType)
 
std::string toStringPrime (PrimeType)
 
std::string toStringCriterion (CriterionType)
 
std::string toStringProblem (ProblemType)
 
std::string toStringNorma (NormaType)
 
std::string toStringCalc (CalcType)
 
std::string toStringPreRed (PreReductionType)
 
std::string toStringOutput (OutputType)
 
std::string toStringPrecision (PrecisionType)
 
Random numbers

All the functions of this module use LFSR258 as an underlying source for pseudo-random numbers.

A free (as in freedom) implementation of this generator can be found at http://simul.iro.umontreal.ca/ as well as the article presenting it. All the functions generating some sort of random number will advance an integer version of LFSR258 by one state and output a transformation of the state to give a double, an int or bits.

double RandU01 ()
 Returns a random number in \([0, 1)\). More...
 
int RandInt (int i, int j)
 Return a uniform pseudo-random integer in \([i, j]\). More...
 
std::uint64_t RandBits (int s)
 Returns the first s pseudo-random bits of the underlying RNG in the form of a s-bit integer. More...
 
void SetSeed (std::uint64_t seed)
 Sets the seed of the generator. More...
 
Mathematical functions

These are complete reimplementation of certain mathematical functions, or wrappers for standard C/C++ functions.

double mysqrt (double x)
 Returns \(\sqrt{x}\) for \(x\ge0\), and \(-1\) for \(x < 0\). More...
 
template<typename T >
double Lg (const T &x)
 Returns the logarithm of \(x\) in base 2. More...
 
double Lg (std::int64_t x)
 Returns the logarithm of \(x\) in base 2. More...
 
template<typename Scal >
Scal abs (Scal x)
 Returns the absolute value of \(x\). More...
 
template<typename T >
std::int64_t sign (const T &x)
 Returns the sign of x. More...
 
template<typename Real >
Real Round (Real x)
 Returns the value of x rounded to the NEAREST integer value. More...
 
std::int64_t Factorial (int t)
 Calculates \(t!\), the factorial of \(t\) and returns it as an std::int64_t. More...
 
Division and modular arithmetic
Remarks
Richard: Pour certaines fonctions, les résultats sont mis dans les premiers arguments de la fonction pour être compatible avec NTL; pour d’autres, ils sont mis dans les derniers arguments pour être compatible avec notre ancienne version de LatMRG en Modula-2. Plutôt détestable. Je crois qu’il faudra un jour réarranger les arguments des fonctions pour qu’elles suivent toutes la même convention que NTL.

This module offers function to perform division and find remainders in a standard way. These functions are usefull in the case where one wants to do divisions or find remainders of operations with negative operands. The reason is that NTL and primitive types do not use the same logic when doing calculations on negative numbers.

Basically, NTL will always floor a division and C++ will always truncate a division (which effectively means the floor function is replaced by a roof function if the answer is a negative number). When calculating the remainder of x/y, both apply the same logic but get a different result because they do not do the same division. In both representations, we have that

\[ y\cdot(x/y) + x%y = x. \]

It turns out that, with negative values, NTL will return an integer with the same sign as y where C++ will return an integer of opposite sign (but both will return the same number modulo y).

template<typename Int >
void Quotient (const Int &a, const Int &b, Int &q)
 Computes a/b, truncates the fractionnal part and puts the result in q. More...
 
template<typename Real >
void Modulo (const Real &a, const Real &b, Real &r)
 Computes the remainder of a/b and stores its positive equivalent mod b in r. More...
 
template<typename Real >
void Divide (Real &q, Real &r, const Real &a, const Real &b)
 Computes the quotient \(q = a/b\) and remainder \(r = a \bmod b\). More...
 
template<typename Real >
void DivideRound (const Real &a, const Real &b, Real &q)
 Computes \(a/b\), rounds the result to the nearest integer and returns the result in \(q\). More...
 
std::int64_t gcd (std::int64_t a, std::int64_t b)
 Returns the value of the greatest common divisor of \(a\) and \(b\) by using Stein's binary GCD algorithm. More...
 
template<typename Int >
void Euclide (const Int &A, const Int &B, Int &C, Int &D, Int &E, Int &F, Int &G)
 This method computes the greater common divisor of A and B with Euclid's algorithm. More...
 
Vectors

These are utilities to manipulate vectors ranging from instanciation to scalar product.

template<typename Real >
void CreateVect (Real *&A, int d)
 Allocates memory to A as an array of Real of dimension d and initializes its elements to 0. More...
 
template<typename Vect >
void CreateVect (Vect &A, int d)
 Creates the vector A of dimensions d+1 and initializes its elements to 0. More...
 
template<typename Real >
void DeleteVect (Real *&A)
 Frees the memory used by the vector A. More...
 
template<typename Vect >
void DeleteVect (Vect &A)
 Frees the memory used by the vector A, destroying all the elements it contains. More...
 
template<typename Real >
void SetZero (Real *A, int d)
 Sets the first d of A to 0. More...
 
template<typename Vect >
void SetZero (Vect &A, int d)
 Sets the first d components of A to 0. More...
 
template<typename Real >
void SetValue (Real *A, int d, const Real &x)
 Sets the first d components of A to the value x. More...
 
template<typename Vect >
std::string toString (const Vect &A, int c, int d, const char *sep=" ")
 Returns a string containing A[c] to A[d-1] formated as [A[c]sep...sepA[d-1]]. More...
 
template<typename Vect >
std::string toString (const Vect &A, int d)
 Returns a string containing the first d components of the vector A as a string. More...
 
template<typename Int , typename Vect1 , typename Vect2 , typename Scal >
void ProdScal (const Vect1 &A, const Vect2 &B, int n, Scal &D)
 Computes the scalar product of vectors A and B truncated to their n first components, then puts the result in D. More...
 
template<typename IntVec >
void Invert (const IntVec &A, IntVec &B, int n)
 Takes an input vector A of dimension n+1 and fill the vector B with the values [-A[n] -A[n-1] ... More...
 
template<typename Vect , typename Scal >
void CalcNorm (const Vect &V, int n, Scal &S, NormType norm)
 Computes the norm norm of vector V trunctated to its n first components, and puts the result in S. More...
 
template<typename Vect >
void CopyVect (Vect &A, const Vect &B, int n)
 Copies the first n components of vector B into vector A. More...
 
template<typename Vect1 , typename Vect2 , typename Scal >
void ModifVect (Vect1 &A, const Vect2 &B, Scal x, int n)
 Adds the first n components of vector B multiplied by x to first n components of vector A. More...
 
template<typename Vect >
void ChangeSign (Vect &A, int n)
 Changes the sign (multiplies by -1) the first n components of vector A. More...
 
std::int64_t GCD2vect (std::vector< std::int64_t > V, int k, int n)
 Computes the greatest common divisor of V[k],...,V[n-1]. More...
 
Matrices
template<typename Real >
void CreateMatr (Real **&A, int d)
 Allocates memory to a square matrix A of dimensions \(d \times d\) and initializes its elements to 0. More...
 
template<typename Real >
void CreateMatr (Real **&A, int line, int col)
 Allocates memory for the matrix A of dimensions \(\text{line} \times \text{col}\) and initializes its elements to 0. More...
 
template<typename IntMat >
void CreateMatr (IntMat &A, int d)
 Resizes the matrix A to a square matrix of dimensions d*d and re-initializes its elements to 0. More...
 
template<typename IntMat >
void CreateMatr (IntMat &A, int line, int col)
 Resizes the matrix A to a matrix of dimensions \(line \times col\) and re-initializes its elements to 0. More...
 
template<typename Real >
void DeleteMatr (Real **&A, int d)
 Frees the memory used by the \(d \times d\) matrix A. More...
 
template<typename Real >
void DeleteMatr (Real **&A, int line, int col)
 Frees the memory used by the matrix A of dimension \(\text{line} \times \text{col}\). More...
 
template<typename IntMat >
void DeleteMatr (IntMat &A)
 Calls the clear() method on A. More...
 
template<typename Matr >
void CopyMatr (Matr &A, const Matr &B, int n)
 Copies the \(n \times n\) submatrix of the first lines and columns of B into matrix A. More...
 
template<typename Matr >
void CopyMatr (Matr &A, const Matr &B, int line, int col)
 Copies the \(\text{line} \times col\) submatrix of the first lines and columns of B into matrix A. More...
 
template<typename MatT >
std::string toStr (const MatT &mat, int d1, int d2)
 Returns a string that is a representation of mat. More...
 
template<typename Int >
bool CheckTriangular (const NTL::matrix< Int > &A, long dim, const Int m)
 Checks that the upper \(\text{dim} \times \text{dim}\) submatrix of A is triangular modulo m. More...
 
template<typename Matr , typename Int >
void Triangularization (Matr &W, Matr &V, int lin, int col, const Int &m)
 Takes a set of generating vectors in the matrix W and iteratively transforms it into an upper triangular lattice basis into the matrix V. More...
 
template<typename Matr , typename Int >
void CalcDual (const Matr &A, Matr &B, int d, const Int &m)
 Takes an upper triangular basis A and computes a dual lattice basis modulo m to this matrix. More...
 
Debugging functions
void MyExit (int status, std::string msg)
 Special exit function. More...
 
Printing functions and operators
template<class K , class T , class C , class A >
std::ostream & operator<< (std::ostream &out, const std::map< K, T, C, A > &x)
 Streaming operator for maps. More...
 
template<class T1 , class T2 >
std::ostream & operator<< (std::ostream &out, const std::pair< T1, T2 > &x)
 Streaming operator for vectors. More...
 
template<class T , class A >
std::ostream & operator<< (std::ostream &out, const std::vector< T, A > &x)
 Streaming operator for vectors. More...
 
template<class K , class C , class A >
std::ostream & operator<< (std::ostream &out, const std::set< K, C, A > &x)
 Streaming operator for sets. More...
 

Variables

const double MAX_LONG_DOUBLE = 9007199254740992.0
 Maximum integer that can be represented exactly as a double: \(2^{53}\). More...
 
const std::int64_t TWO_EXP []
 Table of powers of 2: TWO_EXP[ \(i\)] \(= 2^i\), \(i=0, 1, …, 63\). More...
 

Detailed Description

Lattice namespace.

Enumeration Type Documentation

◆ CalcType

Indicates which type of calculation is considered for the \(P_{\alpha}\) test.

PAL is for the \(P_{\alpha}\) test.
BAL is for the bound on the \(P_{\alpha}\) test.
NORMPAL is for the \(P_{\alpha}\) test PAL, with the result normalized over the BAL bound.
SEEKPAL is for the \(P_{\alpha}\) seek, which searches for good values of the multiplier.

◆ CriterionType

Gives the merit criterion for ranking generators or lattices.

LENGTH: Only using the length of the shortest vector as a criterion. BEYER: the figure of merit is the Beyer quotient \(Q_T\).
SPECTRAL: the figure of merit \(S_T\) is based on the spectral test.
PALPHA: the figure of merit is based on \(P_{\alpha}\).
BOUND_JS: the figure of merit is based on the Joe-Sinescu bound [24].

◆ NormaType

Indicates which normalization is used to compute \(S_t\) in the spectral test, for each dimension \(t\).

BESTLAT: the value used for \(d_t^*\) corresponds to the best lattice.
BESTLAT: the value used for \(d_t^*\) corresponds to the best bound known to us.
LAMINATED: the value used for \(d_t^*\) corresponds to the best laminated lattice.
ROGERS: the value for \(d_t^*\) is obtained from Rogers’ bound on the density of sphere packing.
MINKOWSKI: the value for \(d_t^*\) is obtained from Minkowski’ theoretical bounds on the length of the shortest nonzero vector in the lattice using the \({\mathcal{L}}_2\) norm.
MINKL1: the value for \(d_t^*\) is obtained from the theoretical bounds on the length of the shortest nonzero vector in the lattice using the \({\mathcal{L}}_1\) norm.
NONE: no normalization will be used.

◆ NormType

Indicates which norm is used to measure the length of vectors.

For \(X = (x_1,…,x_t)\),

SUPNORM corresponds to \(\Vert X\Vert= \max(|x_1|,…,|x_t|)\).
L1NORM corresponds to \(\Vert X\Vert= |x_1|+\cdots+|x_t|\).
L2NORM corresponds to \(\Vert X\Vert= (x_1^2+\cdots+x_t^2)^{1/2}\).
ZAREMBANORM corresponds to \(\Vert X\Vert= \max(1, |x_1|)\cdots\max(1, |x_t|)\).

◆ OutputType

Indicates in which form and where the results will be sent.

TERM: the results will appear only on the terminal screen.
RES: the results will be in plain text format and sent to a file with extension .res.
TEX: the results will be in LaTeX format and sent to a file with extension .tex.
GEN: the results will be sent to a file with extension .gen.

◆ PrecisionType

Indicates in which precision the NTL algorithms will be perfoms : DOUBLE – double QUADRUPLE – quad_float (quasi quadruple precision) this is useful when roundoff errors can cause problems EXPONENT – xdouble (extended exponent doubles) this is useful when numbers get too big ARBITRARY – RR (arbitrary precision floating point) this is useful for large precision and magnitudes Generally speaking, the choice DOUBLE will be the fastest, but may be prone to roundoff errors and/or overflow.

◆ PreReductionType

A list of all the possible lattice reduction implemented in LatticeTester.

NORPRERED: no reduction DIETER: Pairwise reduction LLL: LLL reduction BKZ: block Korkine-Zolotarev reduction FULL: shortest vector search.

◆ PrimeType

Indicates whether an integer is prime, probably prime, composite or its status is unknown (or don’t care).

◆ ProblemType

An enum listing the problems that LatticeTester can solve.

Function Documentation

◆ abs()

template<typename Scal >
Scal LatticeTester::abs ( Scal  x)
inline

Returns the absolute value of \(x\).

◆ BernoulliPoly()

double LatticeTester::BernoulliPoly ( int  n,
double  x 
)

Evaluates the Bernoulli polynomial \(B_n(x)\) of degree \(n\) at \(x\).

The first Bernoulli polynomials are:

\begin{align*} B_0(x) &= 1 \\ B_1(x) &= x - 1/2 \\ B_2(x) &= x^2-x+1/6 \\ B_3(x) &= x^3 - 3x^2/2 + x/2 \\ B_4(x) &= x^4-2x^3+x^2-1/30 \\ B_5(x) &= x^5 - 5x^4/2 + 5x^3/3 - x/6 \\ B_6(x) &= x^6-3x^5+5x^4/2-x^2/2+1/42 \\ B_7(x) &= x^7 - 7x^6/2 + 7x^5/2 - 7x^3/6 + x/6 \\ B_8(x) &= x^8-4x^7+14x^6/3 - 7x^4/3 +2x^2/3-1/30. \end{align*}

Only degrees \(n \le 8\) are programmed for now.

◆ CalcDual()

template<typename Matr , typename Int >
void LatticeTester::CalcDual ( const Matr &  A,
Matr &  B,
int  d,
const Int &  m 
)

Takes an upper triangular basis A and computes a dual lattice basis modulo m to this matrix.

For this algorithm to work, A as to be upper triangular and all the coefficients on the diagonal have to divide m.

For B to be a m-dual to A, we have to have that \(AB^t = mI\). It is quite easy to show that, knowing A is upper triangular, B will be a lower triangular matrix with A(i,i)*B(i,i) = m for all i and \( A_i \cdot B_j = 0\) for \(i\neq j\). To get the second condition, we simply have to recursively take for each line

\[B_{i,j} = -\frac{1}{A_{j,j}}\sum_{k=j+1}^i A_{j,k} B_{i,k}.\]

◆ CalcNorm()

template<typename Vect , typename Scal >
void LatticeTester::CalcNorm ( const Vect &  V,
int  n,
Scal &  S,
NormType  norm 
)
inline

Computes the norm norm of vector V trunctated to its n first components, and puts the result in S.

Scal has to be a floating point type.

◆ ChangeSign()

template<typename Vect >
void LatticeTester::ChangeSign ( Vect &  A,
int  n 
)
inline

Changes the sign (multiplies by -1) the first n components of vector A.

◆ CheckTriangular()

template<typename Int >
bool LatticeTester::CheckTriangular ( const NTL::matrix< Int > &  A,
long  dim,
const Int  m 
)

Checks that the upper \(\text{dim} \times \text{dim}\) submatrix of A is triangular modulo m.

This will return true if all the elements under the diagonal are equal to zero modulo m and false otherwise. If m is 0, this function simply verifies that the matrix is triangular.

◆ CopyMatr() [1/2]

template<typename Matr >
void LatticeTester::CopyMatr ( Matr &  A,
const Matr &  B,
int  n 
)
inline

Copies the \(n \times n\) submatrix of the first lines and columns of B into matrix A.

This function does not check for sizes, so A and B both have to be at leat \(n \times n\).

◆ CopyMatr() [2/2]

template<typename Matr >
void LatticeTester::CopyMatr ( Matr &  A,
const Matr &  B,
int  line,
int  col 
)
inline

Copies the \(\text{line} \times col\) submatrix of the first lines and columns of B into matrix A.

This function does not check for sizes, so A and B both have to be at leat \(line \times col\).

◆ CopyVect()

template<typename Vect >
void LatticeTester::CopyVect ( Vect &  A,
const Vect &  B,
int  n 
)
inline

Copies the first n components of vector B into vector A.

◆ CreateMatr() [1/4]

template<typename Real >
void LatticeTester::CreateMatr ( Real **&  A,
int  d 
)
inline

Allocates memory to a square matrix A of dimensions \(d \times d\) and initializes its elements to 0.

◆ CreateMatr() [2/4]

template<typename Real >
void LatticeTester::CreateMatr ( Real **&  A,
int  line,
int  col 
)
inline

Allocates memory for the matrix A of dimensions \(\text{line} \times \text{col}\) and initializes its elements to 0.

◆ CreateMatr() [3/4]

template<typename IntMat >
void LatticeTester::CreateMatr ( IntMat &  A,
int  d 
)
inline

Resizes the matrix A to a square matrix of dimensions d*d and re-initializes its elements to 0.

◆ CreateMatr() [4/4]

template<typename IntMat >
void LatticeTester::CreateMatr ( IntMat &  A,
int  line,
int  col 
)
inline

Resizes the matrix A to a matrix of dimensions \(line \times col\) and re-initializes its elements to 0.

◆ CreateVect() [1/2]

template<typename Real >
void LatticeTester::CreateVect ( Real *&  A,
int  d 
)
inline

Allocates memory to A as an array of Real of dimension d and initializes its elements to 0.

Real has to be a numeric type.

◆ CreateVect() [2/2]

template<typename Vect >
void LatticeTester::CreateVect ( Vect &  A,
int  d 
)
inline

Creates the vector A of dimensions d+1 and initializes its elements to 0.

The type Vect has to have a resize(integer_type) method that sets the size of the instance to the value of the argument.

◆ DeleteMatr() [1/3]

template<typename Real >
void LatticeTester::DeleteMatr ( Real **&  A,
int  d 
)
inline

Frees the memory used by the \(d \times d\) matrix A.

This will not free all the memory allocated to A if A is of greater dimension and it can cause a memory leak.

◆ DeleteMatr() [2/3]

template<typename Real >
void LatticeTester::DeleteMatr ( Real **&  A,
int  line,
int  col 
)
inline

Frees the memory used by the matrix A of dimension \(\text{line} \times \text{col}\).

This will not free all the memory allocated to A if A is of greater dimension and it can cause a memory leak.

◆ DeleteMatr() [3/3]

template<typename IntMat >
void LatticeTester::DeleteMatr ( IntMat &  A)
inline

Calls the clear() method on A.

A has to have a clear() method that frees the memory allocated to it.

◆ DeleteVect() [1/2]

template<typename Real >
void LatticeTester::DeleteVect ( Real *&  A)
inline

Frees the memory used by the vector A.

This calls delete[] on A so trying to access A after using this is unsafe.

◆ DeleteVect() [2/2]

template<typename Vect >
void LatticeTester::DeleteVect ( Vect &  A)
inline

Frees the memory used by the vector A, destroying all the elements it contains.

Vect type has to have a clear() method that deallocates all the elements in the vector.

◆ Digamma()

double LatticeTester::Digamma ( double  x)

Returns the value of the logarithmic derivative of the Gamma function \(\psi(x) = \Gamma'(x) / \Gamma(x)\).

◆ Divide()

template<typename Real >
void LatticeTester::Divide ( Real &  q,
Real &  r,
const Real &  a,
const Real &  b 
)
inline

Computes the quotient \(q = a/b\) and remainder \(r = a \bmod b\).

Truncates \(q\) to the nearest integer towards 0. One always has \(a = qb + r\) and \(|r| < |b|\). This works with std::int64_t, NTL::ZZ and real numbers.

\(a\) \(b\) \(q\) \(r\)
\(5\) 3 1 \(2\)
\(-5\) \(3\) \(-1\) \(-2\)
\(5\) \(-3\) \(-1\) \(2\)
\(-5\) \(-3\) \(1\) \(-2\)

◆ DivideRound()

template<typename Real >
void LatticeTester::DivideRound ( const Real &  a,
const Real &  b,
Real &  q 
)
inline

Computes \(a/b\), rounds the result to the nearest integer and returns the result in \(q\).

This works with std::int64_t, NTL::ZZ and real numbers.

◆ Euclide()

template<typename Int >
void LatticeTester::Euclide ( const Int &  A,
const Int &  B,
Int &  C,
Int &  D,
Int &  E,
Int &  F,
Int &  G 
)

This method computes the greater common divisor of A and B with Euclid's algorithm.

This will store this gcd in G and also the linear combination that permits to get G from A and B. This function should work with std::int64_t and NTL::ZZ.

For \(A\) and \(B\) this will assign to \(C\), \(D\), \(E\), \(F\) and \(G\) values such that:

\begin{align*} C a + D b & = G = \mbox{GCD } (a,b)\\ E a + F b & = 0. \end{align*}

◆ Factorial()

std::int64_t LatticeTester::Factorial ( int  t)

Calculates \(t!\), the factorial of \(t\) and returns it as an std::int64_t.

◆ FourierC1()

double LatticeTester::FourierC1 ( double  x,
std::int64_t  n 
)

Computes and returns the value of the series (see [10])

\[ S(x, n) = \sum_{j=1}^{n} \frac{\cos(2\pi j x)}{j}. \]

Restrictions: \(n>0\) and \(0 \le x \le 1\).

◆ FourierE1()

double LatticeTester::FourierE1 ( double  x,
std::int64_t  n 
)

Computes and returns the value of the series

\[ G(x, n) = \sideset{}{'}\sum_{-n/2<h\le n/2}\; \frac{e^{2\pi i h x}}{|h|}, \]

where the symbol \(\sum^\prime\) means that the term with \(h=0\) is excluded from the sum, and assuming that the imaginary part of \(G(x, n)\) vanishes.

Restrictions: \(n>0\) and \(0 \le x \le 1\).

◆ gcd()

std::int64_t LatticeTester::gcd ( std::int64_t  a,
std::int64_t  b 
)

Returns the value of the greatest common divisor of \(a\) and \(b\) by using Stein's binary GCD algorithm.

◆ GCD2vect()

std::int64_t LatticeTester::GCD2vect ( std::vector< std::int64_t >  V,
int  k,
int  n 
)
inline

Computes the greatest common divisor of V[k],...,V[n-1].

◆ Harmonic()

double LatticeTester::Harmonic ( std::int64_t  n)

Computes the \(n\)-th harmonic number \(H_n = \sum_{j=1}^n 1/j\).

◆ Harmonic2()

double LatticeTester::Harmonic2 ( std::int64_t  n)

Computes the sum

\[ \sideset{}{'}\sum_{-n/2<j\le n/2}\; \frac 1{|j|}, \]

where the symbol \(\sum^\prime\) means that the term with \(j=0\) is excluded from the sum.

◆ Invert()

template<typename IntVec >
void LatticeTester::Invert ( const IntVec &  A,
IntVec &  B,
int  n 
)
inline

Takes an input vector A of dimension n+1 and fill the vector B with the values [-A[n] -A[n-1] ...

-A[1][1]. B is assumed to be of dimension at least n+1.

◆ lFactorial()

std::int64_t LatticeTester::lFactorial ( int  t)

Calculates \(t!\), the factorial of \(t\).

Might throw if t is too large or if std::int64_t can't contain the factorial asked for.

◆ Lg() [1/2]

template<typename T >
double LatticeTester::Lg ( const T &  x)
inline

Returns the logarithm of \(x\) in base 2.

◆ Lg() [2/2]

double LatticeTester::Lg ( std::int64_t  x)
inline

Returns the logarithm of \(x\) in base 2.

◆ ModifVect()

template<typename Vect1 , typename Vect2 , typename Scal >
void LatticeTester::ModifVect ( Vect1 &  A,
const Vect2 &  B,
Scal  x,
int  n 
)
inline

Adds the first n components of vector B multiplied by x to first n components of vector A.

This will modify A. This does wierd type convertion and might not work well if different types are used.

◆ Modulo()

template<typename Real >
void LatticeTester::Modulo ( const Real &  a,
const Real &  b,
Real &  r 
)
inline

Computes the remainder of a/b and stores its positive equivalent mod b in r.

This works with std::int64_t, NTL::ZZ and real valued numbers.

◆ MyExit()

void LatticeTester::MyExit ( int  status,
std::string  msg 
)

Special exit function.

status is the status code to return to the system, msg is the message to print upon exit.

◆ mysqrt()

double LatticeTester::mysqrt ( double  x)
inline

Returns \(\sqrt{x}\) for \(x\ge0\), and \(-1\) for \(x < 0\).

◆ operator<<() [1/4]

template<class K , class T , class C , class A >
std::ostream& LatticeTester::operator<< ( std::ostream &  out,
const std::map< K, T, C, A > &  x 
)

Streaming operator for maps.

Formats a map as: { key1=>val1, ..., keyN=>valN }.

◆ operator<<() [2/4]

template<class T1 , class T2 >
std::ostream& LatticeTester::operator<< ( std::ostream &  out,
const std::pair< T1, T2 > &  x 
)

Streaming operator for vectors.

Formats a pair as: (first,second).

◆ operator<<() [3/4]

template<class T , class A >
std::ostream& LatticeTester::operator<< ( std::ostream &  out,
const std::vector< T, A > &  x 
)

Streaming operator for vectors.

Formats a vector as: [ val1, ..., valN ].

◆ operator<<() [4/4]

template<class K , class C , class A >
std::ostream& LatticeTester::operator<< ( std::ostream &  out,
const std::set< K, C, A > &  x 
)

Streaming operator for sets.

Formats a set as: { val1, ..., valN }.

◆ ProdScal()

template<typename Int , typename Vect1 , typename Vect2 , typename Scal >
void LatticeTester::ProdScal ( const Vect1 &  A,
const Vect2 &  B,
int  n,
Scal &  D 
)
inline

Computes the scalar product of vectors A and B truncated to their n first components, then puts the result in D.

There is a lot to consider when passing types to this function. The best is for Vect1 to be the same type as Vect2 and for Scal to be the same as Int, and that those types are the ones stored in Vect1 and Vect2.

WARNING: This uses so many types without check about them and also assumes all those types can be converted to each other without problem. This is used in many places to compute a floating point norm of vectors with integers values. Take care when using this function.

◆ Quotient()

template<typename Int >
void LatticeTester::Quotient ( const Int &  a,
const Int &  b,
Int &  q 
)
inline

Computes a/b, truncates the fractionnal part and puts the result in q.

This function is overloaded to work as specified on NTL::ZZ integers. Example:

\(a\) \(b\) \(q\)
\(5\) 3 1
\(-5\) \(3\) \(-1\)
\(5\) \(-3\) \(-1\)
\(-5\) \(-3\) \(1\)

◆ RandBits()

std::uint64_t LatticeTester::RandBits ( int  s)

Returns the first s pseudo-random bits of the underlying RNG in the form of a s-bit integer.

It is imperative that \(1 \leq s \leq 64\) because the RNG is 64 bits wide.

◆ RandInt()

int LatticeTester::RandInt ( int  i,
int  j 
)

Return a uniform pseudo-random integer in \([i, j]\).

Note that the numbers \(i\) and \(j\) are part of the possible output. It is important that \(i < j\) because the underlying arithmetic uses unsigned integers to store j-i+1 and that will be undefined behavior.

◆ RandU01()

double LatticeTester::RandU01 ( )

Returns a random number in \([0, 1)\).

The number will have 53 pseudo-random bits.

◆ Round()

template<typename Real >
Real LatticeTester::Round ( Real  x)
inline

Returns the value of x rounded to the NEAREST integer value.

(This does not truncate the integer value as is usual in computer arithmetic.)

◆ SetSeed()

void LatticeTester::SetSeed ( std::uint64_t  seed)

Sets the seed of the generator.

Because of the constraints on the state, seed has to be \( > 2\). If this is never called, a default seed will be used.

◆ SetValue()

template<typename Real >
void LatticeTester::SetValue ( Real *  A,
int  d,
const Real &  x 
)
inline

Sets the first d components of A to the value x.

◆ SetZero() [1/2]

template<typename Real >
void LatticeTester::SetZero ( Real *  A,
int  d 
)
inline

Sets the first d of A to 0.

◆ SetZero() [2/2]

template<typename Vect >
void LatticeTester::SetZero ( Vect &  A,
int  d 
)
inline

Sets the first d components of A to 0.

◆ sign()

template<typename T >
std::int64_t LatticeTester::sign ( const T &  x)
inline

Returns the sign of x.

The sign is 1 if x>0, 0 if x==0 and -1 if x<0.

◆ swap9()

template<typename T >
void LatticeTester::swap9 ( T &  x,
T &  y 
)
inline

Takes references to two variables of a generic type and swaps their content.

This uses the assignment operator, so it might not always work well if this operator's implementation is not thorough.

◆ toStr()

template<typename MatT >
std::string LatticeTester::toStr ( const MatT &  mat,
int  d1,
int  d2 
)

Returns a string that is a representation of mat.

This string will represent the \(d1 \times d2\) submatrix of the first lines and colums of mat.

◆ toString() [1/2]

template<typename Vect >
std::string LatticeTester::toString ( const Vect &  A,
int  c,
int  d,
const char *  sep = " " 
)

Returns a string containing A[c] to A[d-1] formated as [A[c]sep...sepA[d-1]].

In this string, components are separated by string sep. By default, sep is just a whitespace character.

◆ toString() [2/2]

template<typename Vect >
std::string LatticeTester::toString ( const Vect &  A,
int  d 
)

Returns a string containing the first d components of the vector A as a string.

Calls toString(const Vect&, int, int, const char*).

◆ Triangularization()

template<typename Matr , typename Int >
void LatticeTester::Triangularization ( Matr &  W,
Matr &  V,
int  lin,
int  col,
const Int &  m 
)

Takes a set of generating vectors in the matrix W and iteratively transforms it into an upper triangular lattice basis into the matrix V.

W and V have to have more lines than lin and more columns than col since this algorithm will only operate on the upper lin*col matrix of W. All the computations will be done modulo m, which means that you must know a rescalling factor for the vector system to call this function. After the execution, W will be a matrix containing irrelevant information and V will contain an upper triangular basis.

For more details please look at [5]. This algorithm basically implements what is written at the end of the article, that is the matrix W contains the set of vectors that is used and modified at each step to get a new vector from the basis.

Variable Documentation

◆ MAX_LONG_DOUBLE

const double LatticeTester::MAX_LONG_DOUBLE = 9007199254740992.0

Maximum integer that can be represented exactly as a double: \(2^{53}\).

◆ TWO_EXP

const std::int64_t LatticeTester::TWO_EXP
Initial value:
= {
1, 2, 4, 8, 16, 32, 64, 128, 256, 512, 1024, 2048, 4096, 8192, 16384, 32768,
65536, 131072, 262144, 524288, 1048576, 2097152, 4194304, 8388608, 16777216,
33554432, 67108864, 134217728, 268435456, 536870912, 1073741824, 2147483648,
4294967296, 8589934592, 17179869184, 34359738368, 68719476736, 137438953472,
274877906944, 549755813888, 1099511627776, 2199023255552, 4398046511104,
8796093022208, 17592186044416, 35184372088832, 70368744177664,
140737488355328, 281474976710656, 562949953421312, 1125899906842624,
2251799813685248, 4503599627370496, 9007199254740992, 18014398509481984,
36028797018963968, 72057594037927936, 144115188075855872, 288230376151711744,
576460752303423488, 1152921504606846976, 2305843009213693952,
(std::int64_t)4611686018427387904U, (std::int64_t)9223372036854775808U }

Table of powers of 2: TWO_EXP[ \(i\)] \(= 2^i\), \(i=0, 1, …, 63\).