RISC JKU

nullspace

This module is about computing the right kernel of matrices with polynomial entries. It provides various general purpose tools for constructing solvers.

By a solver, we mean a function of the following form:

nullspace.solver(mat, degrees=[], infolevel=0)

INPUT:

  • mat – a matrix \(A\) over some polynomial ring
  • degrees – list of nonnegative integers, providing, for each variable \(x\), a bound on the degree with which \(x\) is going to appear in the output. This data is optional. If set to [], the solver has to figure out the degrees by itself.
  • infolevel – a nonnegative integer indicating the desired verbosity (default=0), or a pair (u, v) where u indicates the desired verbosity and v specifies how many leading spaces should be added to the printed lines (default=0).

OUTPUT:

  • list of vectors \(v\) such that \(A*v=0\)

Depending on the application from which a polynomial matrix arises, the preferred way of computing its nullspace may be different. The purpose of this module is not to provide a single solver which is good for every possible input, but instead, it provides a collection of functions by which solvers can be composed.

For example, gauss is a function for constructing a solver using fraction free Gaussian elimination. It works for most polynomial domains.

sage: my_solver = gauss()
sage: A = MatrixSpace(ZZ['x'], 4, 5).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

sage: A = MatrixSpace(ZZ['x', 'y'], 4, 5).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

sage: A = MatrixSpace(GF(1093)['x', 'y'], 4, 5).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

Several other functions create solvers which reduce the nullspace problem to one or more nullspace problems over simpler domains which are then solved by another solver. This solver can be chosen by the user.

For example, kronecker is a function for constructing a solver for matrices of multivariate polynomials. It requires as parameter a solver for matrices of univariate polynomials, for instance a solver created by gauss.

sage: my_solver = kronecker(gauss())
sage: A = MatrixSpace(ZZ['x', 'y'], 4, 5).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

For large random examples, the solver created by kronecker(gauss()) is likely to be significantly faster than the solver created by gauss(). Even faster, for matrices of polynomials with integer coefficients, is a solver using chinese remaindering combined with kronecker substitution combined with gaussian elimination:

sage: my_solver = cra(kronecker(gauss()))
sage: A = MatrixSpace(ZZ['x', 'y'], 4, 5).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

Alternatively:

sage: my_solver = kronecker(cra(gauss()))
sage: A = MatrixSpace(ZZ['x', 'y'], 4, 5).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

Here is the same example with a variant of the same solver (not necessarily faster).

sage: my_solver = cra(kronecker(lagrange(sage_native, start_point=5), newton(sage_native)), max_modulus = 2^16, proof=True, ncpus=4)
sage: A = MatrixSpace(ZZ['x', 'y'], 4, 5).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

A particular solver will typically only be applicable to matrices with entries in a domain for which it was designed for. Below is an overview of the functions provided in this module, together with the supported input domains, and the input domains of the subsolvers they require. The notation \(K[x,...]\) refers to univariate or multivariate polynomial ring, understanding the same reading (unvariate vs. multivariate) in corresponding rows of the 2nd and 3rd column.

method input domain requires subsolver for
cra \(K[x,...]\) where \(K\) is \(ZZ\), \(QQ\), or \(GF(p)\) \(GF(p)[x,...]\)
galois \(QQ(alpha)[x,...]\) \(GF(p)[x,...]\)
clear \(K(x,...)\) \(K[x,...]\)
clear \(K[x,...]\) where \(K\) is the fraction field of \(R\) \(R[x,...]\)
compress \(K[x,...]\) or \(K(x,...)\) same domain and \(GF(p)\)
kronecker \(K[x,...]\) \(K[x]\) and \(GF(p)[x]\)
gauss \(K[x,...]\) None
wiedemann \(K[x,...]\) or \(K(x,...)\) None
lagrange \(K[x]\) or \(K(x)\) where \(K\) is a field \(K\)
hermite \(K[x]\) where \(K\) is a field None
newton \(K[x]\) where \(K\) is a field \(K\)
merge \(K[x,...][y,...]\) \(K[x,...,y,...]\)
quick_check \(K[x,...]\) where \(K\) is \(ZZ\), \(QQ\), or \(GF(p)\) same domain and \(GF(p)\)
sage_native \(K[x,...]\) or \(K(x,...)\) or \(K\) None

