Pari/GP Reference Documentation HOME Contents - Global index - GP keyboard shortcuts

Associative and central simple algebras


algabsdim   algadd   algalgtobasis   algaut   algb   algbasis   algbasistoalg   algcenter   algcentralproj   algchar   algcharpoly   algdecomposition   algdegree   algdim   algdisc   algdivl   algdivr   alggroup   alghasse   alghassef   alghassei   algindex   alginit   alginv   alginvbasis   algisassociative   algiscommutative   algisdivision   algisdivl   algisinv   algisramified   algissemisimple   algissimple   algissplit   alglathnf   algleftmultable   algmul   algmultable   algneg   algnorm   algpoleval   algpow   algprimesubalg   algquotient   algradical   algramifiedplaces   algrandom   algrelmultable   algsimpledec   algsplittingdata   algsplittingfield   algsplittingmatrix   algsqr   algsub   algsubalg   algtableinit   algtensor   algtrace   algtype  
 

This section collects functions related to associative algebras and central simple algebras over number fields. Let A be a finite-dimensional unitary associative algebra over a field K. We say that A is central if the center of A is K, and that A is simple if it has no nontrivial two-sided ideals.

We provide functions to manipulate 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 manipulate central simple algebras over number fields. We represent them either by the left multiplication table on a basis over the center, 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. An algebra is semisimple if its Jacobson radical is trivial. A semisimple algebra is isomorphic to a direct sum of simple algebras. The dimension of a central simple algebra A over K is always a square d^2, and the integer d is called the degree of the algebra A over K. A central simple algebra A over a field K is always isomorphic to M_d(D) for some integer d and some central division algebra D of degree e : the integer e is called the index of A.

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 = 0}^{d-1}x^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 homomorphism (L/K,σ,b) → M_d(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 ∈ (1/d)ℤ/ℤ ⊂ ℚ/ℤ is called the Hasse invariant of A.

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. 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 = 0}^{d-1}x^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 = 0}^{d-1}x^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 b_0,...,b_{nd^2-1} where b_i = x^{i/(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 e_0,...,e_{d^2-1} be a basis of A over K and assume that the left multiplication table of A on (e_i) is integral. Then the submodule 𝒪 = ⨁ _{i = 0}^{d^2-1}ℤ_K e_i is an order in A, called the natural order. Let ω_0,...,ω_{n-1} be a ℤ-basis of ℤ_K. The natural basis of 𝒪 is b_0,...,b_{nd^2-1} where b_i = ω_{(i mod n)}e_{i/n}.

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, 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 = 0}^{d-1}x^iL. For a central simple algebra defined by a multiplication table over its center K on a basis (e_i), the algebraic representation is a t_COL with coefficients in K representing the element on the basis (e_i).

Warning. The coefficients in the decomposition A = ⨁ _{i = 0}^{d-1}x^iL are not the same as those in the decomposition A = ⨁ _{i = 0}^{d-1}Lx^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 behaviour 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.

algabsdim HOME   TOP

Given an algebra al output by alginit or by algtableinit, returns the dimension of al over its prime subfield (ℚ or 𝔽_p).

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

The library syntax is long algabsdim(GEN al).

algadd HOME   TOP

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 HOME   TOP

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 HOME   TOP

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 HOME   TOP

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 HOME   TOP

Given an 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 HOME   TOP

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 HOME   TOP

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 I_2, 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) \\ = (I_2)
  %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 HOME   TOP

Given a table algebra al output by algtableinit and a t_VEC z = [z_1,...,z_n] of orthogonal central idempotents, returns a t_VEC [al_1,...,al_n] of algebras such that al_i = z_i al. If maps = 1, each al_i 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⨁ 𝔽_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 HOME   TOP

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 HOME   TOP

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, 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, 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

Also accepts a square matrix with coefficients in al.

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

algdecomposition HOME   TOP

al being a table algebra output by algtableinit, returns [J,[al_1,...,al_n]] where J is a basis of the Jacobson radical of al and al_1,...,al_n are the simple factors of the semisimple algebra al/J.

The library syntax is GEN alg_decomposition(GEN al).

algdegree HOME   TOP

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 HOME   TOP

Given a central simple algebra al output by alginit, returns the dimension of al over its center. Given a table algebra al output by algtableinit, returns the dimension of al over its prime subfield (ℚ or 𝔽_p).

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

The library syntax is long algdim(GEN al).

algdisc HOME   TOP

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 = algabsdim(A);
  ? 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 HOME   TOP

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^{-1}y. 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 HOME   TOP

Given two elements x and y in al, return xy^{-1}. Also accepts matrices with coefficients in al.

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

alggroup HOME   TOP

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

alghasse HOME   TOP

Given a central simple algebra al output by alginit and a prime ideal or an integer between 1 and r_1+r_2, 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 HOME   TOP

