Associative and central simple algebras

This section collects functions related to associative algebras and central simple algebras (CSA) over number fields.


Algebra definitions

Let A be a finite-dimensional unital associative algebra over a field K. The algebra A is central if its center is K and it is simple if it has no nontrivial two-sided ideals.

We provide functions to handle associative algebras of finite dimension over ℚ or 𝔽p. We represent them by the left multiplication table on a basis over the prime subfield; the function algtableinit creates the object representing an associative algebra. We also provide functions to handle central simple algebras over a number field K. We represent them either by the left multiplication table on a basis over the center K or by a cyclic algebra (see below); the function alginit creates the object representing a central simple algebra.

The set of elements of an algebra A that annihilate every simple left A-module is a two-sided ideal, called the Jacobson radical of A. If the Jacobson radical is trivial, the algebra is semisimple: it is isomorphic to a direct product of simple algebras. The dimension of a CSA over its center K is always a square d^2 and the integer d is called the degree of the algebra over K. A CSA over a field K is always isomorphic to Mk(D) for some integer k and some central division algebra D of degree e: the integer e is the index of the algebra.

Let L/K be a cyclic extension of degree d, let σ be a generator of Gal(L/K) and let b ∈ K*. Then the cyclic algebra (L/K,σ,b) is the algebra ⨁ i = 0d-1x^iL with x^d = b and ℓ x = xσ(ℓ) for all ℓ ∈ L. The algebra (L/K,σ,b) is a central simple K-algebra of degree d, and it is an L-vector space. Left multiplication is L-linear and induces a K-algebra isomorphism (L/K,σ,b) ⨂ K L → Md(L).

Let K be a nonarchimedean local field with uniformizer π, and let L/K be the unique unramified extension of degree d. Then every central simple algebra A of degree d over K is isomorphic to (L/K, Frob, π^h) for some integer h. The element h/d ∈ ℚ/ℤ is called the Hasse invariant of A.


Orders in algebras

Let A be an algebra of finite dimension over ℚ. An order in A is a finitely generated ℤ-submodule 𝒪 such that ℚ𝒪 = A, that is also a subring with unit. By default the data computed by alginit contains a ℤ-basis of a maximal order 𝒪0. We define natural orders in central simple algebras defined by a cyclic algebra or by a multiplication table over the center. Let A = (L/K,σ,b) = ⨁ i = 0d-1x^iL be a cyclic algebra over a number field K of degree n with ring of integers ℤK. Let ℤL be the ring of integers of L, and assume that b is integral. Then the submodule 𝒪 = ⨁ i = 0d-1x^iℤL is an order in A, called the natural order. Let ω0,...,ωnd-1 be a ℤ-basis of ℤL. The natural basis of 𝒪 is b0,...,bnd^2-1 where bi = xi/(nd)ω(i mod nd). Now let A be a central simple algebra of degree d over a number field K of degree n with ring of integers ℤK. Let e0,...,ed^2-1 be a basis of A over K and assume that the left multiplication table of A on (ei) is integral. Then the submodule 𝒪 = ⨁ i = 0d^2-1K ei is an order in A, called the natural order. Let ω0,...,ωn-1 be a ℤ-basis of ℤK. The natural basis of 𝒪 is b0,...,bnd^2-1 where bi = ω(i mod n)ei/n.


Lattices in algebras

We also provide functions to handle full lattices in algebras over ℚ. A full lattice J ⊂ A is represented by a 2-component t_VEC [I,t] representing J = tI, where

* I is an integral nonsingular upper-triangular matrix representing a sublattice of 𝒪0 expressed on the integral basis, and

* t ∈ ℚ > 0 is a t_INT or t_FRAC.

For the sake of efficiency you should use matrices I that are primitive and in Hermite Normal Form; this makes the representation unique. No GP function uses this property, but all GP functions return lattices in this form. The prefix for lattice functions is alglat.


GP conventions for algebras

As with number fields, we represent elements of central simple algebras in two ways, called the algebraic representation and the basis representation, and you can convert betweeen the two with the functions algalgtobasis and algbasistoalg. In every central simple algebra object, we store a ℤ-basis of an order 𝒪0, and the basis representation is simply a t_COL with coefficients in ℚ expressing the element in that basis. If no maximal order was computed by alginit, then 𝒪0 is the natural order. If a maximal order was computed, then 𝒪0 is a maximal order containing the natural order. For a cyclic algebra A = (L/K,σ,b), the algebraic representation is a t_COL with coefficients in L representing the element in the decomposition A = ⨁ i = 0d-1x^iL. For a central simple algebra defined by a multiplication table over its center K on a basis (ei), the algebraic representation is a t_COL with coefficients in K representing the element on the basis (ei).

Warning. The coefficients in the decomposition A = ⨁ i = 0d-1x^iL are not the same as those in the decomposition A = ⨁ i = 0d-1Lx^i! The i-th coefficients are related by conjugating by x^i, which on L amounts to acting by σ^i.

Warning. For a central simple algebra over ℚ defined by a multiplication table, we cannot distinguish between the basis and the algebraic representations from the size of the vectors. The behavior is then to always interpret the column vector as a basis representation if the coefficients are t_INT or t_FRAC, and as an algebraic representation if the coefficients are t_POL or t_POLMOD.


algadd(al, x, y)

Given two elements x and y in al, computes their sum x+y in the algebra al.

  ? A = alginit(nfinit(y),[-1,1]);
  ? algadd(A,[1,0]~,[1,2]~)
  %2 = [2, 2]~

Also accepts matrices with coefficients in al.

The library syntax is GEN algadd(GEN al, GEN x, GEN y).


algalgtobasis(al, x)

Given an element x in the central simple algebra al output by alginit, transforms it to a column vector on the integral basis of al. This is the inverse function of algbasistoalg.

  ? A = alginit(nfinit(y^2-5),[2,y]);
  ? algalgtobasis(A,[y,1]~)
  %2 = [0, 2, 0, -1, 2, 0, 0, 0]~
  ? algbasistoalg(A,algalgtobasis(A,[y,1]~))
  %3 = [Mod(Mod(y, y^2 - 5), x^2 - 2), 1]~

The library syntax is GEN algalgtobasis(GEN al, GEN x).


algaut(al)

Given a cyclic algebra al = (L/K,σ,b) output by alginit, returns the automorphism σ.

  ? nf = nfinit(y);
  ? p = idealprimedec(nf,7)[1];
  ? p2 = idealprimedec(nf,11)[1];
  ? A = alginit(nf,[3,[[p,p2],[1/3,2/3]],[0]]);
  ? algaut(A)
  %5 = -1/3*x^2 + 1/3*x + 26/3

The library syntax is GEN algaut(GEN al).


algb(al)

Given a cyclic algebra al = (L/K,σ,b) output by alginit, returns the element b ∈ K.

  nf = nfinit(y);
  ? p = idealprimedec(nf,7)[1];
  ? p2 = idealprimedec(nf,11)[1];
  ? A = alginit(nf,[3,[[p,p2],[1/3,2/3]],[0]]);
  ? algb(A)
  %5 = Mod(-77, y)

The library syntax is GEN algb(GEN al).


algbasis(al)

Given a central simple algebra al output by alginit, returns a ℤ-basis of the order 𝒪0 stored in al with respect to the natural order in al. It is a maximal order if one has been computed.

  A = alginit(nfinit(y), [-1,-1]);
  ? algbasis(A)
  %2 =
  [1 0 0 1/2]
  
  [0 1 0 1/2]
  
  [0 0 1 1/2]
  
  [0 0 0 1/2]