AUTHOR:

  • Manuel Kauers (2012-09-16)
nullspace.clear(subsolver)

Constructs a solver which clears denominators in a given matrix over \(FF(R)[x..]\) or \(FF(R[x..])\) and then calls a subsolver on a matrix over \(R[x..]\)

INPUT:

  • subsolver – a solver for matrices over \(R[x..]\)

OUTPUT:

  • a solver for matrices over \(FF(R)[x..]\) or \(FF(R[x..])\)

EXAMPLE:

sage: A = MatrixSpace(ZZ['x'].fraction_field(), 4, 5).random_element()
sage: my_solver = clear(gauss())
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

sage: A = MatrixSpace(QQ['x'], 4, 5).random_element()
sage: my_solver = clear(gauss())
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

For turning a matrix over \(QQ(x)\) to one over \(ZZ[x]\), you need to clear twice: once the denominator of the rational functions, and then once more the denominators of the coefficients

sage: A = MatrixSpace(QQ['x'].fraction_field(), 4, 5).random_element()
sage: my_solver = clear(clear(gauss()))
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM: clears denominators and then applies the subsolver.

Warning

For unstructured matrices, clearing denominators may significantly increase the size of the system. In such situations, consider using nullspace.lagrange .

nullspace.compress(subsolver, presolver=<function sage_native at 0xde49b54>, modulus=8388608)

Constructs a solver which throws away unnecessary rows and columns and then applies the subsolver

INPUT:

  • subsolver – a solver for matrices over \(R[x..]\)
  • presolver – a solver for matrices over \(GF(p)\), defaults to sage_native
  • modulus – modulus used for the precomputation in the homomorphic image if \(R.characteristic()==0\). If this number is not a prime, it will be replaced by the next smaller integer which is prime.

OUTPUT:

  • a solver for matrices over \(R[x..]\)

EXAMPLE:

sage: A = MatrixSpace(ZZ['x'], 7, 4).random_element()*MatrixSpace(ZZ['x'], 4, 5).random_element()
sage: my_solver = compress(gauss())
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0, 0, 0, 0)

ALGORITHM:

  1. Compute the homomorphic image of the solution space using the presolver. Discard all those columns of the matrix where the corresponding row of every solution vector contains a zero.
  2. Starting from the heaviest row down to the lightest (measuring the weight of a row by counting the number of nonzero entries and the number of monomials), check whether the solution space changes if the row is dropped. Keep on dropping rows until the difference between the number of rows and the number of columns agrees with the dimension of the solution space.
  3. Apply the subsolver to the resulting matrix. Fill in zeros into the output according to the columns deleted in step one. Then return the result.
nullspace.cra(subsolver, max_modulus=8388608, proof=False, ncpus=1)

Creates a subsolver based on chinese remaindering for matrices over \(K[x]\) or \(K[x,y,..]\) where \(K\) is \(ZZ\) or \(QQ\) or \(GF(p)\).

INPUT:

  • subsolver – a solver for matrices over \(GF(p)[x]\) (if the original matrix contains univariate polynomails) or \(GF(p)[x,y,...]\) (if it contains multivariate polynomials)
  • max_modulus – a positive integer. The solver will iterate over the primes less than this number in decreasing order. Defaults to the largest word size integer for which we expect Sage to use hardware arithmetic.
  • proof – a boolean value. If set to False (default), a termination is only tested in a homomorphic image, which saves much time but may, with a very low probability, lead to a wrong output.
  • ncpus – number of cpus that may be used in parallel by the solver (default=1).

OUTPUT:

  • a solver for matrices with entries in \(K[x,...]\) based on a subsolver for matrices with entries in \(GF(p)[x,...]\).

EXAMPLES:

sage: A = MatrixSpace(ZZ['x', 'y'], 4, 7).random_element(degree=3)
sage: my_solver = cra(kronecker(gauss()))
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM:

  1. If the coefficient domain is a finite field, the problem is delegated to the subsolver and we return whatever we obtain from it.
  2. For various word-size primes \(p\), reduce the coefficients in mat modulo \(p\) and call the subsolver on the resulting matrix.
  3. Using the Chinese Reminder Algorithm and rational reconstruction, combine these results to a solution candidate with integer coefficients.
  4. If this solution candidate is correct, return it and stop. If proof is set to True, this check is performed rigorously, otherwise (default) only modulo some prime.
  5. If the solution candidate is not correct, consider some more primes and try again.