Given a central simple algebra al output by alginit, returns a t_VEC [PR, h_f] describing the local Hasse invariants at the finite places of the center: PR is a t_VEC of primes and h_f is a t_VECSMALL of integers modulo the degree d of al.

  ? 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 HOME   TOP

Given a central simple algebra al output by alginit, returns a t_VECSMALL h_i of r_1 integers modulo the degree d of al, where r_1 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 HOME   TOP

Return 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 M_d(D). If pl is set, it should be a prime ideal of K or an integer between 1 and r_1+r_2, 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 HOME   TOP

Initialize 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 (e_1 = 1,..., e_n); the table is given as a t_VEC of n matrices in M_n(K), giving the left multiplication by the basis elements e_i, in the given basis. Assumes that e_1 = 1, that the multiplication table is integral, and that K[e_1,...,e_n] describes a central simple algebra over K.

  { m_i = [0,-1,0, 0;
           1, 0,0, 0;
           0, 0,0,-1;
           0, 0,1, 0];
    m_j = [0, 0,-1,0;
           0, 0, 0,1;
           1, 0, 0,0;
           0,-1, 0,0];
    m_k = [0, 0, 0, 0;
           0, 0,-1, 0;
           0, 1, 0, 0;
           1, 0, 0,-1];
    A = alginit(nfinit(y), [matid(4), m_i,m_j,m_k],  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,h_f], h_i] is a triple containing an integer d > 1, a pair [PR, h_f] describing the Hasse invariants at finite places, and h_i 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 h_f 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 h_v is given by h_i[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 h_v 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 M_d(K).

In all cases, this function computes a maximal order for the algebra by default, which may require a lot of time. Setting flag = 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 same splitting field L in nfinit format.

* 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 algord.

* A basis of the natural order expressed on the basis of 𝒪_0, accessed with alginvord.

* 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 flag) where v is a variable number.

alginv HOME   TOP

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 HOME   TOP

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 HOME   TOP

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 M_n(K), giving the left multiplications by the basis elements e_1,..., e_n (structure constants). We check whether the first basis element e_1 is 1 and e_i(e_je_k) = (e_ie_j)e_k 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 HOME   TOP

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 HOME   TOP

Given a central simple algebra al output by alginit, test whether al is a division algebra. If pl is set, it should be a prime ideal of K or an integer between 1 and r_1+r_2, and in that case test 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 HOME   TOP

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 HOME   TOP

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 HOME   TOP

Given a central simple algebra al output by alginit, test 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 r_1+r_2, and in that case test 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 HOME   TOP

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
  ? m_i=[0,-1,0,0;1,0,0,0;0,0,0,-1;0,0,1,0]; \\ quaternion algebra (-1,-1)
  ? m_j=[0,0,-1,0;0,0,0,1;1,0,0,0;0,-1,0,0];
  ? m_k=[0,0,0,-1;0,0,-1,0;0,1,0,0;1,0,0,0];
  ? mt = [matid(4), m_i, m_j, m_k];
  ? A = algtableinit(mt);
  ? algissemisimple(A)
  %9 = 1

The library syntax is GEN algissemisimple(GEN al).

algissimple HOME   TOP

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
  ? m_i=[0,-1,0,0;1,0,0,0;0,0,0,-1;0,0,1,0];
  ? m_j=[0,0,-1,0;0,0,0,1;1,0,0,0;0,-1,0,0];
  ? m_k=[0,0,0,-1;0,0,b,0;0,1,0,0;1,0,0,0];
  ? mt = [matid(4), m_i, m_j, m_k];
  ? 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 sum F_4+F_2
  ? algissimple(A)
  %13 = 0

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

algissplit HOME   TOP

Given a central simple algebra al output by alginit, test 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 r_1+r_2, and in that case test 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).

alglathnf HOME   TOP

Given an algebra al and a square invertible matrix m with size the dimension of al, returns the lattice generated by the columns of m.

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

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

algleftmultable HOME   TOP

Given an element x in al, computes its left multiplication table. If x is given in basis form, returns its multiplication table on the integral basis; if x is given in algebraic form, returns its multiplication table on the basis corresponding to the algebraic form of elements of al. In every case, if x is a t_COL of length n, then the output is a n x n t_MAT. Also accepts a square matrix with coefficients in al.

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

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

algmul HOME   TOP

Given two elements x and y in al, computes their product x*y 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 HOME   TOP

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 e_1 = 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 HOME   TOP

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 HOME   TOP

Given an element x in al, computes its norm. If al is a table algebra output by algtableinit, returns the absolute norm of x, which is an element of 𝔽_p of ℚ; if al is a central simple algebra output by alginit, 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

Also accepts a square matrix with coefficients in al.

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

algpoleval HOME   TOP