The library syntax is GEN algbasis(GEN al).


algbasistoalg(al, x)

Given an element x in the central simple algebra al output by alginit, transforms it to its algebraic representation in al. This is the inverse function of algalgtobasis.

  ? A = alginit(nfinit(y^2-5),[2,y]);
  ? z = algbasistoalg(A,[0,1,0,0,2,-3,0,0]~);
  ? liftall(z)
  %3 = [(-1/2*y - 2)*x + (-1/4*y + 5/4), -3/4*y + 7/4]~
  ? algalgtobasis(A,z)
  %4 = [0, 1, 0, 0, 2, -3, 0, 0]~

The library syntax is GEN algbasistoalg(GEN al, GEN x).


algcenter(al)

If al is a table algebra output by algtableinit, returns a basis of the center of the algebra al over its prime field (ℚ or 𝔽p). If al is a central simple algebra output by alginit, returns the center of al, which is stored in al.

A simple example: the 2 x 2 upper triangular matrices over ℚ, generated by I2, a = [0,1;0,0] and b = [0,0;0,1], such that a^2 = 0, ab = a, ba = 0, b^2 = b: the diagonal matrices form the center.

  ? mt = [matid(3),[0,0,0;1,0,1;0,0,0],[0,0,0;0,0,0;1,0,1]];
  ? A = algtableinit(mt);
  ? algcenter(A) \\ = (I2)
  %3 =
  [1]
  
  [0]
  
  [0]

An example in the central simple case:

  ? nf = nfinit(y^3-y+1);
  ? A = alginit(nf, [-1,-1]);
  ? algcenter(A).pol
  %3 = y^3 - y + 1

The library syntax is GEN algcenter(GEN al).


algcentralproj(al, z, {maps = 0})

Given a table algebra al output by algtableinit and a t_VEC z = [z1,...,zn] of orthogonal central idempotents, returns a t_VEC [al1,...,aln] of algebras such that ali = zi al. If maps = 1, each ali is a t_VEC [quo,proj,lift] where quo is the quotient algebra, proj is a t_MAT representing the projection onto this quotient and lift is a t_MAT representing a lift.

A simple example: 𝔽2 x 𝔽4, generated by 1 = (1,1), e = (1,0) and x such that x^2+x+1 = 0. We have e^2 = e, x^2 = x+1 and ex = 0.

  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,2);
  ? e = [0,1,0]~;
  ? e2 = algsub(A,[1,0,0]~,e);
  ? [a,a2] = algcentralproj(A,[e,e2]);
  ? algdim(a)
  %6 = 1
  ? algdim(a2)
  %7 = 2

The library syntax is GEN alg_centralproj(GEN al, GEN z, long maps).


algchar(al)

Given an algebra al output by alginit or algtableinit, returns the characteristic of al.

  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,13);
  ? algchar(A)
  %3 = 13

The library syntax is GEN algchar(GEN al).


algcharpoly(al, b, {v = 'x}, {abs = 0})

Given an element b in al, returns its characteristic polynomial as a polynomial in the variable v. If al is a table algebra output by algtableinit or if abs = 1, returns the absolute characteristic polynomial of b, which is an element of 𝔽p[v] or ℚ[v]; if al is a central simple algebra output by alginit and abs = 0, returns the reduced characteristic polynomial of b, which is an element of K[v] where K is the center of al.

  ? al = alginit(nfinit(y), [-1,-1]); \\ (-1,-1)Q
  ? algcharpoly(al, [0,1]~)
  %2 = x^2 + 1
  ? algcharpoly(al, [0,1]~,,1)
  %3 = x^4 + 2*x^2 + 1
  ? nf = nfinit(y^2-5);
  ? al = alginit(nf,[-1,y]);
  ? a = [y,1+x]~*Mod(1,y^2-5)*Mod(1,x^2+1);
  ? P = lift(algcharpoly(al,a))
  %7 = x^2 - 2*y*x + (-2*y + 5)
  ? algcharpoly(al,a,,1)
  %8 = x^8 - 20*x^6 - 80*x^5 + 110*x^4 + 800*x^3 + 1500*x^2 - 400*x + 25
  ? lift(P*subst(P,y,-y)*Mod(1,y^2-5))^2
  %9 = x^8 - 20*x^6 - 80*x^5 + 110*x^4 + 800*x^3 + 1500*x^2 - 400*x + 25

Also accepts a square matrix with coefficients in al.

The library syntax is GEN algcharpoly(GEN al, GEN b, long v = -1, long abs) where v is a variable number.


algdegree(al)

Given a central simple algebra al output by alginit, returns the degree of al.

  ? nf = nfinit(y^3-y+1);
  ? A = alginit(nf, [-1,-1]);
  ? algdegree(A)
  %3 = 2

The library syntax is long algdegree(GEN al).


algdim(al, {abs = 0})

If al is a table algebra output by algtableinit or if abs = 1, returns the dimension of al over its prime subfield (ℚ or 𝔽p). If al is a central simple algebra output by alginit and abs = 0, returns the dimension of al over its center.

  ? nf = nfinit(y^3-y+1);
  ? A = alginit(nf, [-1,-1]);
  ? algdim(A)
  %3 = 4
  ? algdim(A,1)
  %4 = 12

The library syntax is long algdim(GEN al, long abs).


algdisc(al)

Given a central simple algebra al output by alginit, computes the discriminant of the order 𝒪0 stored in al, that is the determinant of the trace form \rm{Tr} : 𝒪0 x 𝒪0 → ℤ.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-3,1-y]);
  ? [PR,h] = alghassef(A)
  %3 = [[[2, [2, 0]~, 1, 2, 1], [3, [3, 0]~, 1, 2, 1]], Vecsmall([0, 1])]
  ? n = algdegree(A);
  ? D = algdim(A,1);
  ? h = vector(#h, i, n - gcd(n,h[i]));
  ? n^D * nf.disc^(n^2) * idealnorm(nf, idealfactorback(nf,PR,h))^n
  %4 = 12960000
  ? algdisc(A)
  %5 = 12960000

The library syntax is GEN algdisc(GEN al).


algdivl(al, x, y)

Given two elements x and y in al, computes their left quotient x\y in the algebra al: an element z such that xz = y (such an element is not unique when x is a zerodivisor). If x is invertible, this is the same as x-1y. Assumes that y is left divisible by x (i.e. that z exists). Also accepts matrices with coefficients in al.

The library syntax is GEN algdivl(GEN al, GEN x, GEN y).


algdivr(al, x, y)

Given two elements x and y in al, returns xy-1. Also accepts matrices with coefficients in al.

The library syntax is GEN algdivr(GEN al, GEN x, GEN y).


alggroup(gal, {p = 0})

Initializes the group algebra K[G] over K = ℚ (p omitted) or 𝔽p where G is the underlying group of the galoisinit structure gal. The input gal is also allowed to be a t_VEC of permutations that is closed under products.

Example:

  ? K = nfsplitting(x^3-x+1);
  ? gal = galoisinit(K);
  ? al = alggroup(gal);
  ? algissemisimple(al)
  %4 = 1
  ? G = [Vecsmall([1,2,3]), Vecsmall([1,3,2])];
  ? al2 = alggroup(G, 2);
  ? algissemisimple(al2)
  %8 = 0

The library syntax is GEN alggroup(GEN gal, GEN p = NULL).


alggroupcenter(gal, {p = 0}, {&cc})

Initializes the center Z(K[G]) of the group algebra K[G] over K = ℚ (p = 0 or omitted) or 𝔽p where G is the underlying group of the galoisinit structure gal. The input gal is also allowed to be a t_VEC of permutations that is closed under products. Sets cc to a t_VEC [elts,conjclass,rep,flag] where elts is a sorted t_VEC containing the list of elements of G, conjclass is a t_VECSMALL of the same length as elts containing the index of the conjugacy class of the corresponding element (an integer between 1 and the number of conjugacy classes), and rep is a t_VECSMALL of length the number of conjugacy classes giving for each conjugacy class the index in elts of a representative of this conjugacy class. Finally flag is 1 if and only if the permutation representation of G is transitive, in which case the i-th element of elts is characterized by g[1] = i; this is always the case when gal is a galoisinit structure. The basis of Z(K[G]) as output consists of the indicator functions of the conjugacy classes in the ordering given by cc. Example:

  ? K = nfsplitting(x^4+x+1);
  ? gal = galoisinit(K); \\ S4
  ? al = alggroupcenter(gal,,&cc);
  ? algiscommutative(al)
  %4 = 1
  ? #cc[3] \\ number of conjugacy classes of S4
  %5 = 5
  ? gal = [Vecsmall([1,2,3]),Vecsmall([1,3,2])]; \\ C2
  ? al = alggroupcenter(gal,,&cc);
  ? cc[3]
  %8 = Vecsmall([1, 2])
  ? cc[4]
  %9 = 0

The library syntax is GEN alggroupcenter(GEN gal, GEN p = NULL, GEN *cc = NULL).


alghasse(al, pl)

Given a central simple algebra al output by alginit and a prime ideal or an integer between 1 and r1+r2, returns a t_FRAC h : the local Hasse invariant of al at the place specified by pl.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,y]);
  ? alghasse(A, 1)
  %3 = 1/2
  ? alghasse(A, 2)
  %4 = 0
  ? alghasse(A, idealprimedec(nf,2)[1])
  %5 = 1/2
  ? alghasse(A, idealprimedec(nf,5)[1])
  %6 = 0