nullspace.galois(subsolver, max_modulus=8388608, proof=False)

Creates a subsolver based on chinese remaindering for matrices over \(K[x]\) or \(K[x,y,..]\) where \(K\) is a single algebraic extension of \(QQ\)

INPUT:

  • subsolver – a solver for matrices over \(GF(p)[x]\) (if the original matrix contains univariate polynomails) or \(GF(p)[x,y,...]\) (if it contains multivariate polynomials)
  • max_modulus – a positive integer. The solver will iterate over the primes less than this number in decreasing order. Defaults to the largest word size integer for which we expect Sage to use hardware arithmetic.
  • proof – a boolean value. If set to False (default), a termination is only tested in a homomorphic image, which saves much time but may, with a very low probability, lead to a wrong output.

OUTPUT:

  • a solver for matrices with entries in \(QQ(\alpha)[x,...]\) based on a subsolver for matrices with entries in \(GF(p)[x,...]\).

EXAMPLES:

sage: R.<x> = QQ['x']
sage: K.<a> = NumberField(x^3-2, 'a')
sage: A = MatrixSpace(K['x'], 4, 5).random_element(degree=3)
sage: my_solver = galois(gauss())
sage: ##V = my_solver(A) 
sage: ##A*V[0]
##(0, 0, 0, 0)

ALGORITHM:

  1. If the coefficient domain is a finite field, the problem is delegated to the subsolver and we return whatever we obtain from it.
  2. For various word-size primes \(p\), reduce the coefficients in mat modulo \(p\) and call the subsolver on the resulting matrix. Primes are chosen in such a way that the minimal polynomial of the generator of the number field splits into linear factors modulo \(p\).
  3. The generator of the number field is mapped in turn to each of the roots of the linear factors of the minimal polynomial mod \(p\), the results are then combined by interpolation to a polynomial over \(GF(p)\), and chinese remaindering and rational reconstruction are used to lift it back to the original number field.
  4. If this solution candidate is correct, return it and stop. If proof is set to True, this check is performed rigorously, otherwise (default) only modulo some new prime.
  5. If the solution candidate is not correct, consider some more primes and try again.
nullspace.gauss(pivot=<function _pivot at 0xde499cc>, ncpus=1, fun=None)

Creates a solver based on fraction free gaussian elimination.

INPUT:

  • pivot – a function which takes as input: a matrix mat (as list of list of ring elements), four integers \(r\), \(n\), \(c\), \(m\) specifying the index range mat[r:n][c:m], and the zero element of the ring, and returns and a pair \((i, j)\) such that mat[i][j]!=zero and \((i, j)\) is a good choice for a pivot. The function returns None iff there are no nonzero elements in the given range of the matrix.

    The default function chooses \((i,j)\) such that mat has many zeros in row \(i\) and column \(j\), and mat[i][j] has a small number of terms.

  • ncpus – maximum number of cpus that may be used in parallel by this solver

  • fun – if different from None, at the beginning of each iteration of the outer loops of forward and backward elimination, the solver calls fun(mat, idx), where mat is the current matrix (as list of lists of elements) and idx is a counter. This functionality is intended for analyzing the elimination process, e.g., for inspecting how well the solver succees in maintaining the sparsity of the matrix. The solver assumes that the function won’t modify the matrix.

OUTPUT:

  • A solver based on fraction free gaussian elimination.

EXAMPLES:

sage: A = MatrixSpace(GF(1093)['x'], 4, 7).random_element(degree=3)
sage: my_solver = gauss()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM: fraction-free gaussian elimination with heuristic content removal and Markoviz pivot search.

nullspace.hermite(early_termination=True)

Creates a solver which computes a nullspace basis of minimal degree.

