Skip to main content

Overview

Giac provides comprehensive row reduction algorithms for solving linear systems, computing reduced row echelon form (RREF), and related operations. The core function is mrref with various algorithm options.

Algorithm Types

enum matrix_algorithms {
  RREF_GAUSS_JORDAN = 0,
  RREF_GUESS = 1,
  RREF_BAREISS = 2,
  RREF_MODULAR = 3,
  RREF_PADIC = 4,
  RREF_LAGRANGE = 5
};
RREF_GAUSS_JORDAN
0
Standard Gauss-Jordan elimination
RREF_GUESS
1
Automatic algorithm selection (recommended)
RREF_BAREISS
2
Bareiss algorithm (fraction-free)
RREF_MODULAR
3
Modular algorithm
RREF_PADIC
4
p-adic algorithm
RREF_LAGRANGE
5
Lagrange interpolation

Core Row Reduction

mrref

Compute reduced row echelon form with extensive options.
int mrref(const matrice & a, matrice & res, vecteur & pivots, gen & det,
          int l, int lmax, int c, int cmax,
          int fullreduction, int dont_swap_below, bool convert_internal,
          int algorithm, int rref_or_det_or_lu, GIAC_CONTEXT);

int mrref(const matrice & a, matrice & res, std::vector<int> & permutation,
          vecteur & pivots, gen & det,
          int l, int lmax, int c, int cmax,
          int fullreduction, int dont_swap_below, bool convert_internal,
          int algorithm, int rref_or_det_or_lu, GIAC_CONTEXT);

bool mrref(const matrice & a, matrice & res, vecteur & pivots, gen & det, GIAC_CONTEXT);

Parameters

a
const matrice &
Input matrix
res
matrice &
Output reduced matrix
pivots
vecteur &
Output vector of pivot positions
det
gen &
Output determinant (if computing)
l
int
Starting row (0-indexed)
lmax
int
Ending row (exclusive)
c
int
Starting column (0-indexed)
cmax
int
Ending column (exclusive)
fullreduction
int
  • 0: Reduction below diagonal only
  • non-zero: Full reduction (above and below diagonal)
dont_swap_below
int
For rows < this value, search for pivot in the row instead of column (no row swaps)
convert_internal
bool
Convert to rational fractions internally
algorithm
int
Algorithm selection (see matrix_algorithms enum)
rref_or_det_or_lu
int
  • 0: RREF computation
  • 1: Determinant computation
  • 2: LU decomposition
  • 3: LU without permutation
return
int
  • 0: Failure
  • 1: Success
  • 2: Success with inversion (no need to remove identity)

Simple Usage

matrice A = ...; // Input matrix
matrice R;
vecteur pivots;
gen det;

bool success = mrref(A, R, pivots, det, contextptr);
Giac Command: _rref

modrref

Modular row reduction.
bool modrref(const matrice & a, matrice & res, vecteur & pivots, gen & det,
             int l, int lmax, int c, int cmax,
             int fullreduction, int dont_swap_below, 
             const gen & modulo, bool ckprime, int rref_or_det_or_lu);

bool modrref(const matrice & a, matrice & res, vecteur & pivots, gen & det,
             const gen & modulo);
modulo
const gen &
Modulo value for arithmetic
ckprime
bool
Check if modulo is prime
return
bool
true if reduction successful

Specialized Row Reduction

smallmodrref

Efficient integer modular RREF for small moduli.
void smallmodrref(int nthreads, std::vector<std::vector<int>> & N,
                  vecteur & pivots, std::vector<int> & permutation,
                  std::vector<int> & maxrankcols, longlong & idet,
                  int l, int lmax, int c, int cmax,
                  int fullreduction, int dont_swap_below, int modulo,
                  int rref_or_det_or_lu, bool reset,
                  smallmodrref_temp_t * workptr, bool allow_block, int carac);
nthreads
int
Number of threads for parallel computation
N
std::vector<std::vector<int>> &
Input/output integer matrix
modulo
int
Modulo value (must be < 2^31)
allow_block
bool
Enable block algorithm optimization

doublerref

Floating-point RREF.
void doublerref(matrix_double & N, vecteur & pivots,
                std::vector<int> & permutation, std::vector<int> & maxrankcols,
                double & idet, int l, int lmax, int c, int cmax,
                int fullreduction, int dont_swap_below,
                int rref_or_det_or_lu, double eps);
N
matrix_double &
Input/output double-precision matrix
eps
double
Tolerance for zero detection

Linear Combinations

linear_combination

Compute linear combination of two vectors.
void linear_combination(const gen & c1, const vecteur & v1,
                       const gen & c2, const vecteur & v2,
                       const gen & c, const gen & cinv,
                       vecteur & v, double eps, int cstart = 0);
c1
const gen &
Coefficient for first vector
v1
const vecteur &
First vector
c2
const gen &
Coefficient for second vector
v2
const vecteur &
Second vector
c
const gen &
Scaling factor
cinv
const gen &
Inverse of c (for optimization)
v
vecteur &
Output result vector
eps
double
Tolerance: values with ||v|| < eps are replaced by 0
cstart
int
default:"0"
Starting index for operation

modlinear_combination

Modular linear combination.
void modlinear_combination(vecteur & v1, const gen & c2, const vecteur & v2,
                          const gen & modulo, int cstart, int cend = 0);

void modlinear_combination(std::vector<int> & v1, int c2,
                          const std::vector<int> & v2, int modulo,
                          int cstart, int cend, bool pseudo);