The library syntax is GEN alghasse(GEN al, GEN pl).


alghassef(al)

Given a central simple algebra al output by alginit, returns a t_VEC [PR, hf] describing the local Hasse invariants at the finite places of the center: PR is a t_VEC of primes and hf is a t_VECSMALL of integers modulo the degree d of al. The Hasse invariant of al at the primes outside PR is 0, but PR can include primes at which the Hasse invariant is 0.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,2*y-1]);
  ? [PR,hf] = alghassef(A);
  ? PR
  %4 = [[19, [10, 2]~, 1, 1, [-8, 2; 2, -10]], [2, [2, 0]~, 1, 2, 1]]
  ? hf
  %5 = Vecsmall([1, 0])

The library syntax is GEN alghassef(GEN al).


alghassei(al)

Given a central simple algebra al output by alginit, returns a t_VECSMALL hi of r1 integers modulo the degree d of al, where r1 is the number of real places of the center: the local Hasse invariants of al at infinite places.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,y]);
  ? alghassei(A)
  %3 = Vecsmall([1, 0])

The library syntax is GEN alghassei(GEN al).


algindex(al, {pl})

Returns the index of the central simple algebra A over K (as output by alginit), that is the degree e of the unique central division algebra D over K such that A is isomorphic to some matrix algebra Mk(D). If pl is set, it should be a prime ideal of K or an integer between 1 and r1+r2, and in that case return the local index at the place pl instead.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,y]);
  ? algindex(A, 1)
  %3 = 2
  ? algindex(A, 2)
  %4 = 1
  ? algindex(A, idealprimedec(nf,2)[1])
  %5 = 2
  ? algindex(A, idealprimedec(nf,5)[1])
  %6 = 1
  ? algindex(A)
  %7 = 2

The library syntax is long algindex(GEN al, GEN pl = NULL).


alginit(B, C, {v}, {maxord = 1})

Initializes the central simple algebra defined by data B, C and variable v, as follows.

* (multiplication table) B is the base number field K in nfinit form, C is a "multiplication table" over K. As a K-vector space, the algebra is generated by a basis (e1 = 1,..., en); the table is given as a t_VEC of n matrices in Mn(K), giving the left multiplication by the basis elements ei, in the given basis. Assumes that e1 = 1, that the multiplication table is integral, and that ( ⨁ i = 1^nK ei,C) describes a central simple algebra over K.

  { mi = [0,-1,0, 0;
           1, 0,0, 0;
           0, 0,0,-1;
           0, 0,1, 0];
    mj = [0, 0,-1,0;
           0, 0, 0,1;
           1, 0, 0,0;
           0,-1, 0,0];
    mk = [0, 0, 0, -1;
           0, 0,-1, 0;
           0, 1, 0, 0;
           1, 0, 0, 0];
    A = alginit(nfinit(y), [matid(4), mi,mj,mk], , 0); }

represents (in a complicated way) the quaternion algebra (-1,-1)_ℚ. See below for a simpler solution.