INPUT:

  • early_termination – a boolean value. If set to True (default), the calculation is aborted as soon as the first solution vector has been found. If set to False, the calculation continues up to some (potentially rather pessimistic) worst case bound on the possible degrees of the solution vectors. If degree information is supplied to the solver, the early_termination setting is ignored.

OUTPUT:

  • A solver based on Hermite-Pade approximation

EXAMPLES:

sage: A = MatrixSpace(GF(1093)['x'], 4, 7).random_element(degree=3)
sage: my_solver = hermite()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM: Hermite-Pade approximation

nullspace.kronecker(subsolver, presolver=None)

Creates a solver for matrices of multivariate polynomials over some domain \(K\), based on a given solver for matrices over univariate polynomials

INPUT:

  • subsolver – a solver for univariate polynomial matrices over \(K\)
  • presolver – a solver for univariate polynomial matrices over prime fields. If None is given, the presolver will be set to subsolver

OUTPUT:

  • a solver for matrices with entries in \(K[x,y,...]\) where K is some domain such that the given subsolver can solve matrices in \(K[x]\). If the solver is called without degree information about the solution vectors, the presolver is called on various matrices with entries in \(GF(p)[x]\) to determine the degrees. In this case if \(K\) is not a prime field, its elements must allow for coercion to prime field elements.

EXAMPLE:

sage: A = MatrixSpace(GF(1093)['x','y'], 4, 7).random_element(degree=3)
sage: mysolver = kronecker(gauss())
sage: V = mysolver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM:

  1. If applied to a matrix of univariate polynomials, the function will delegate the whole problem to the subsolver, and return its result.
  2. When no degree information is specified, the presolver will be applied to homomorphic images of the matrix with all the variables except one set to some ground field element and reduce the coefficients modulo a word size prime (unless the coefficient field already is a prime field). The degrees seen in these result are taken as the degrees of the respective variables in the final result.
  3. By Kronecker-substitution, reduce the task to a nullspace-problem for a matrix over K[x], apply the subsolver to it, undo the Kronecker-substitution, and return the result.
nullspace.lagrange(subsolver, start_point=10, ncpus=1)

Creates a solver for matrices of univariate polynomials or rational functions over some field \(K\), based on a given solver for matrices over \(K\), using evaluation+interpolation.

INPUT:

  • subsolver – a solver for matrices over \(K\)
  • start_point – first evaluation point to be used
  • ncpus – maximum number of cpus that may be used in parallel by the solver (default=1).

OUTPUT:

  • a solver for matrices with entries in \(K[x]\) or \(K(x)\) where K is some field such that the given subsolver can solve matrices with entries in \(K\).

EXAMPLE:

sage: A = MatrixSpace(GF(1093)['x'], 4, 7).random_element(degree=3)
sage: my_solver = lagrange(sage_native)
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM:

  1. For x replaced by start_point, start_point+1, start_point+2, ..., compute the nullspace by the given subsolver.
  2. Use (fast) interpolation and rational reconstruction to combine the images to a nullspace basis for the original matrix over K[x]

Note

  • The solver raises a ValueError if the ground field is too small.
  • It is assumed that the subsolver returns a normal form of the nullspace basis, so that results of calls for different evaluations have some common preimage in the matrices over \(K(x)\)
nullspace.merge(subsolver)

Constructs a solver which first merges towers of polynomial or rational function extensions into a single one and then applies the subsolver.

INPUT:

  • subsolver – a solver for matrices over the target ring

OUTPUT:

  • a solver for matrices over base rings which are obtained from some ring \(R\) by extending twice by polynomials or rational functions, i.e., \(R=B[x..][y..]\) or \(R=F(B[x..])[y..]\) or \(R=F(B[x..][y..])\) or \(R=F(F(B[x..])[y..])\). In the first case, the target ring will be \(B[x..,y..]\), in the remaining cases it will be \(F(B[x..,y..])\). If \(R\) is not of one of the four forms listed above, the matrix is passed to the subsolver without modification.

EXAMPLES:

sage: my_solver = merge(kronecker(gauss()))
sage: A = MatrixSpace(ZZ['x']['y'], 3, 4).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0)
sage: my_solver = merge(clear(kronecker(gauss())))
sage: A = MatrixSpace(ZZ['x'].fraction_field()['y'], 3, 4).random_element()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0)
nullspace.newton(subsolver, inverse=<function <lambda> at 0xde49f7c>)