Given an element b in al and a polynomial T in K[X], computes T(b) in al.

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

algpow HOME   TOP

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 HOME   TOP

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 HOME   TOP

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

algradical HOME   TOP

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), generated by (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 ℚ, generated by I_2, 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 HOME   TOP

Given a central simple algebra al output by alginit, return 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 r_1 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 HOME   TOP

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 HOME   TOP

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;
  ? {m_i = [0,a,0,0;
            1,0,0,0;
            0,0,0,a;
            0,0,1,0];}
  ? {m_j = [0, 0,b, 0;
            0, 0,0,-b;
            1, 0,0, 0;
            0,-1,0, 0];}
  ? {m_k = [0, 0,0,-a*b;
            0, 0,b,   0;
            0,-a,0,   0;
            1, 0,0,   0];}
  ? mt = [matid(4), m_i, m_j, m_k];
  ? A = alginit(nf,mt,'x);
  ? M = algrelmultable(A);
  ? M[2] == m_i
  %8 = 1
  ? M[3] == m_j
  %9 = 1
  ? M[4] == m_k
  %10 = 1

The library syntax is GEN algrelmultable(GEN al).

algsimpledec HOME   TOP

al being the output of algtableinit representing a semisimple algebra, returns a t_VEC [al_1,al_2,...,al_n] such that al is isomorphic to the direct sum of the simple algebras al_i. When flag = 1, each component is instead a t_VEC [al_i,proj_i,lift_i] where proj_i and lift_i are matrices respectively representing the projection map on the i-th factor and a section of it. 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 ℚ.

Warning. The images of the lift_i are not guaranteed to form a direct sum.

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

algsplittingdata HOME   TOP

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;
  ? {m_i = [0,a,0,0;
            1,0,0,0;
            0,0,0,a;
            0,0,1,0];}
  ? {m_j = [0, 0,b, 0;
            0, 0,0,-b;
            1, 0,0, 0;
            0,-1,0, 0];}
  ? {m_k = [0, 0,0,-a*b;
            0, 0,b,   0;
            0,-a,0,   0;
            1, 0,0,   0];}
  ? mt = [matid(4), m_i, m_j, m_k];
  ? 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 HOME   TOP

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;
  {m_i = [0,a,0,0;
         1,0,0,0;
         0,0,0,a;
         0,0,1,0];}
  {m_j = [0, 0,b, 0;
         0, 0,0,-b;
         1, 0,0, 0;
         0,-1,0, 0];}
  {m_k = [0, 0,0,-a*b;
         0, 0,b,   0;
         0,-a,0,   0;
         1, 0,0,   0];}
  mt = [matid(4), m_i, m_j, m_k];
  A = alginit(nf,mt,'x);
  algsplittingfield(A).pol
  %8 = x^2 - y

The library syntax is GEN algsplittingfield(GEN al).

algsplittingmatrix HOME   TOP

A central simple algebra al output by alginit contains data describing an isomorphism φ : A ⨂ _K L → M_d(L), where d is the degree of the algebra and L is an extension of L with [L:K] = d. Returns the matrix φ(x).

  ? A = alginit(nfinit(y), [-1,-1]);
  ? algsplittingmatrix(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)]

Also accepts matrices with coefficients in al.

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

algsqr HOME   TOP

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 HOME   TOP

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 HOME   TOP

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

  ? 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)
  %5 = 2

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

algtableinit HOME   TOP

Initialize 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 (e_1 = 1, e_2,..., e_n); the table is given as a t_VEC of n matrices in M_n(K), giving the left multiplication by the basis elements e_i, in the given basis. Assumes that e_1 = 1, that K e_1⨁ ...⨁ K e_n] 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 I_2, 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) \\ = (I_2)
  %7 =
  [1]
  
  [0]
  
  [0]

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

algtensor HOME   TOP

Given two algebras al1 and al2, computes their tensor product. For table algebras output by algtableinit, the flag maxord is ignored. For central simple algebras output by alginit, 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).

algtrace HOME   TOP

Given an element x in al, computes its trace. If al is a table algebra output by algtableinit, returns the absolute trace of x, which is an element of 𝔽_p or ℚ; if al is the output of alginit, 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

Also accepts a square matrix with coefficients in al.

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

algtype HOME   TOP

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;
  ?  {m_i = [0,a,0,0;
             1,0,0,0;
             0,0,0,a;
             0,0,1,0];}
  ?  {m_j = [0, 0,b, 0;
             0, 0,0,-b;
             1, 0,0, 0;
             0,-1,0, 0];}
  ?  {m_k = [0, 0,0,-a*b;
             0, 0,b,   0;
             0,-a,0,   0;
             1, 0,0,   0];}
  ?  mt = [matid(4), m_i, m_j, m_k];
  ?  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).