* (cyclic algebra) B is an rnf structure attached to a cyclic number field extension L/K of degree d, C is a t_VEC [sigma,b] with 2 components: sigma is a t_POLMOD representing an automorphism generating Gal(L/K), b is an element in K*. This represents the cyclic algebra (L/K,σ,b). Currently the element b has to be integral.

   ? Q = nfinit(y); T = polcyclo(5, 'x); F = rnfinit(Q, T);
   ? A = alginit(F, [Mod(x^2,T), 3]);

defines the cyclic algebra (L/ℚ, σ, 3), where L = ℚ(ζ5) and σ:ζζ^2 generates Gal(L/ℚ).

* (quaternion algebra, special case of the above) B is an nf structure attached to a number field K, C = [a,b] is a vector containing two elements of K* with a not a square in K, returns the quaternion algebra (a,b)K. The variable v ('x by default) must have higher priority than the variable of K.pol and is used to represent elements in the splitting field L = K[x]/(x^2-a).

   ? Q = nfinit(y); A = alginit(Q, [-1,-1]);  \\  (-1,-1)_ℚ

* (algebra/K defined by local Hasse invariants) B is an nf structure attached to a number field K, C = [d, [PR,hf], hi] is a triple containing an integer d > 1, a pair [PR, hf] describing the Hasse invariants at finite places, and hi the Hasse invariants at archimedean (real) places. A local Hasse invariant belongs to (1/d)ℤ/ℤ ⊂ ℚ/ℤ, and is given either as a t_FRAC (lift to (1/d)ℤ), a t_INT or t_INTMOD modulo d (lift to ℤ/dℤ); a whole vector of local invariants can also be given as a t_VECSMALL, whose entries are handled as t_INTs. PR is a list of prime ideals (prid structures), and hf is a vector of the same length giving the local invariants at those maximal ideals. The invariants at infinite real places are indexed by the real roots K.roots: if the Archimedean place v is attached to the j-th root, the value of hv is given by hi[j], must be 0 or 1/2 (or d/2 modulo d), and can be nonzero only if d is even.

By class field theory, provided the local invariants hv sum to 0, up to Brauer equivalence, there is a unique central simple algebra over K with given local invariants and trivial invariant elsewhere. In particular, up to isomorphism, there is a unique such algebra A of degree d.

We realize A as a cyclic algebra through class field theory. The variable v ('x by default) must have higher priority than the variable of K.pol and is used to represent elements in the (cyclic) splitting field extension L/K for A.

   ? nf = nfinit(y^2+1);
   ? PR = idealprimedec(nf,5); #PR
   %2 = 2
   ? hi = [];
   ? hf = [PR, [1/3,-1/3]];
   ? A = alginit(nf, [3,hf,hi]);
   ? algsplittingfield(A).pol
   %6 = x^3 - 21*x + 7

* (matrix algebra, toy example) B is an nf structure attached to a number field K, C = d is a positive integer. Returns a cyclic algebra isomorphic to the matrix algebra Md(K).

In all cases, this function computes a maximal order for the algebra by default, which may require a lot of time. Setting maxord = 0 prevents this computation.

The pari object representing such an algebra A is a t_VEC with the following data:

* A splitting field L of A of the same degree over K as A, in rnfinit format, accessed with algsplittingfield.

* The Hasse invariants at the real places of K, accessed with alghassei.

* The Hasse invariants of A at the finite primes of K that ramify in the natural order of A, accessed with alghassef.

* A basis of an order 𝒪0 expressed on the basis of the natural order, accessed with algbasis.

* A basis of the natural order expressed on the basis of 𝒪0, accessed with alginvbasis.

* The left multiplication table of 𝒪0 on the previous basis, accessed with algmultable.

* The characteristic of A (always 0), accessed with algchar.

* The absolute traces of the elements of the basis of 𝒪0.

* If A was constructed as a cyclic algebra (L/K,σ,b) of degree d, a t_VEC [σ,σ^2,...,σd-1]. The function algaut returns σ.

* If A was constructed as a cyclic algebra (L/K,σ,b), the element b, accessed with algb.

* If A was constructed with its multiplication table mt over K, the t_VEC of t_MAT mt, accessed with algrelmultable.

* If A was constructed with its multiplication table mt over K, a t_VEC with three components: a t_COL representing an element of A generating the splitting field L as a maximal subfield of A, a t_MAT representing an L-basis ℬ of A expressed on the ℤ-basis of 𝒪0, and a t_MAT representing the ℤ-basis of 𝒪0 expressed on ℬ. This data is accessed with algsplittingdata.

The library syntax is GEN alginit(GEN B, GEN C, long v = -1, long maxord) where v is a variable number.


alginv(al, x)

Given an element x in al, computes its inverse x-1 in the algebra al. Assumes that x is invertible.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? alginv(A,[1,1,0,0]~)
  %2 = [1/2, 1/2, 0, 0]~

Also accepts matrices with coefficients in al.

The library syntax is GEN alginv(GEN al, GEN x).


alginvbasis(al)

Given an central simple algebra al output by alginit, returns a ℤ-basis of the natural order in al with respect to the order 𝒪0 stored in al.

  A = alginit(nfinit(y), [-1,-1]);
  ? alginvbasis(A)
  %2 =
  [1 0 0 -1]
  
  [0 1 0 -1]
  
  [0 0 1 -1]
  
  [0 0 0  2]

The library syntax is GEN alginvbasis(GEN al).


algisassociative(mt, p = 0)

Returns 1 if the multiplication table mt is suitable for algtableinit(mt,p), 0 otherwise. More precisely, mt should be a t_VEC of n matrices in Mn(K), giving the left multiplications by the basis elements e1,..., en (structure constants). We check whether the first basis element e1 is 1 and ei(e_jek) = (e_iej)ek for all i,j,k.

   ? mt = [matid(3),[0,0,0;1,0,1;0,0,0],[0,0,0;0,0,0;1,0,1]];
   ? algisassociative(mt)
   %2 = 1

May be used to check a posteriori an algebra: we also allow mt as output by algtableinit (p is ignored in this case).

The library syntax is GEN algisassociative(GEN mt, GEN p).


algiscommutative(al)

al being a table algebra output by algtableinit or a central simple algebra output by alginit, tests whether the algebra al is commutative.

  ? mt = [matid(3),[0,0,0;1,0,1;0,0,0],[0,0,0;0,0,0;1,0,1]];
  ? A = algtableinit(mt);
  ? algiscommutative(A)
  %3 = 0
  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,2);
  ? algiscommutative(A)
  %6 = 1

The library syntax is GEN algiscommutative(GEN al).


algisdivision(al, {pl})

Given a central simple algebra al output by alginit, tests whether al is a division algebra. If pl is set, it should be a prime ideal of K or an integer between 1 and r1+r2, and in that case tests whether al is locally a division algebra at the place pl instead.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,y]);
  ? algisdivision(A, 1)
  %3 = 1
  ? algisdivision(A, 2)
  %4 = 0
  ? algisdivision(A, idealprimedec(nf,2)[1])
  %5 = 1
  ? algisdivision(A, idealprimedec(nf,5)[1])
  %6 = 0
  ? algisdivision(A)
  %7 = 1

The library syntax is GEN algisdivision(GEN al, GEN pl = NULL).


algisdivl(al, x, y, {&z})

Given two elements x and y in al, tests whether y is left divisible by x, that is whether there exists z in al such that xz = y, and sets z to this element if it exists.

  ? A = alginit(nfinit(y), [-1,1]);
  ? algisdivl(A,[x+2,-x-2]~,[x,1]~)
  %2 = 0
  ? algisdivl(A,[x+2,-x-2]~,[-x,x]~,&z)
  %3 = 1
  ? z
  %4 = [Mod(-2/5*x - 1/5, x^2 + 1), 0]~

Also accepts matrices with coefficients in al.

The library syntax is GEN algisdivl(GEN al, GEN x, GEN y, GEN *z = NULL).


algisinv(al, x, {&ix})

Given an element x in al, tests whether x is invertible, and sets ix to the inverse of x.

  ? A = alginit(nfinit(y), [-1,1]);
  ? algisinv(A,[-1,1]~)
  %2 = 0
  ? algisinv(A,[1,2]~,&ix)
  %3 = 1
  ? ix
  %4 = [Mod(Mod(-1/3, y), x^2 + 1), Mod(Mod(2/3, y), x^2 + 1)]~

Also accepts matrices with coefficients in al.

The library syntax is GEN algisinv(GEN al, GEN x, GEN *ix = NULL).


algisramified(al, {pl})

Given a central simple algebra al output by alginit, tests whether al is ramified, i.e. not isomorphic to a matrix algebra over its center. If pl is set, it should be a prime ideal of K or an integer between 1 and r1+r2, and in that case tests whether al is locally ramified at the place pl instead.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,y]);
  ? algisramified(A, 1)
  %3 = 1
  ? algisramified(A, 2)
  %4 = 0
  ? algisramified(A, idealprimedec(nf,2)[1])
  %5 = 1
  ? algisramified(A, idealprimedec(nf,5)[1])
  %6 = 0
  ? algisramified(A)
  %7 = 1

The library syntax is GEN algisramified(GEN al, GEN pl = NULL).


algissemisimple(al)

al being a table algebra output by algtableinit or a central simple algebra output by alginit, tests whether the algebra al is semisimple.

  ? mt = [matid(3),[0,0,0;1,0,1;0,0,0],[0,0,0;0,0,0;1,0,1]];
  ? A = algtableinit(mt);
  ? algissemisimple(A)
  %3 = 0
  ? mi=[0,-1,0,0;1,0,0,0;0,0,0,-1;0,0,1,0]; \\ quaternion algebra (-1,-1)
  ? mj=[0,0,-1,0;0,0,0,1;1,0,0,0;0,-1,0,0];
  ? mk=[0,0,0,-1;0,0,-1,0;0,1,0,0;1,0,0,0];
  ? mt = [matid(4), mi, mj, mk];
  ? A = algtableinit(mt);
  ? algissemisimple(A)
  %9 = 1