Constructs a solver based on x-adic lifting for matrices with entries over \(K[x]\) where \(K\) is a field.

INPUT:

  • subsolver – a solver for matrices over \(K\)
  • inverse – a function for computing the inverse of a given regular square matrix with entries in \(K\). Defaults to lambda mat : mat.inverse()

OUTPUT:

  • a solver for matrices with entries in \(K[x]\) based on a subsolver for matrices with entries in \(K\).

EXAMPLES:

sage: A = MatrixSpace(GF(1093)['x'], 4, 7).random_element(degree=3)
sage: my_solver = newton(sage_native)
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM:

  1. First call the subsolver on a homomorphic image to get an idea about the rank of the matrix and the shape of the reduced echolon form of the solution matrix.
  2. Based on this information, rewrite the original homogeneous system \(A\cdot X=0\) into an inhomogeneous system \(U\cdot Y+V=0\) where \(U\) is invertible.
  3. Compute the inverse \(U(x=0)^{-1}\) of \(U(x=0)\) using the given inverse function.
  4. Using Newton-iteration, lift \(U(x=0)^{-1}\) to a polynomial matrix \(W\) with \(U\cdot W=1\bmod x^{2k}\) where \(k\) exceeds the specified degree bound (or a worst-case bound, if no degree bound was specified).
  5. Compute \(W\cdot V\), rewrite it as nullspace basis mod \(x^{2k}\), apply rational reconstruction, clear denominators, and return the result.

Note

It is assumed that the subsolver returns a nullspace basis which is normalized so that when the basis vectors are viewed as the columns of a matrix, this matrix contains every possible unit vector as row

nullspace.quick_check(subsolver, modsolver=<function sage_native at 0xde49b54>, modulus=8388608)

Constructs a solver which first tests in a homomorphic image whether the nullspace is empty and applies the subsolver only if it is not.

INPUT:

  • subsolver – a solver
  • modsolver – a solver for matrices over \(GF(p)\), defaults to sage_native
  • modulus – modulus used for the precomputation in the homomorphic image if \(R.characteristic()==0\). If this number is not a prime, it will be replaced by the next smaller integer which is a prime.

OUTPUT:

  • a solver for matrices over \(K[x..]\) or \(K(x..)\) with \(K\) one of \(ZZ\), \(QQ\), \(GF(p)\), and which can be handled by the given subsolver.

EXAMPLE:

sage: A = MatrixSpace(ZZ['x'], 4, 5).random_element()
sage: my_solver = quick_check(gauss())
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)
sage: V = my_solver(A.transpose())
sage: len(V)
0
nullspace.sage_native(mat, degrees=[], infolevel=0)

Computes the nullspace of the given matrix using Sage’s native method.

INPUT:

  • mat – any matrix
  • degrees – ignored
  • infolevel – a nonnegative integer indicating the desired verbosity.

OUTPUT:

  • a list of vectors that form a basis of the right kernel of mat

EXAMPLES:

sage: A = MatrixSpace(GF(1093)['x'], 4, 7).random_element(degree=3)
sage: V = sage_native(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM: just calls right_kernel_matrix() on mat and converts its output to a list of vectors

nullspace.wiedemann()

Constructs a solver using Wiedemann’s algorithm

INPUT:

  • none.

OUTPUT:

  • a solver for matrices over \(R[x..]\)

EXAMPLE:

sage: A = MatrixSpace(ZZ['x'], 4, 5).random_element()
sage: my_solver = wiedemann()
sage: V = my_solver(A)
sage: A*V[0]
(0, 0, 0, 0)

ALGORITHM: Wiedemann’s algorithm.

Note

  1. If the matrix is not a square matrix, it will be left-multiplied by its transpose before the algorithm is started.
  2. If the matrix is a square matrix, it need not actually be a Sage matrix object. Instead, it can be any Python object A which provides the following functions: A.dimensions(), A.parent(), and A*v for a vector v.
  3. The solver returns at most one solution vector, even if the nullspace has higher dimension. If it returns the empty list, it only means that the nullspace is probably empty.

Previous topic

generalized_series