v1
vecteur &
First vector (modified in place)
c2
const gen & | int
Coefficient for second vector
v2
const vecteur &
Second vector
modulo
const gen & | int
Modulo value
cstart
int
Starting index
cend
int
default:"0"
Ending index (0 means end of vector)

System Solving

p-adic Linear System Solving

padic_linsolve

Solve a linear system using p-adic methods.
int padic_linsolve(const matrice & a, const vecteur & b, vecteur & res,
                   gen & p, gen & det_mod_p, gen & h2,
                   unsigned reconstruct = 0, int maxtry = 4);
a
const matrice &
Coefficient matrix (integer entries)
b
const vecteur &
Right-hand side vector (integer entries)
res
vecteur &
Output solution vector
p
gen &
Prime for p-adic computation
det_mod_p
gen &
Determinant mod p
h2
gen &
Auxiliary parameter
reconstruct
unsigned
default:"0"
Rational reconstruction parameter
maxtry
int
default:"4"
Maximum number of attempts
return
int
  • 0: No invertible element found
  • -1: Determinant is zero (no solution)
  • 1: Success
Giac Command: _padic_linsolve

padic_linsolve_prepare

Prepare for solving non-Cramer systems.
int padic_linsolve_prepare(const matrice & a, gen & p,
                           std::vector<int> & ranklines,
                           std::vector<int> & rankcols,
                           matrice & asub, matrice & ainv,
                           vecteur & compat, vecteur & kernel);
ranklines
std::vector<int> &
Output: row indices for maximal rank submatrix
rankcols
std::vector<int> &
Output: column indices for maximal rank submatrix
asub
matrice &
Output: maximal rank submatrix
ainv
matrice &
Output: inverse of asub mod p
compat
vecteur &
Output: compatibility conditions
kernel
vecteur &
Output: kernel basis
return
int
Rank of matrix, or -1 on failure

padic_linsolve_solve

Solve using prepared system.
bool padic_linsolve_solve(const matrice & a, const gen & p,
                          const std::vector<int> & ranklines,
                          const std::vector<int> & rankcols,
                          const matrice & asub, const matrice & ainv,
                          const vecteur & compat, const vecteur & b,
                          vecteur & sol);

Indefinite System Solving

bool solve_indef(matrice & A, const vecteur * b, vecteur & x,
                 int * p, int * n, int * z, GIAC_CONTEXT);
A
matrice &
System matrix (LDL factorization)
b
const vecteur *
Right-hand side
x
vecteur &
Output solution
p
int *
Number of positive eigenvalues
n
int *
Number of negative eigenvalues
z
int *
Number of zero eigenvalues
return
bool
true if solve successful

Utility Functions

mdividebypivot

Divide matrix rows by their pivot elements.
void mdividebypivot(matrice & a, int lastcol = -1);
void mdividebypivot(matrice & a, int lastcol, GIAC_CONTEXT);
a
matrice &
Matrix to modify (in-place)
lastcol
int
default:"-1"
  • -1: Divide last column
  • -2: Do not divide last column
  • ≥0: Stop dividing at this column

add_identity

Append identity matrix to the right.
void add_identity(matrice & arref);
void add_identity(std::vector<std::vector<int>> & arref);
arref
matrice &
Matrix to augment (modified in-place)

remove_identity

Remove identity matrix from augmented result.
bool remove_identity(matrice & res);
bool remove_identity(std::vector<std::vector<int>> & res, int modulo);
res
matrice &
Augmented matrix (modified in-place)
return
bool
true if identity was found and removed

Quadratic Forms (gauss.h)

gauss

Gaussian reduction of quadratic forms.
vecteur gauss(const gen & q, const vecteur & x, GIAC_CONTEXT);
vecteur gauss(const gen & q, const vecteur & x, vecteur & D,
              vecteur & U, vecteur & P, GIAC_CONTEXT);
q
const gen &
Quadratic form (symbolic expression)
x
const vecteur &
Vector of variables
D
vecteur &
Output diagonal coefficients
U
vecteur &
Output transformation matrix
P
vecteur &
Output change of variables
return
vecteur
Reduced quadratic form data
Giac Command: _gauss

quad

Find matrix of a quadratic form.
vecteur quad(int & b, const gen & q, const vecteur & x, GIAC_CONTEXT);
b
int &
Output: form degree (2 if purely quadratic, 0/1/3 if contains terms of other degrees)
q
const gen &
Symbolic quadratic expression
x
const vecteur &
Variables
return
vecteur
Matrix representation of quadratic form

q2a and a2q

Convert between quadratic forms and matrices.
gen _q2a(const gen & args, GIAC_CONTEXT);  // quadratic form to matrix
gen _a2q(const gen & args, GIAC_CONTEXT);  // matrix to quadratic form
Giac Commands: _q2a, _a2q

Conic Sections

conique_reduite

Reduce a conic equation to standard form.
int conique_reduite(const gen & equation_conique,
                    const gen & pointsurconique,
                    const vecteur & nom_des_variables,
                    gen & x0, gen & y0, vecteur & V1, vecteur & V2,
                    gen & propre, gen & equation_reduite,
                    vecteur & param_curves, gen & ratparam,
                    bool numeric, GIAC_CONTEXT,
                    gen * aptr = 0, gen * bptr = 0);
return
int
  • 0: Not a conic
  • Less than 0: Not a proper conic
  • 2: Parabola
  • 3: Ellipse
  • 4: Hyperbola
Giac Command: _conique_reduite

See Also

Build docs developers (and LLMs) love