The library syntax is GEN algissemisimple(GEN al).


algissimple(al, {ss = 0})

al being a table algebra output by algtableinit or a central simple algebra output by alginit, tests whether the algebra al is simple. If ss = 1, assumes that the algebra al is semisimple without testing it.

  ? mt = [matid(3),[0,0,0;1,0,1;0,0,0],[0,0,0;0,0,0;1,0,1]];
  ? A = algtableinit(mt); \\ matrices [*,*; 0,*]
  ? algissimple(A)
  %3 = 0
  ? algissimple(A,1) \\ incorrectly assume that A is semisimple
  %4 = 1
  ? mi=[0,-1,0,0;1,0,0,0;0,0,0,-1;0,0,1,0];
  ? mj=[0,0,-1,0;0,0,0,1;1,0,0,0;0,-1,0,0];
  ? mk=[0,0,0,-1;0,0,b,0;0,1,0,0;1,0,0,0];
  ? mt = [matid(4), mi, mj, mk];
  ? A = algtableinit(mt); \\ quaternion algebra (-1,-1)
  ? algissimple(A)
  %10 = 1
  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,2); \\ direct product F4 x F2
  ? algissimple(A)
  %13 = 0

The library syntax is GEN algissimple(GEN al, long ss).


algissplit(al, {pl})

Given a central simple algebra al output by alginit, tests whether al is split, i.e. isomorphic to a matrix algebra over its center. If pl is set, it should be a prime ideal of K or an integer between 1 and r1+r2, and in that case tests whether al is locally split at the place pl instead.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,y]);
  ? algissplit(A, 1)
  %3 = 0
  ? algissplit(A, 2)
  %4 = 1
  ? algissplit(A, idealprimedec(nf,2)[1])
  %5 = 0
  ? algissplit(A, idealprimedec(nf,5)[1])
  %6 = 1
  ? algissplit(A)
  %7 = 0

The library syntax is GEN algissplit(GEN al, GEN pl = NULL).


alglatadd(al, lat1, lat2, {&ptinter})

Given an algebra al and two lattices lat1 and lat2 in al, computes the sum lat1 + lat2. If ptinter is present, set it to the intersection lat1 ∩ lat2.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? lat1 = alglathnf(al,[1,1,0,0,0,0,0,0]~);
  ? lat2 = alglathnf(al,[1,0,1,0,0,0,0,0]~);
  ? latsum = alglatadd(al,lat1,lat2,&latinter);
  ? matdet(latsum[1])
  %5 = 4
  ? matdet(latinter[1])
  %6 = 64

The library syntax is GEN alglatadd(GEN al, GEN lat1, GEN lat2, GEN *ptinter = NULL).


alglatcontains(al, lat, x, {&ptc})

Given an algebra al, a lattice lat and x in al, tests whether x ∈ lat. If ptc is present, sets it to the t_COL of coordinates of x in the basis of lat.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? a1 = [1,-1,0,1,2,0,1,2]~;
  ? lat1 = alglathnf(al,a1);
  ? alglatcontains(al,lat1,a1,&c)
  %4 = 1
  ? c
  %5 = [-1, -2, -1, 1, 2, 0, 1, 1]~

The library syntax is GEN alglatcontains(GEN al, GEN lat, GEN x, GEN *ptc = NULL).


alglatelement(al, lat, c)

Given an algebra al, a lattice lat and a t_COL c, returns the element of al whose coordinates on the ℤ-basis of lat are given by c.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? a1 = [1,-1,0,1,2,0,1,2]~;
  ? lat1 = alglathnf(al,a1);
  ? c = [1..8]~;
  ? elt = alglatelement(al,lat1,c);
  ? alglatcontains(al,lat1,elt,&c2)
  %6 = 1
  ? c==c2
  %7 = 1

The library syntax is GEN alglatelement(GEN al, GEN lat, GEN c).


alglathnf(al, m, {d = 0})

Given an algebra al and a matrix m with columns representing elements of al, returns the lattice L generated by the columns of m. If provided, d must be a rational number such that L contains d times the natural basis of al. The argument m is also allowed to be a t_VEC of t_MAT, in which case m is replaced by the concatenation of the matrices, or a t_COL, in which case m is replaced by its left multiplication table as an element of al.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? a = [1,1,-1/2,1,1/3,-1,1,1]~;
  ? mt = algtomatrix(al,a,1);
  ? lat = alglathnf(al,mt);
  ? lat[2]
  %5 = 1/6

The library syntax is GEN alglathnf(GEN al, GEN m, GEN d).


alglatindex(al, lat1, lat2)

Given an algebra al and two lattices lat1 and lat2 in al, computes the generalized index of lat1 relative to lat2, i.e. |lat2/lat1∩ lat2|/|lat1/lat1∩ lat2|.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? lat1 = alglathnf(al,[1,1,0,0,0,0,0,0]~);
  ? lat2 = alglathnf(al,[1,0,1,0,0,0,0,0]~);
  ? alglatindex(al,lat1,lat2)
  %4 = 1
  ? lat1==lat2
  %5 = 0

The library syntax is GEN alglatindex(GEN al, GEN lat1, GEN lat2).


alglatinter(al, lat1, lat2, {&ptsum})

Given an algebra al and two lattices lat1 and lat2 in al, computes the intersection lat1∩ lat2. If ptsum is present, sets it to the sum lat1 + lat2.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? lat1 = alglathnf(al,[1,1,0,0,0,0,0,0]~);
  ? lat2 = alglathnf(al,[1,0,1,0,0,0,0,0]~);
  ? latinter = alglatinter(al,lat1,lat2,&latsum);
  ? matdet(latsum[1])
  %5 = 4
  ? matdet(latinter[1])
  %6 = 64

The library syntax is GEN alglatinter(GEN al, GEN lat1, GEN lat2, GEN *ptsum = NULL).


alglatlefttransporter(al, lat1, lat2)

Given an algebra al and two lattices lat1 and lat2 in al, computes the left transporter from lat1 to lat2, i.e. the set of x ∈ al such that x.lat1 ⊂ lat2.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? lat1 = alglathnf(al,[1,-1,0,1,2,0,5,2]~);
  ? lat2 = alglathnf(al,[0,1,-2,-1,0,0,3,1]~);
  ? tr = alglatlefttransporter(al,lat1,lat2);
  ? a = alglatelement(al,tr,[0,0,0,0,0,0,1,0]~);
  ? alglatsubset(al,alglatmul(al,a,lat1),lat2)
  %6 = 1
  ? alglatsubset(al,alglatmul(al,lat1,a),lat2)
  %7 = 0

The library syntax is GEN alglatlefttransporter(GEN al, GEN lat1, GEN lat2).


alglatmul(al, lat1, lat2)

