Pari/GP Reference Documentation | Contents - Global index - GP keyboard shortcuts |
Algebra definitions GP conventions for algebras Lattices in algebras Orders in algebras algadd algalgtobasis algaut algb algbasis algbasistoalg algcenter algcentralproj algchar algcharpoly algdegree algdim algdisc algdivl algdivr alggroup alggroupcenter alghasse alghassef alghassei algindex alginit alginv alginvbasis algisassociative algiscommutative algisdivision algisdivl algisinv algisramified algissemisimple algissimple algissplit alglatadd alglatcontains alglatelement alglathnf alglatindex alglatinter alglatlefttransporter alglatmul alglatrighttransporter alglatsubset algmakeintegral algmul algmultable algneg algnorm algpoleval algpow algprimesubalg algquotient algradical algramifiedplaces algrandom algrelmultable algsimpledec algsplit algsplittingdata algsplittingfield algsqr algsub algsubalg algtableinit algtensor algtomatrix algtrace algtype | |
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 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 d2 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-1xiL with xd = 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.
Let H be the Hamilton quaternion algebra, that is the 4-dimensional
algebra over ℝ with basis 1,i,j,ij and multiplication given
by i2 = j2 = -1 and ji = -ij, which is also the cyclic
algebra (ℂ/ℝ,z
| |
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
| |
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 * I is an integral nonsingular upper-triangular matrix representing a sublattice of 𝒪0 expressed on the integral basis, and
* t ∈ ℚ > 0 is a
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
| |
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
Warning. The coefficients in the decomposition A = ⨁ i = 0d-1xiL are not the same as those in the decomposition A = ⨁ i = 0d-1Lxi! The i-th coefficients are related by conjugating by xi, 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
An element of the Hamilton quaternion algebra H can be represented as a
| |
algadd({al}, x, y) | |
Given two elements x and y in al (Hamilton quaternions if omitted), computes their sum x+y in the algebra al.
? A = alginit(nfinit(y),[-1,1]); ? algadd(A,[1,x]~,[1,2,3,4]~) % = [2, 1, 1, 6]~ ? algadd(,sqrt(2+I),[-1,0,1,2]~) % = [0.4553466902, 0.3435607497, 1, 2]~ Also accepts matrices with coefficients in al.
If x and y are given in the same format, then one should simply use
The library syntax is
| |
algalgtobasis(al, x) | |
Given an element x in the central simple algebra al output
by
? 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
| |
algaut(al) | |
Given a cyclic algebra al = (L/K,σ,b) output by
? 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
| |
algb(al) | |
Given a cyclic algebra al = (L/K,σ,b) output by
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
| |
algbasis(al) | |
Given a central simple algebra al output by
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
| |
algbasistoalg(al, x) | |
Given an element x in the central simple algebra al output
by
? 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
| |
algcenter(al) | |
If al is a table algebra output by
A simple example: the 2 x 2 upper triangular matrices over ℚ,
generated by I2, a =
? 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
| |
algcentralproj(al, z, {maps = 0}) | |
Given a table algebra al output by A simple example: 𝔽2 x 𝔽4, generated by 1 = (1,1), e = (1,0) and x such that x2+x+1 = 0. We have e2 = e, x2 = 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
| |
algchar(al) | |
Given an algebra al output by
? 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
| |
algcharpoly({al}, b, {v = 'x}, {abs = 0}) | |
Given an element b in al (Hamilton quaternions if omitted), returns
its characteristic polynomial as a polynomial in the variable v. If al
is a table algebra output by
? 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 ? algcharpoly(,[sqrt(2),-1,0,Pi]~) \\ Hamilton quaternions %10 = x^2 - 2.8284271247*x + 12.8696044010 Also accepts a square matrix with coefficients in al.
The library syntax is
| |
algdegree(al) | |
Given a central simple algebra al output by
? nf = nfinit(y^3-y+1); ? A = alginit(nf, [-1,-1]); ? algdegree(A) %3 = 2
The library syntax is
| |
algdim(al, {abs = 0}) | |
If al is a table algebra output by
? nf = nfinit(y^3-y+1); ? A = alginit(nf, [-1,-1]); ? algdim(A) %3 = 4 ? algdim(A,1) %4 = 12 ? C = alginit(I,0); \\ complex numbers as a real algebra ? algdim(C,1) %6 = 2
The library syntax is
| |
algdisc(al) | |
Given a central simple algebra al output by
? 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
| |
algdivl({al}, x, y) | |
Given two elements x and y in al (Hamilton quaternions if omitted), 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 square matrices with coefficients in al.
? A = alginit(nfinit(y),[-1,1]); ? x = [1,1]~; algisinv(A,x) % = 0 ? z = algmul(A,x,algrandom(A,2)) % = [0, 0, 0, 8]~ ? algdivl(A,x,z) % = [4, 4, 0, 0]~
The library syntax is
| |
algdivr({al}, x, y) | |
Given two elements x and y in al (Hamilton quaternions if omitted), returns xy-1. Also accepts square matrices with coefficients in al.
The library syntax is
| |
alggroup(gal, {p = 0}) | |
Initializes the group algebra K[G] over K = ℚ (p omitted) or 𝔽p
where G is the underlying group of the 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
| |
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
? 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
| |
alghasse(al, {pl}) | |
Given a central simple algebra al output by
? 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 ? H = alginit(1.,1/2); \\ Hamilton quaternion algebra ? alghasse(H) %8 = 1/2
The library syntax is
| |
alghassef(al) | |
Given a central simple algebra al output by
? 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
| |
alghassei(al) | |
Given a central simple algebra al output by
? nf = nfinit(y^2-5); ? A = alginit(nf, [-1,y]); ? alghassei(A) %3 = Vecsmall([1, 0])
The library syntax is
| |
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
| |
alginit(B, C, {v}, {flag = 3}) | |
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
{ 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
? 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 σ:ζ
* (quaternion algebra, special case of the above) B is an
? Q = nfinit(y); A = alginit(Q, [-1,-1]); \\ (-1,-1)ℚ
* (algebra/K defined by local Hasse invariants)
B is an 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
(
? 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
* (algebras over ℝ) If B is a In all cases over a number field, this function factors various discriminants and computes a maximal order for the algebra by default, which may require a lot of time. This can be controlled by flag, whose binary digits mean: * 1: compute a maximal order. * 2: fully factor the discriminants instead of using a lazy factorisation. If this digit of flag is set to 0, the local Hasse invariants are not computed. If only a partial factorisation is known, the computed order is only guaranteed to be maximal at the known prime factors.
The pari object representing such an algebra A is a
* A splitting field L of A of the same degree over K as A, in
* The Hasse invariants at the real places of K, accessed with
* The Hasse invariants of A at the finite primes of K that ramify in
the natural order of A, accessed with
* A basis of an order 𝒪0 expressed on the basis of the natural
order, accessed with
* A basis of the natural order expressed on the basis of 𝒪0,
accessed with
* The left multiplication table of 𝒪0 on the previous basis,
accessed with
* The characteristic of A (always 0), accessed with * 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
* If A was constructed as a cyclic algebra (L/K,σ,b), the
element b, accessed with
* If A was constructed with its multiplication table mt over K,
the
* If A was constructed with its multiplication table mt over K,
a
The library syntax is
| |
alginv({al}, x) | |
Given an element x in al (Hamilton quaternions if omitted), 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]~ ? alginv(,[1,0,Pi,sqrt(2)]~) \\ Hamilton quaternions %3 = [0.0777024661, 0, -0.2441094967, -0.1098878814]~ Also accepts square matrices with coefficients in al.
The library syntax is
| |
alginvbasis(al) | |
Given an central simple algebra al output by
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
| |
algisassociative(mt, p = 0) | |
Returns 1 if the multiplication table
? 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
The library syntax is
| |
algiscommutative(al) | |
al being a table algebra output by
? 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
| |
algisdivision(al, {pl}) | |
Given a central simple algebra al output by
? 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
| |
algisdivl({al}, x, y, {&z}) | |
Given two elements x and y in al (Hamilton quaternions if omitted), 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 square matrices with coefficients in al.
The library syntax is
| |
algisinv({al}, x, {&ix}) | |
Given an element x in al (Hamilton quaternions if omitted), 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 square matrices with coefficients in al.
The library syntax is
| |
algisramified(al, {pl}) | |
Given a central simple algebra al output by
? 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
| |
algissemisimple(al) | |
al being a table algebra output by
? 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
| |
algissimple(al, {ss = 0}) | |
al being a table algebra output by
? 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
| |
algissplit(al, {pl}) | |
Given a central simple algebra al output by
? 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
| |
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
| |
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
? 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
| |
alglatelement(al, lat, c) | |
Given an algebra al, a lattice lat and a
? 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
| |
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
? 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
| |
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
| |
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
| |
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
| |
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
| |
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
| |
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
| |
algmakeintegral(mt, {maps = 0}) | |
mt being a multiplication table over ℚ in the same format as the
input of
? 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
| |
algmul({al}, x, y) | |
Given two elements x and y in al (Hamilton quaternions if omitted), 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, 3, 5, -4]~ ? algmul(,[1,2,3,4]~,sqrt(I)) \\ Hamilton quaternions % = [-0.7071067811, 2.1213203435, 4.9497474683, 0.7071067811]~ Also accepts matrices with coefficients in al.
The library syntax is
| |
algmultable(al) | |
Returns a multiplication table of al over its
prime subfield (ℚ or 𝔽p) or over ℝ for real algebras, as a
? 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] ? H = alginit(1.,1/2); \\ Hamilton quaternions ? algmultable(H)[3] \\ multiplication by j %7 = [0 0 -1 0] [0 0 0 1] [1 0 0 0] [0 -1 0 0]
The library syntax is
| |
algneg({al}, x) | |
Given an element x in al, computes its opposite -x in the algebra al (Hamilton quaternions if omitted).
? 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
| |
algnorm({al}, x, {abs = 0}) | |
Given an element x in al (Hamilton quaternions if omitted),
computes its norm. If al is a table algebra output by
? 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 ? algnorm(,[0,sqrt(3),0,sqrt(2)]~) \\ Hamilton quaternions %10 = 5.0000000000 Also accepts a square matrix with coefficients in al.
The library syntax is
| |
algpoleval({al}, T, b) | |
Given an element b in al (Hamilton quaternions if omitted) and a
polynomial T in K[X], computes T(b) in al. Also accepts as input a
? 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 ? algpoleval(,polcyclo(8),[1,0,0,1]~/sqrt(2)) \\ Hamilton quaternions %8 = [0.E-38, 0, 0, 0.E-38]~
The library syntax is
| |
algpow({al}, x, n) | |
Given an element x in al (Hamilton quaternions if omitted) and an integer n, computes the power xn in the algebra al.
? A = alginit(nfinit(y), [-1,-1]); ? algpow(A,[1,1,0,0]~,7) %2 = [8, -8, 0, 0]~ ? algpow(,[1,2,3,sqrt(3)]~,-3) \\ Hamilton quaternions % = [-0.0095664563, 0.0052920822, 0.0079381233, 0.0045830776]~ Also accepts a square matrix with coefficients in al.
The library syntax is
| |
algprimesubalg(al) | |
al being the output of
? 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
| |
algquotient(al, I, {maps = 0}) | |
al being a table algebra output by
? 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
| |
algradical(al) | |
al being a table algebra output by Here is an example with A = ℚ[x]/(x2), 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 =
? 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
| |
algramifiedplaces(al) | |
Given a central simple algebra al output by
? nf = nfinit(y^2-5); ? A = alginit(nf, [-1,y]); ? algramifiedplaces(A) %3 = [1, [2, [2, 0]~, 1, 2, 1]]
The library syntax is
| |
algrandom({al}, b) | |
Given an algebra al and a nonnegative integer b, returns a random element in al with coefficients in [-b,b].
? al = alginit(nfinit(y),[-1,-1]); ? algrandom(al,3) % = [2, 0, 3, -1]~
If al is an algebra over ℝ (Hamilton quaternions if omitted) and
b is a positive
? algrandom(,1.) % = [-0.1806334680, -0.2810504190, 0.5011479961, 0.9498643737]~
The library syntax is
| |
algrelmultable(al) | |
Given a central simple algebra al output by
? 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
| |
algsimpledec(al, {maps = 0}) | |
al being the output of
The library syntax is
| |
algsplit(al, {v = 'x}) | |
If al is a table algebra over 𝔽p output by
* map is a
* mapi is an N x N matrix with 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,
? 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
| |
algsplittingdata(al) | |
Given a central simple algebra al output by * an element t of al such that L = K(t) is a maximal subfield of al;
* a matrix
* a matrix
? 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
| |
algsplittingfield(al) | |
Given a central simple algebra al output by
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
| |
algsqr({al}, x) | |
Given an element x in al (Hamilton quaternions if omitted), computes its square x2 in the algebra al.
? A = alginit(nfinit(y), [-1,-1]); ? algsqr(A,[1,0,2,0]~) %2 = [-3, 0, 4, 0]~ ? algsqr(,[0,0,0,Pi]~) \\ Hamilton quaternions %3 = [-9.8696044010, 0, 0, 0]~ Also accepts a square matrix with coefficients in al.
The library syntax is
| |
algsub({al}, x, y) | |
Given two elements x and y in al (Hamilton quaternions if omitted), 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.
If x and y are given in the same format, then one should simply use
The library syntax is
| |
algsubalg(al, B) | |
al being a table algebra output by 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
| |
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 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
* The characteristic of A, accessed with
* The multiplication table of A, accessed with * The traces of the elements of the basis.
A simple example: the 2 x 2 upper triangular matrices over ℚ,
generated by I2, a =
? 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
| |
algtensor(al1, al2, {flag = 3}) | |
Given two algebras al1 and al2, computes their tensor
product. flag has the same meaning as in 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
| |
algtomatrix({al}, x, {abs = 0}) | |
Given an element x in al (Hamilton quaternions if omitted),
returns the image of x under a homomorphism to a matrix algebra. If
al is a table algebra output by
? 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] ? algtomatrix(,[1,2,3,4]~) \\ Hamilton quaternions %4 = [1 + 2*I -3 - 4*I] [3 - 4*I 1 - 2*I] ? algtomatrix(,I,1) %5 = [0 -1 0 0] [1 0 0 0] [0 0 0 -1] [0 0 1 0] Also accepts matrices with coefficients in al.
The library syntax is
| |
algtrace({al}, x, {abs = 0}) | |
Given an element x in al (Hamilton quaternions if omitted),
computes its trace. If al is a table algebra output by
? 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 ? algtrace(,[1.,2,3,4]~) \\ Hamilton quaternions %10 = 2.0000000000 ? algtrace(,[1.,2,3,4]~,0) %11 = 4.0000000000 Also accepts a square matrix with coefficients in al.
The library syntax is
| |
algtype(al) | |
Given an algebra al output by * 0: not a valid algebra.
* 1: table algebra output by
* 2: central simple algebra output by
* 3: central simple algebra output by * 4: division algebra over ℝ (ℝ, ℂ or Hamilton quaternion algebra H).
? 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 ? H = alginit(1.,1/2); ? algtype(H) %16 = 4
The library syntax is
| |