Given an algebra al and two lattices lat1 and lat2 in al, computes the lattice generated by the products of elements of lat1 and lat2. One of lat1 and lat2 is also allowed to be an element of al; in this case, computes the product of the element and the lattice.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? a1 = [1,-1,0,1,2,0,1,2]~;
  ? a2 = [0,1,2,-1,0,0,3,1]~;
  ? lat1 = alglathnf(al,a1);
  ? lat2 = alglathnf(al,a2);
  ? lat3 = alglatmul(al,lat1,lat2);
  ? matdet(lat3[1])
  %7 = 29584
  ? lat3 == alglathnf(al, algmul(al,a1,a2))
  %8 = 0
  ? lat3 == alglatmul(al, lat1, a2)
  %9 = 0
  ? lat3 == alglatmul(al, a1, lat2)
  %10 = 0

The library syntax is GEN alglatmul(GEN al, GEN lat1, GEN lat2).


alglatrighttransporter(al, lat1, lat2)

Given an algebra al and two lattices lat1 and lat2 in al, computes the right transporter from lat1 to lat2, i.e. the set of x ∈ al such that lat1.x ⊂ lat2.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? lat1 = alglathnf(al,matdiagonal([1,3,7,1,2,8,5,2]));
  ? lat2 = alglathnf(al,matdiagonal([5,3,8,1,9,8,7,1]));
  ? tr = alglatrighttransporter(al,lat1,lat2);
  ? a = alglatelement(al,tr,[0,0,0,0,0,0,0,1]~);
  ? alglatsubset(al,alglatmul(al,lat1,a),lat2)
  %6 = 1
  ? alglatsubset(al,alglatmul(al,a,lat1),lat2)
  %7 = 0

The library syntax is GEN alglatrighttransporter(GEN al, GEN lat1, GEN lat2).


alglatsubset(al, lat1, lat2, {&ptindex})

Given an algebra al and two lattices lat1 and lat2 in al, tests whether lat1 ⊂ lat2. If it is true and ptindex is present, sets it to the index of lat1 in lat2.

  ? al = alginit(nfinit(y^2+7), [-1,-1]);
  ? lat1 = alglathnf(al,[1,1,0,0,0,0,0,0]~);
  ? lat2 = alglathnf(al,[1,0,1,0,0,0,0,0]~);
  ? alglatsubset(al,lat1,lat2)
  %4 = 0
  ? latsum = alglatadd(al,lat1,lat2);
  ? alglatsubset(al,lat1,latsum,&index)
  %6 = 1
  ? index
  %7 = 4

The library syntax is GEN alglatsubset(GEN al, GEN lat1, GEN lat2, GEN *ptindex = NULL).


algmakeintegral(mt, {maps = 0})

mt being a multiplication table over ℚ in the same format as the input of algtableinit, computes an integral multiplication table mt2 for an isomorphic algebra. When maps = 1, returns a t_VEC [mt2,S,T] where S and T are matrices respectively representing the map from the algebra defined by mt to the one defined by mt2 and its inverse.

  ? mt = [matid(2),[0,-1/4;1,0]];
  ? algtableinit(mt);
    ***   at top-level: algtableinit(mt)
    ***                 ^ —  —  —  —  — -
    *** algtableinit: domain error in algtableinit: denominator(mt) != 1
  ? mt2 = algmakeintegral(mt);
  ? al = algtableinit(mt2);
  ? algisassociative(al)
  %4 = 1
  ? [mt2, S, T] = algmakeintegral(mt,1);
  ? S
  %6 =
  [1   0]
  
  [0 1/4]
  ? T
  %7 =
  [1 0]
  
  [0 4]
  ? vector(#mt, i, S * (mt * T[,i]) * T) == mt2
  %8 = 1

The library syntax is GEN algmakeintegral(GEN mt, long maps).


algmul(al, x, y)

Given two elements x and y in al, computes their product xy in the algebra al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algmul(A,[1,1,0,0]~,[0,0,2,1]~)
  %2 = [2, 3, 5, -4]~

Also accepts matrices with coefficients in al.

The library syntax is GEN algmul(GEN al, GEN x, GEN y).


algmultable(al)

Returns a multiplication table of al over its prime subfield (ℚ or 𝔽p), as a t_VEC of t_MAT: the left multiplication tables of basis elements. If al was output by algtableinit, returns the multiplication table used to define al. If al was output by alginit, returns the multiplication table of the order 𝒪0 stored in al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? M = algmultable(A);
  ? #M
  %3 = 4
  ? M[1]  \\ multiplication by e1 = 1
  %4 =
  [1 0 0 0]
  
  [0 1 0 0]
  
  [0 0 1 0]
  
  [0 0 0 1]
  
  ? M[2]
  %5 =
  [0 -1  1  0]
  
  [1  0  1  1]
  
  [0  0  1  1]
  
  [0  0 -2 -1]

The library syntax is GEN algmultable(GEN al).


algneg(al, x)

Given an element x in al, computes its opposite -x in the algebra al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algneg(A,[1,1,0,0]~)
  %2 = [-1, -1, 0, 0]~

Also accepts matrices with coefficients in al.

The library syntax is GEN algneg(GEN al, GEN x).


algnorm(al, x, {abs = 0})

Given an element x in al, computes its norm. If al is a table algebra output by algtableinit or if abs = 1, returns the absolute norm of x, which is an element of 𝔽p of ℚ; if al is a central simple algebra output by alginit and abs = 0 (default), returns the reduced norm of x, which is an element of the center of al.

  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,19);
  ? algnorm(A,[0,-2,3]~)
  %3 = 18
  ? nf = nfinit(y^2-5);
  ? B = alginit(nf,[-1,y]);
  ? b = [x,1]~;
  ? n = algnorm(B,b)
  %7 = Mod(-y + 1, y^2 - 5)
  ? algnorm(B,b,1)
  %8 = 16
  ? nfeltnorm(nf,n)^algdegree(B)
  %9 = 16

Also accepts a square matrix with coefficients in al.

The library syntax is GEN algnorm(GEN al, GEN x, long abs).


algpoleval(al, T, b)

Given an element b in al and a polynomial T in K[X], computes T(b) in al. Also accepts as input a t_VEC [b,mb] where mb is the left multiplication table of b.

  ? nf = nfinit(y^2-5);
  ? al = alginit(nf,[y,-1]);
  ? b = [1..8]~;
  ? pol = algcharpoly(al,b,,1);
  ? algpoleval(al,pol,b)==0
  %5 = 1
  ? mb = algtomatrix(al,b,1);
  ? algpoleval(al,pol,[b,mb])==0
  %7 = 1

The library syntax is GEN algpoleval(GEN al, GEN T, GEN b).


algpow(al, x, n)

Given an element x in al and an integer n, computes the power x^n in the algebra al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algpow(A,[1,1,0,0]~,7)
  %2 = [8, -8, 0, 0]~

Also accepts a square matrix with coefficients in al.

The library syntax is GEN algpow(GEN al, GEN x, GEN n).


algprimesubalg(al)

al being the output of algtableinit representing a semisimple algebra of positive characteristic, returns a basis of the prime subalgebra of al. The prime subalgebra of al is the subalgebra fixed by the Frobenius automorphism of the center of al. It is abstractly isomorphic to a product of copies of 𝔽p.

  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,2);
  ? algprimesubalg(A)
  %3 =
  [1 0]
  
  [0 1]
  
  [0 0]

The library syntax is GEN algprimesubalg(GEN al).


algquotient(al, I, {maps = 0})

al being a table algebra output by algtableinit and I being a basis of a two-sided ideal of al represented by a matrix, returns the quotient al/I. When maps = 1, returns a t_VEC [al/I,proj,lift] where proj and lift are matrices respectively representing the projection map and a section of it.

  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,2);
  ? AQ = algquotient(A,[0;1;0]);
  ? algdim(AQ)
  %4 = 2

The library syntax is GEN alg_quotient(GEN al, GEN I, long maps).


algradical(al)

al being a table algebra output by algtableinit, returns a basis of the Jacobson radical of the algebra al over its prime field (ℚ or 𝔽p).

Here is an example with A = ℚ[x]/(x^2), with the basis (1,x):

  ? mt = [matid(2),[0,0;1,0]];
  ? A = algtableinit(mt);
  ? algradical(A) \\ = (x)
  %3 =
  [0]
  
  [1]

Another one with 2 x 2 upper triangular matrices over ℚ, with basis I2, a = [0,1;0,0] and b = [0,0;0,1], such that a^2 = 0, ab = a, ba = 0, b^2 = b:

  ? mt = [matid(3),[0,0,0;1,0,1;0,0,0],[0,0,0;0,0,0;1,0,1]];
  ? A = algtableinit(mt);
  ? algradical(A) \\ = (a)
  %6 =
  [0]
  
  [1]
  
  [0]

The library syntax is GEN algradical(GEN al).


algramifiedplaces(al)

Given a central simple algebra al output by alginit, returns a t_VEC containing the list of places of the center of al that are ramified in al. Each place is described as an integer between 1 and r1 or as a prime ideal.

  ? nf = nfinit(y^2-5);
  ? A = alginit(nf, [-1,y]);
  ? algramifiedplaces(A)
  %3 = [1, [2, [2, 0]~, 1, 2, 1]]

The library syntax is GEN algramifiedplaces(GEN al).


algrandom(al, b)

Given an algebra al and an integer b, returns a random element in al with coefficients in [-b,b].

The library syntax is GEN algrandom(GEN al, GEN b).


algrelmultable(al)

Given a central simple algebra al output by alginit defined by a multiplication table over its center (a number field), returns this multiplication table.

  ? nf = nfinit(y^3-5); a = y; b = y^2;
  ? {mi = [0,a,0,0;
            1,0,0,0;
            0,0,0,a;
            0,0,1,0];}
  ? {mj = [0, 0,b, 0;
            0, 0,0,-b;
            1, 0,0, 0;
            0,-1,0, 0];}
  ? {mk = [0, 0,0,-a*b;
            0, 0,b,   0;
            0,-a,0,   0;
            1, 0,0,   0];}
  ? mt = [matid(4), mi, mj, mk];
  ? A = alginit(nf,mt,'x);
  ? M = algrelmultable(A);
  ? M[2] == mi
  %8 = 1
  ? M[3] == mj
  %9 = 1
  ? M[4] == mk
  %10 = 1

The library syntax is GEN algrelmultable(GEN al).


algsimpledec(al, {maps = 0})

al being the output of algtableinit, returns a t_VEC [J,[al1,al2,...,aln]] where J is a basis of the Jacobson radical of al and al/J is isomorphic to the direct product of the simple algebras ali. When maps = 1, each ali is replaced with a t_VEC [ali,proji,lifti] where proji and lifti are matrices respectively representing the projection map alali and a section of it. Modulo J, the images of the lifti form a direct sum in al/J, so that the images of 1 ∈ ali under lifti are central primitive idempotents of al/J. The factors are sorted by increasing dimension, then increasing dimension of the center. This ensures that the ordering of the isomorphism classes of the factors is deterministic over finite fields, but not necessarily over ℚ.

The library syntax is GEN algsimpledec(GEN al, long maps).


algsplit(al, {v = 'x})

If al is a table algebra over 𝔽p output by algtableinit that represents a simple algebra, computes an isomorphism between al and a matrix algebra Md(𝔽p^n) where N = nd^2 is the dimension of al. Returns a t_VEC [map,mapi], where:

* map is a t_VEC of N matrices of size d x d with t_FFELT coefficients using the variable v, representing the image of the basis of al under the isomorphism.

* mapi is an N x N matrix with t_INT coefficients, representing the image in al by the inverse isomorphism of the basis (bi) of Md(𝔽p[α]) (where α has degree n over 𝔽p) defined as follows: let Ei,j be the matrix having all coefficients 0 except the (i,j)-th coefficient equal to 1, and define bi3+n(i2+di1)+1 = Ei1+1,i2+1 αi3, where 0 ≤ i1,i2 < d and 0 ≤ i3 < n.

Example:

  ? al0 = alginit(nfinit(y^2+7), [-1,-1]);
  ? al = algtableinit(algmultable(al0), 3); \\ isomorphic to M2(F9)
  ? [map,mapi] = algsplit(al, 'a);
  ? x = [1,2,1,0,0,0,0,0]~; fx = map*x
  %4 =
  [2*a 0]
  
  [  0 2]
  ? y = [0,0,0,0,1,0,0,1]~; fy = map*y
  %5 =
  [1   2*a]
  
  [2 a + 2]
  ? map*algmul(al,x,y) == fx*fy
  %6 = 1
  ? map*mapi[,6]
  %7 =
  [0 0]
  
  [a 0]

Warning. If al is not simple, algsplit(al) can trigger an error, but can also run into an infinite loop. Example:

  ? al = alginit(nfinit(y),[-1,-1]); \\ ramified at 2
  ? al2 = algtableinit(algmultable(al),2); \\ maximal order modulo 2
  ? algsplit(al2); \\ not semisimple, infinite loop

The library syntax is GEN algsplit(GEN al, long v = -1) where v is a variable number.


algsplittingdata(al)

Given a central simple algebra al output by alginit defined by a multiplication table over its center K (a number field), returns data stored to compute a splitting of al over an extension. This data is a t_VEC [t,Lbas,Lbasinv] with 3 components:

* an element t of al such that L = K(t) is a maximal subfield of al;

* a matrix Lbas expressing a L-basis of al (given an L-vector space structure by multiplication on the right) on the integral basis of al;

* a matrix Lbasinv expressing the integral basis of al on the previous L-basis.

  ? nf = nfinit(y^3-5); a = y; b = y^2;
  ? {mi = [0,a,0,0;
            1,0,0,0;
            0,0,0,a;
            0,0,1,0];}
  ? {mj = [0, 0,b, 0;
            0, 0,0,-b;
            1, 0,0, 0;
            0,-1,0, 0];}
  ? {mk = [0, 0,0,-a*b;
            0, 0,b,   0;
            0,-a,0,   0;
            1, 0,0,   0];}
  ? mt = [matid(4), mi, mj, mk];
  ? A = alginit(nf,mt,'x);
  ? [t,Lb,Lbi] = algsplittingdata(A);
  ? t
  %8 = [0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0]~;
  ? matsize(Lb)
  %9 = [12, 2]
  ? matsize(Lbi)
  %10 = [2, 12]

The library syntax is GEN algsplittingdata(GEN al).


algsplittingfield(al)

Given a central simple algebra al output by alginit, returns an rnf structure: the splitting field of al that is stored in al, as a relative extension of the center.

  nf = nfinit(y^3-5);
  a = y; b = y^2;
  {mi = [0,a,0,0;
         1,0,0,0;
         0,0,0,a;
         0,0,1,0];}
  {mj = [0, 0,b, 0;
         0, 0,0,-b;
         1, 0,0, 0;
         0,-1,0, 0];}
  {mk = [0, 0,0,-a*b;
         0, 0,b,   0;
         0,-a,0,   0;
         1, 0,0,   0];}
  mt = [matid(4), mi, mj, mk];
  A = alginit(nf,mt,'x);
  algsplittingfield(A).pol
  %8 = x^2 - y

The library syntax is GEN algsplittingfield(GEN al).


algsqr(al, x)

Given an element x in al, computes its square x^2 in the algebra al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algsqr(A,[1,0,2,0]~)
  %2 = [-3, 0, 4, 0]~

Also accepts a square matrix with coefficients in al.

The library syntax is GEN algsqr(GEN al, GEN x).


algsub(al, x, y)

Given two elements x and y in al, computes their difference x-y in the algebra al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algsub(A,[1,1,0,0]~,[1,0,1,0]~)
  %2 = [0, 1, -1, 0]~

Also accepts matrices with coefficients in al.

The library syntax is GEN algsub(GEN al, GEN x, GEN y).


algsubalg(al, B)

al being a table algebra output by algtableinit and B being a basis of a subalgebra of al represented by a matrix, computes an algebra al2 isomorphic to B.

Returns [al2,B2] where B2 is a possibly different basis of the subalgebra al2, with respect to which the multiplication table of al2 is defined.

  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,2);
  ? B = algsubalg(A,[1,0; 0,0; 0,1]);
  ? algdim(A)
  %4 = 3
  ? algdim(B[1])
  %5 = 2
  ? m = matcompanion(x^4+1);
  ? mt = [m^i | i <- [0..3]];
  ? al = algtableinit(mt);
  ? B = [1,0;0,0;0,1/2;0,0];
  ? al2 = algsubalg(al,B);
  ? algdim(al2[1])
  ? al2[2]
  %13 =
  [1 0]
  
  [0 0]
  
  [0 1]
  
  [0 0]

The library syntax is GEN algsubalg(GEN al, GEN B).


algtableinit(mt, {p = 0})

Initializes the associative algebra over K = ℚ (p omitted) or 𝔽p defined by the multiplication table mt. As a K-vector space, the algebra is generated by a basis (e1 = 1, e2,..., en); the table is given as a t_VEC of n matrices in Mn(K), giving the left multiplication by the basis elements ei, in the given basis. Assumes that e1 = 1, that K e1⨁ ...⨁ K en] describes an associative algebra over K, and in the case K = ℚ that the multiplication table is integral. If the algebra is already known to be central and simple, then the case K = 𝔽p is useless, and one should use alginit directly.

The point of this function is to input a finite dimensional K-algebra, so as to later compute its radical, then to split the quotient algebra as a product of simple algebras over K.

The pari object representing such an algebra A is a t_VEC with the following data:

* The characteristic of A, accessed with algchar.

* The multiplication table of A, accessed with algmultable.

* The traces of the elements of the basis.

A simple example: the 2 x 2 upper triangular matrices over ℚ, generated by I2, a = [0,1;0,0] and b = [0,0;0,1], such that a^2 = 0, ab = a, ba = 0, b^2 = b:

  ? mt = [matid(3),[0,0,0;1,0,1;0,0,0],[0,0,0;0,0,0;1,0,1]];
  ? A = algtableinit(mt);
  ? algradical(A) \\ = (a)
  %6 =
  [0]
  
  [1]
  
  [0]
  ? algcenter(A) \\ = (I2)
  %7 =
  [1]
  
  [0]
  
  [0]

The library syntax is GEN algtableinit(GEN mt, GEN p = NULL).


algtensor(al1, al2, {maxord = 1})

Given two algebras al1 and al2, computes their tensor product. Computes a maximal order by default. Prevent this computation by setting maxord = 0.

Currently only implemented for cyclic algebras of coprime degree over the same center K, and the tensor product is over K.

The library syntax is GEN algtensor(GEN al1, GEN al2, long maxord).


algtomatrix(al, x, {abs = 1})

Given an element x in al, returns the image of x under a homomorphism to a matrix algebra. If al is a table algebra output by algtableinit or if abs = 1, returns the left multiplication table on the integral basis; if al is a central simple algebra and abs = 0, returns φ(x) where φ : A ⨂ K L → Md(L) (where d is the degree of the algebra and L is an extension of L with [L:K] = d) is an isomorphism stored in al. Also accepts a square matrix with coefficients in al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algtomatrix(A,[0,0,0,2]~)
  %2 =
  [Mod(x + 1, x^2 + 1) Mod(Mod(1, y)*x + Mod(-1, y), x^2 + 1)]
  
  [Mod(x + 1, x^2 + 1)                   Mod(-x + 1, x^2 + 1)]
  ? algtomatrix(A,[0,1,0,0]~,1)
  %2 =
  [0 -1  1  0]
  
  [1  0  1  1]
  
  [0  0  1  1]
  
  [0  0 -2 -1]
  ? algtomatrix(A,[0,x]~,1)
  %3 =
  [-1  0 0 -1]
  
  [-1  0 1  0]
  
  [-1 -1 0 -1]
  
  [ 2  0 0  1]

Also accepts matrices with coefficients in al.

The library syntax is GEN algtomatrix(GEN al, GEN x, long abs).


algtrace(al, x, {abs = 0})

Given an element x in al, computes its trace. If al is a table algebra output by algtableinit or if abs = 1, returns the absolute trace of x, which is an element of 𝔽p or ℚ; if al is the output of alginit and abs = 0 (default), returns the reduced trace of x, which is an element of the center of al.

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algtrace(A,[5,0,0,1]~)
  %2 = 11
  ? algtrace(A,[5,0,0,1]~,1)
  %3 = 22
  ? nf = nfinit(y^2-5);
  ? A = alginit(nf,[-1,y]);
  ? a = [1+x+y,2*y]~*Mod(1,y^2-5)*Mod(1,x^2+1);
  ? t = algtrace(A,a)
  %7 = Mod(2*y + 2, y^2 - 5)
  ? algtrace(A,a,1)
  %8 = 8
  ? algdegree(A)*nfelttrace(nf,t)
  %9 = 8

Also accepts a square matrix with coefficients in al.

The library syntax is GEN algtrace(GEN al, GEN x, long abs).


algtype(al)

Given an algebra al output by alginit or by algtableinit, returns an integer indicating the type of algebra:

* 0: not a valid algebra.

* 1: table algebra output by algtableinit.

* 2: central simple algebra output by alginit and represented by a multiplication table over its center.

* 3: central simple algebra output by alginit and represented by a cyclic algebra.

  ? algtype([])
  %1 = 0
  ? mt = [matid(3), [0,0,0; 1,1,0; 0,0,0], [0,0,1; 0,0,0; 1,0,1]];
  ? A = algtableinit(mt,2);
  ? algtype(A)
  %4 = 1
  ? nf = nfinit(y^3-5);
  ?  a = y; b = y^2;
  ?  {mi = [0,a,0,0;
             1,0,0,0;
             0,0,0,a;
             0,0,1,0];}
  ?  {mj = [0, 0,b, 0;
             0, 0,0,-b;
             1, 0,0, 0;
             0,-1,0, 0];}
  ?  {mk = [0, 0,0,-a*b;
             0, 0,b,   0;
             0,-a,0,   0;
             1, 0,0,   0];}
  ?  mt = [matid(4), mi, mj, mk];
  ?  A = alginit(nf,mt,'x);
  ? algtype(A)
  %12 = 2
  ? A = alginit(nfinit(y), [-1,-1]);
  ? algtype(A)
  %14 = 3

The library syntax is long algtype(GEN al).