libPARI - Algebraic Number Theory
None of the following routines thoroughly check their input: they distinguish between bona fide structures as output by PARI routines, but designing perverse data will easily fool them. To give an example, a square matrix will be interpreted as an ideal even though the Z-module generated by its columns may not be an Z_K
-module (i.e. the expensive nfisideal
routine will not be called).
long
nftyp(GEN x)
. Returns the type of number field structure stored in x
, typ_NF
, typ_BNF
, or typ_BNR
. Other answers are possible, meaning x
is not a number field structure.
GEN
get_nf(GEN x, long *t)
. Extract an nf structure from x
if possible and return it, otherwise return NULL
. Sets t
to the nftyp
of x
in any case.
GEN
get_bnf(GEN x, long *t)
. Extract a bnf
structure from x
if possible and return it, otherwise return NULL
. Sets t
to the nftyp
of x
in any case.
GEN
get_nfpol(GEN x, GEN *nf)
try to extract an nf structure from x
, and sets *nf
to NULL
(failure) or to the nf. Returns the (monic, integral) polynomial defining the field.
GEN
get_bnfpol(GEN x, GEN *bnf, GEN *nf)
try to extract a bnf and an nf structure from x
, and sets *bnf
and *nf
to NULL
(failure) or to the corresponding structure. Returns the (monic, integral) polynomial defining the field.
GEN
checknf(GEN x)
if an nf structure can be extracted from x
, return it; otherwise raise an exception. The more general get_nf
is often more flexible.
GEN
checkbnf(GEN x)
if an bnf structure can be extracted from x
, return it; otherwise raise an exception. The more general get_bnf
is often more flexible.
GEN
checkbnf_i(GEN bnf)
same as checkbnf
but return NULL
instead of raising an exception.
void
checkbnr(GEN bnr)
Raise an exception if the argument is not a bnr structure.
GEN
checknf_i(GEN nf)
same as checknf
but return NULL
instead of raising an exception.
void
checkbnrgen(GEN bnr)
Raise an exception if the argument is not a bnr structure, complete with explicit generators for the ray class group. This is normally useless and checkbnr
should be instead, unless you are absolutely certain that the generators will be needed at a later point, and you are about to embark in a costly intermediate computation. PARI functions do check that generators are present in bnr before accessing them: they will raise an error themselves; many functions that may require them, e.g. bnrconductor
, often do not actually need them.
void
checkrnf(GEN rnf)
Raise an exception if the argument is not an rnf structure.
int
checkrnf_i(GEN rnf)
same as checkrnf
but return 0
on failure and 1
on success.
void
checkbid(GEN bid)
Raise an exception if the argument is not a bid structure.
GEN
checkbid_i(GEN bid)
same as checkbid
but return NULL
instead of raising an exception and return bid
on success.
GEN
checkgal(GEN x)
if a galoisinit structure can be extracted from x
, return it; otherwise raise an exception.
void
checksqmat(GEN x, long N)
check whether x
is a square matrix of dimension N
. May be used to check for ideals if N
is the field degree.
void
checkprid(GEN pr)
Raise an exception if the argument is not a prime ideal structure.
int
checkprid_i(GEN pr)
same as checkprid
but return 0
instead of raising an exception and return 1
on success.
int
is_nf_factor(GEN F)
return 1
if F
is an ideal factorization and 0
otherwise.
int
is_nf_extfactor(GEN F)
return 1
if F
is an extended ideal factorization (allowing 0
or negative exponents) and 0
otherwise.
GEN
get_prid(GEN ideal)
return the underlying prime ideal structure if one can be extracted from ideal
(ideal or extended ideal), and return NULL
otherwise.
void
checkabgrp(GEN v)
Raise an exception if the argument is not an abelian group structure, i.e. a t_VEC
with either 2
or 3
entries: [N,
cyc]
or [N,
cyc,
gen]
.
GEN
abgrp_get_no(GEN x)
extract the cardinality N
from an abelian group structure.
GEN
abgrp_get_cyc(GEN x)
extract the elementary divisors cyc from an abelian group structure.
GEN
abgrp_get_gen(GEN x)
extract the generators gen from an abelian group structure.
void
checkmodpr(GEN modpr)
Raise an exception if the argument is not a modpr
structure (from nfmodprinit
).
GEN
get_modpr(GEN x)
return x
if it is a modpr
structure and NULL
otherwise.
GEN
checknfelt_mod(GEN nf, GEN x, const char *s)
given an nf structure nf
and a t_POLMOD
x
, return the attached polynomial representative (shallow) if x
and nf
are compatible. Raise an exception otherwise. Set s
to the name of the caller for a meaningful error message.
void
check_ZKmodule(GEN x, const char *s)
check whether x
looks like Z_K
-module (a pair [A,I]
, where A
is a matrix and I
is a list of ideals; A
has as many columns as I
has elements. Otherwise raises an exception. Set s
to the name of the caller for a meaningful error message.
long
idealtyp(GEN *ideal, GEN *fa)
The input is ideal
, a pointer to an ideal (or extended ideal), which is usually modified, fa
being set as a side-effect. Returns the type of the underlying ideal among id_PRINCIPAL
(a number field element), id_PRIME
(a prime ideal) id_MAT
(an ideal in matrix form).
If ideal
pointed to an ideal, set fa
to NULL
, and possibly simplify ideal
(for instance the zero ideal is replaced by gen_0
). If it pointed to an extended ideal, replace ideal
by the underlying ideal and set fa
to the factorization matrix component.
nf
structureThese functions expect a true nf argument attached to a number field K =
Q[x]/(T)
, e.g. a bnf will not work. Let n = [K:
Q]
be the field degree.
GEN
nf_get_pol(GEN nf)
returns the polynomial T
(monic, in Z[x]
).
long
nf_get_varn(GEN nf)
returns the variable number of the number field defining polynomial.
long
nf_get_r1(GEN nf)
returns the number of real places r_1
.
long
nf_get_r2(GEN nf)
returns the number of complex places r_2
.
void
nf_get_sign(GEN nf, long *r1, long *r2)
sets r_1
and r_2
to the number of real and complex places respectively. Note that r_1+2r_2
is the field degree.
long
nf_get_degree(GEN nf)
returns the number field degree, n = r_1 + 2r_2
.
GEN
nf_get_disc(GEN nf)
returns the field discriminant.
GEN
nf_get_index(GEN nf)
returns the index of T
, i.e. the index of the order generated by the power basis (1,x,...,x^{n-1})
in the maximal order of K
.
GEN
nf_get_zk(GEN nf)
returns a basis (w_1,w_2,...,w_n)
for the maximal order of K
. Those are polynomials in Q[x]
of degree < n
; it is guaranteed that w_1 = 1
.
GEN
nf_get_invzk(GEN nf)
returns the matrix (m_{i,j})\in M_n(
Z)
giving the power basis (x^i)
in terms of the (w_j)
, i.e such that x^{j-1} =
sum_{i = 1}^n m_{i,j} w_i
for all 1 <= j <= n
; since w_1 = 1 = x^0
, we have m_{i,1} =
delta_{i,1}
for all i
. The conversion functions in the algtobasis
family essentially amount to a left multiplication by this matrix.
GEN
nf_get_roots(GEN nf)
returns the r_1
real roots of the polynomial defining the number fields: first the r_1
real roots (as t_REAL
s), then the r_2
representatives of the pairs of complex conjugates.
GEN
nf_get_allroots(GEN nf)
returns all the complex roots of T
: first the r_1
real roots (as t_REAL
s), then the r_2
pairs of complex conjugates.
GEN
nf_get_M(GEN nf)
returns the (r_1+r_2) x n
matrix M
giving the embeddings of K
: M[i,j]
contains w_j(
alpha_i)
, where alpha_i
is the i
-th element of nf_get_roots(nf)
. In particular, if v
is an n
-th dimensional t_COL
representing the element sum_{i = 1}^n v[i] w_i
of K
, then RgM_RgC_mul(M,v)
represents the embeddings of v
.
GEN
nf_get_G(GEN nf)
returns a n x n
real matrix G
such that Gv.Gv = T_2(v)
, where v
is an n
-th dimensional t_COL
representing the element sum_{i = 1}^n v[i] w_i
of K
and T_2
is the standard Euclidean form on K\otimes
R, i.e. T_2(v) =
sum_{
sigma} |
sigma(v)|^2
, where sigma runs through all n
complex embeddings of K
.
GEN
nf_get_roundG(GEN nf)
returns a rescaled version of G
, rounded to nearest integers, specifically RM_round_maxrank
(G)
.
GEN
nf_get_ramified_primes(GEN nf)
returns the vector of ramified primes.
GEN
nf_get_Tr(GEN nf)
returns the matrix of the Trace quadratic form on the basis (w_1,...,w_n)
: its (i,j)
entry is Tr w_i w_j
.
GEN
nf_get_diff(GEN nf)
returns the primitive part of the inverse of the above Trace matrix.
long
nf_get_prec(GEN nf)
returns the precision (in words) to which the nf was computed.
bnf
structureThese functions expect a true bnf argument, e.g. a bnr will not work.
GEN
bnf_get_nf(GEN bnf)
returns the underlying nf.
GEN
bnf_get_clgp(GEN bnf)
returns the class group in bnf, which is a 3
-component vector [h,
cyc,
gen]
.
GEN
bnf_get_cyc(GEN bnf)
returns the elementary divisors of the class group (cyclic components) [d_1,..., d_k]
, where d_k | ... | d_1
.
GEN
bnf_get_gen(GEN bnf)
returns the generators [g_1,...,g_k]
of the class group. Each g_i
has order d_i
, and the full module of relations between the g_i
is generated by the d_ig_i = 0
.
GEN
bnf_get_no(GEN bnf)
returns the class number.
GEN
bnf_get_reg(GEN bnf)
returns the regulator.
GEN
bnf_get_logfu(GEN bnf)
returns (complex floating point approximations to) the logarithms of the complex embeddings of our system of fundamental units.
GEN
bnf_get_fu(GEN bnf)
returns the fundamental units. Raise an error if the bnf does not contain units in algebraic form.
GEN
bnf_get_fu_nocheck(GEN bnf)
as bnf_get_fu
without checking whether units are present. Do not use this unless you initialize the bnf yourself!
GEN
bnf_get_tuU(GEN bnf)
returns a generator of the torsion part of Z_K^*
.
long
bnf_get_tuN(GEN bnf)
returns the order of the torsion part of Z_K^*
, i.e. the number of roots of unity in K
.
bnr
structureThese functions expect a true bnr argument.
GEN
bnr_get_bnf(GEN bnr)
returns the underlying bnf.
GEN
bnr_get_nf(GEN bnr)
returns the underlying nf.
GEN
bnr_get_clgp(GEN bnr)
returns the ray class group.
GEN
bnr_get_no(GEN bnr)
returns the ray class number.
GEN
bnr_get_cyc(GEN bnr)
returns the elementary divisors of the ray class group (cyclic components) [d_1,..., d_k]
, where d_k | ... | d_1
.
GEN
bnr_get_gen(GEN bnr)
returns the generators [g_1,...,g_k]
of the ray class group. Each g_i
has order d_i
, and the full module of relations between the g_i
is generated by the d_ig_i = 0
. Raise a generic error if the bnr does not contain the ray class group generators.
GEN
bnr_get_gen_nocheck(GEN bnr)
as bnr_get_gen
without checking whether generators are present. Do not use this unless you initialize the bnr yourself!
GEN
bnr_get_bid(GEN bnr)
returns the bid attached to the bnr modulus.
GEN
bnr_get_mod(GEN bnr)
returns the modulus attached to the bnr.
rnf
structureThese functions expect a true rnf argument, attached to an extension L/K
, K =
Q[y]/(T)
, L = K[x]/(P)
.
long
rnf_get_degree(GEN rnf)
returns the relative degree [L:K]
.
long
rnf_get_absdegree(GEN rnf)
returns the absolute degree [L:
Q]
.
long
rnf_get_nfdegree(GEN rnf)
returns the degree of the base field [K:
Q]
.
GEN
rnf_get_nf(GEN rnf)
returns the base field K
, an nf structure.
GEN
rnf_get_nfpol(GEN rnf)
returns the polynomial T
defining the base field K
.
long
rnf_get_nfvarn(GEN rnf)
returns the variable y
attached to the base field K
.
void
rnf_get_nfzk(GEN rnf, GEN *b, GEN *cb)
returns the integer basis of the base field K
.
GEN
rnf_get_pol(GEN rnf)
returns the relative polynomial defining L/K
.
long
rnf_get_varn(GEN rnf)
returns the variable x
attached to L
.
GEN
rnf_get_zk(GEN nf)
returns the relative integer basis generating Z_L
as a Z_K
-module, as a pseudo-matrix (A,I)
in HNF.
GEN
rnf_get_disc(GEN rnf)
is the output [
d, s]
of rnfdisc
.
GEN
rnf_get_idealdisc(GEN rnf)
is the ideal discriminant d from rnfdisc
.
GEN
rnf_get_index(GEN rnf)
is the index ideal f
GEN
rnf_get_polabs(GEN rnf)
returns an absolute polynomial defining L/
Q.
GEN
rnf_get_alpha(GEN rnf)
a root alpha of the polynomial defining the base field, modulo polabs
(cf. rnfequation
)
GEN
rnf_get_k(GEN rnf)
a small integer k
such that theta =
beta + k
alpha is a root of polabs
, where beta is a root of pol
and alpha a root of the polynomial defining the base field, as in rnf_get_alpha
(cf. also rnfequation
).
GEN
rnf_get_invzk(GEN rnf)
contains A^{-1}
, where (A,I)
is the chosen pseudo-basis for Z_L
over Z_K
.
GEN
rnf_get_map(GEN rnf)
returns technical data attached to the map K\to L
. Currently, this contains data from rnfequation
, as well as the polynomials T
and P
.
bid
structureThese functions expect a true bid argument, attached to a modulus I = I_0 I_ oo
in a number field K
.
GEN
bid_get_mod(GEN bid)
returns the modulus attached to the bid.
GEN
bid_get_grp(GEN bid)
returns the Abelian group attached to (
Z_K/I)^*
.
GEN
bid_get_ideal(GEN bid)
return the finite part I_0
of the bid modulus (an integer ideal).
GEN
bid_get_arch(GEN bid)
return the Archimedean part I_ oo
of the bid modulus as a vector of real places in vec01 format, see "Label se:signatures".
GEN
bid_get_archp(GEN bid)
return the Archimedean part I_ oo
of the bid modulus, as a vector of real places in indices format see "Label se:signatures".
GEN
bid_get_fact(GEN bid)
returns the ideal factorization I_0 =
prod_i
p_i^{e_i}
.
bid_get_ideal
(
bid)
, via idealfactor
.
GEN
bid_get_no(GEN bid)
returns the cardinality of the group (
Z_K/I)^*
.
GEN
bid_get_cyc(GEN bid)
returns the elementary divisors of the group (
Z_K/I)^*
(cyclic components) [d_1,..., d_k]
, where d_k | ... | d_1
.
GEN
bid_get_gen(GEN bid)
returns the generators of (
Z_K/I)^*
contained in bid. Raise a generic error if bid does not contain generators.
GEN
bid_get_gen_nocheck(GEN bid)
as bid_get_gen
without checking whether generators are present. Do not use this unless you initialize the bid yourself!
GEN
bid_get_sprk(GEN bid)
return a list of structures attached to the (
Z_K/
p^e)^*
where p^e
divides I_0
exactly.
GEN
bid_get_sarch(GEN bid)
return the structure attached to (
Z_K/I_ oo )^*
, by nfarchstar
.
GEN
bid_get_U(GEN bid)
return the matrix with integral coefficients relating the local generators (from chinese remainders) to the global SNF generators (bid.gen
).
GEN
bid_get_ind(GEN bid)
return a t_VECSMALL
v
of indices used while converting from local generators to the global generators: v[i]
is the number of generators used to describe (
Z_K /
prod_{j < i}
p_j^{e_j})^*
.
If the required data is not part of the structure, it is computed then inserted, and the new value is returned.
These functions expect a bnf
argument:
GEN
bnf_build_cycgen(GEN bnf)
the bnf contains generators [g_1,...,g_k]
of the class group, each with order d_i
. Then g_i^{d_i} = (x_i)
is a principal ideal. This function returns the x_i
as a factorization matrix (famat
) giving the element in factored form as a product of S
-units.
GEN
bnf_build_matalpha(GEN bnf)
the class group was computed using a factorbase S
of prime ideals p_i
, i <= r
. They satisfy relations of the form prod_j
p_i^{e_{i,j}} = (
alpha_j)
, where the e_{i,j}
are given by the matrices bnf[1]
(W
, singling out a minimal set of generators in S
) and bnf[2]
(B
, expressing the rest of S
in terms of the singled out generators). This function returns the alpha_j
in factored form as a product of S
-units.
GEN
bnf_build_units(GEN bnf)
returns a minimal set of generators for the unit group. The first element is a torsion unit, the others have infinite order.
These functions expect a rnf
argument:
GEN
rnf_build_nfabs(GEN rnf, long prec)
given a rnf structure attched to L/K
, (compute and) return an nf structure attached to L
at precision prec
.
void
rnfcomplete(GEN rnf)
as rnf_build_nfabs
using the precision of K
for prec
.
GEN
rnf_zkabs(GEN rnf)
returns a Z-basis in HNF for Z_L
as a pair [T, v]
, where T
is rnf_get_polabs
(
rnf)
and v
a vector of elements lifted from Q[X]/(T)
. Note that rnf_build_nfabs
essentially applies nfinit
to the output of this function.
GEN
nfnewprec(GEN x, long prec)
. Raise an exception if x
is not a number field structure (nf, bnf or bnr). Otherwise, sets its accuracy to prec
and return the new structure. This is mostly useful with prec
larger than the accuracy to which x
was computed, but it is also possible to decrease the accuracy of x
(truncating relevant components, which may speed up later computations). This routine may modify the original x
(see below).
This routine is straightforward for nf structures, but for the other ones, it requires all principal ideals corresponding to the bnf relations in algebraic form (they are originally only available via floating point approximations). This in turn requires many calls to bnfisprincipal0
, which is often slow, and may fail if the initial accuracy was too low. In this case, the routine will not actually fail but recomputes a bnf from scratch!
Since this process may be very expensive, the corresponding data is cached (as a clone) in the original x
so that later precision increases become very fast. In particular, the copy returned by nfnewprec
also contains this additional data.
GEN
bnfnewprec(GEN x, long prec)
. As nfnewprec
, but extracts a bnf structure form x
before increasing its accuracy, and returns only the latter.
GEN
bnrnewprec(GEN x, long prec)
. As nfnewprec
, but extracts a bnr structure form x
before increasing its accuracy, and returns only the latter.
GEN
nfnewprec_shallow(GEN nf, long prec)
GEN
bnfnewprec_shallow(GEN bnf, long prec)
GEN
bnrnewprec_shallow(GEN bnr, long prec)
Shallow functions underlying the above, except that the first argument must now have the corresponding number field type. I.e. one cannot call nfnewprec_shallow(nf, prec)
if nf
is actually a bnf.
The number field K =
Q[X]/(T)
is represented by an nf
(or bnf
or bnr
structure). An algebraic number belonging to K
is given as
@3* a t_INT
, t_FRAC
or t_POL
(implicitly modulo T
), or
@3* a t_POLMOD
(modulo T
), or
@3* a t_COL
v
of dimension N = [K:
Q]
, representing the element in terms of the computed integral basis (e_i)
, as
sum(i = 1, N, v[i] * nf.zk[i])
The preferred forms are t_INT
and t_COL
of t_INT
. Routines can handle denominators but it is much more efficient to remove denominators first (Q_remove_denom
) and take them into account at the end.
@3Safe routines. The following routines do not assume that their nf
argument is a true nf (it can be any number field type, e.g. a bnf), and accept number field elements in all the above forms. They return their result in t_COL
form.
GEN
nfadd(GEN nf, GEN x, GEN y)
returns x+y
.
GEN
nfsub(GEN nf, GEN x, GEN y)
returns x-y
.
GEN
nfdiv(GEN nf, GEN x, GEN y)
returns x / y
.
GEN
nfinv(GEN nf, GEN x)
returns x^{-1}
.
GEN
nfmul(GEN nf,GEN x,GEN y)
returns xy
.
GEN
nfpow(GEN nf,GEN x,GEN k)
returns x^k
, k
is in Z.
GEN
nfpow_u(GEN nf,GEN x, ulong k)
returns x^k
, k >= 0
.
GEN
nfsqr(GEN nf,GEN x)
returns x^2
.
long
nfval(GEN nf, GEN x, GEN pr)
returns the valuation of x
at the maximal ideal p attached to the prid pr
. Returns LONG_MAX
if x
is 0
.
GEN
nfnorm(GEN nf,GEN x)
absolute norm of x
.
GEN
nftrace(GEN nf,GEN x)
absolute trace of x
.
GEN
nfpoleval(GEN nf, GEN pol, GEN a)
evaluate the t_POL
pol
(with coefficients in nf
) on the algebraic number a
(also in nf
).
GEN
FpX_FpC_nfpoleval(GEN nf, GEN pol, GEN a, GEN p)
evaluate the FpX
pol
on the algebraic number a
(also in nf
).
The following three functions implement trivial functionality akin to Euclidean division for which we currently have no real use. Of course, even if the number field is actually Euclidean, these do not in general implement a true Euclidean division.
GEN
nfdiveuc(GEN nf, GEN a, GEN b)
returns the algebraic integer closest to x / y
. Functionally identical to ground( nfdiv(nf,x,y) )
.
GEN
nfdivrem(GEN nf, GEN a, GEN b)
returns the vector [q,r]
, where
q = nfdiveuc(nf, a, b);
r = nfsub(nf, a, nfmul(nf,q,b)); \\ or r = nfmod(nf,a,b);
GEN
nfmod(GEN nf, GEN a, GEN b)
returns r
such that
q = nfdiveuc(nf, a, b);
r = nfsub(nf, a, nfmul(nf,q,b));
GEN
nf_to_scalar_or_basis(GEN nf, GEN x)
let x
be a number field element. If it is a rational scalar, i.e. can be represented by a t_INT
or t_FRAC
, return the latter. Otherwise returns its basis representation (nfalgtobasis
). Shallow function.
GEN
nf_to_scalar_or_alg(GEN nf, GEN x)
let x
be a number field element. If it is a rational scalar, i.e. can be represented by a t_INT
or t_FRAC
, return the latter. Otherwise returns its lifted t_POLMOD
representation (lifted nfbasistoalg
). Shallow function.
GEN
RgX_to_nfX(GEN nf, GEN x)
let x
be a t_POL
whose coefficients are number field elements; apply nf_to_scalar_or_basis
to each coefficient and return the resulting new polynomial. Shallow function.
GEN
RgM_to_nfM(GEN nf, GEN x)
let x
be a t_MAT
whose coefficients are number field elements; apply nf_to_scalar_or_basis
to each coefficient and return the resulting new matrix. Shallow function.
GEN
RgC_to_nfC(GEN nf, GEN x)
let x
be a t_COL
or t_VEC
whose coefficients are number field elements; apply nf_to_scalar_or_basis
to each coefficient and return the resulting new t_COL
. Shallow function.
@3Unsafe routines. The following routines assume that their nf
argument is a true nf (e.g. a bnf is not allowed) and their argument are restricted in various ways, see the precise description below.
GEN
nfinvmodideal(GEN nf, GEN x, GEN A)
given an algebraic integer x
and a non-zero integral ideal A
in HNF, returns a y
such that xy = 1
modulo A
.
GEN
nfpowmodideal(GEN nf, GEN x, GEN n, GEN ideal)
given an algebraic integer x
, an integer n
, and a non-zero integral ideal A
in HNF, returns an algebraic integer congruent to x^n
modulo A
.
GEN
nfmuli(GEN nf, GEN x, GEN y)
returns x x y
assuming that both x
and y
are either t_INT
s or ZV
s of the correct dimension.
GEN
nfsqri(GEN nf, GEN x)
returns x^2
assuming that x
is a t_INT
or a ZV
of the correct dimension.
GEN
nfC_nf_mul(GEN nf, GEN v, GEN x)
given a t_VEC
or t_COL
v
of elements of K
in t_INT
, t_FRAC
or t_COL
form, multiply it by the element x
(arbitrary form). This is faster than multiplying coordinatewise since pre-computations related to x
(computing the multiplication table) are done only once. The components of the result are in most cases t_COL
s but are allowed to be t_INT
s or t_FRAC
s. Shallow function.
GEN
nfC_multable_mul(GEN v, GEN mx)
same as nfC_nf_mul
, where the argument x
is replaced by its multiplication table mx
.
GEN
zkC_multable_mul(GEN v, GEN x)
same as nfC_nf_mul
, where v
is a vector of algebraic integers, x
is an algebraic integer, and x
is replaced by zk_multable
(x)
.
GEN
zk_multable(GEN nf, GEN x)
given a ZC
x
(implicitly representing an algebraic integer), returns the ZM
giving the multiplication table by x
. Shallow function (the first column of the result points to the same data as x
).
GEN
zk_inv(GEN nf, GEN x)
given a ZC
x
(implicitly representing an algebraic integer), returns the QC
giving the inverse x^{-1}
. Return NULL
if x
is 0
. Not memory clean but safe for gerepileupto
.
GEN
zkmultable_inv(GEN mx)
as zk_inv
, where the argument given is zk_multable
(x)
.
GEN
zkmultable_capZ(GEN mx)
given a non-zero zkmultable mx attached to x \in
Z_K
, return the positive generator of (x)
cap Z.
GEN
zk_scalar_or_multable(GEN nf, GEN x)
given a t_INT
or ZC
x
, returns a t_INT
equal to x
if the latter is a scalar (t_INT
or ZV_isscalar
(x)
is 1) and zk_multable
(
nf,x)
otherwise. Shallow function.
The following routines implement multiplication in a commutative R
-algebra, generated by (e_1 = 1,..., e_n)
, and given by a multiplication table M
: elements in the algebra are n
-dimensional t_COL
s, and the matrix M
is such that for all 1 <= i,j <= n
, its column with index (i-1)n + j
, say (c_k)
, gives e_i.e_j =
sum c_k e_k
. It is assumed that e_1
is the neutral element for the multiplication (a convenient optimization, true in practice for all multiplications we needed to implement). If x
has any other type than t_COL
where an algebra element is expected, it is understood as x e_1
.
GEN
multable(GEN M, GEN x)
given a column vector x
, representing the quantity sum_{i = 1}^N x_i e_i
, returns the multiplication table by x
. Shallow function.
GEN
ei_multable(GEN M, long i)
returns the multiplication table by the i
-th basis element e_i
. Shallow function.
GEN
tablemul(GEN M, GEN x, GEN y)
returns x.y
.
GEN
tablesqr(GEN M, GEN x)
returns x^2
.
GEN
tablemul_ei(GEN M, GEN x, long i)
returns x.e_i
.
GEN
tablemul_ei_ej(GEN M, long i, long j)
returns e_i.e_j
.
GEN
tablemulvec(GEN M, GEN x, GEN v)
given a vector v
of elements in the algebra, returns the x.v[i]
.
The following routines implement naive linear algebra using the black box field mechanism:
GEN
nfM_det(GEN nf, GEN M)
GEN
nfM_inv(GEN nf, GEN M)
GEN
nfM_mul(GEN nf, GEN A, GEN B)
GEN
nfM_nfC_mul(GEN nf, GEN A, GEN B)
Computational algebraic theory performs extensively linear algebra on Z-modules with a natural multiplicative structure (K^*
, fractional ideals in K
, Z_K^*
, ideal class group), thereby raising elements to horrendously large powers. A seemingly innocuous elementary linear algebra operation like C_i\leftarrow C_i - 10000 C_1
involves raising entries in C_1
to the 10000
-th power. Understandably, it is often more efficient to keep elements in factored form rather than expand every such expression. A factorization matrix (or famat) is a two column matrix, the first column containing elements (arbitrary objects which may be repeated in the column), and the second one contains exponents (t_INT
s, allowed to be 0). By abuse of notation, the empty matrix cgetg(1, t_MAT)
is recognized as the trivial factorization (no element, no exponent).
Even though we think of a famat with columns g
and e
as one meaningful object when fully expanded as prod g[i]^{e[i]}
, famats are basically about concatenating information to keep track of linear algebra: the objects stored in a famat need not be operation-compatible, they will not even be compared to each other (with one exception: famat_reduce
). Multiplying two famats just concatenates their elements and exponents columns. In a context where a famat is expected, an object x
which is not of type t_MAT
will be treated as the factorization x^1
. The following functions all return famats:
GEN
famat_mul(GEN f, GEN g)
f
, g
are famat, or objects whose type is not t_MAT
(understood as f^1
or g^1
). Returns fg
. The empty factorization is the neutral element for famat multiplication.
GEN
famat_mul_shallow(GEN f, GEN g)
shallow version of famat_mul
.
GEN
famat_pow(GEN f, GEN n)
n
is a t_INT
. If f
is a t_MAT
, assume it is a famat and return f^n
(multiplies the exponent column by n
). Otherwise, understand it as an element and returns the 1
-line famat f^n
.
GEN
famat_pow_shallow(GEN f, GEN n)
shallow version of famat_pow
.
GEN
famat_mulpow_shallow(GEN f, GEN g, GEN e)
famat
corresponding to f.g^e
. Shallow function.
GEN
famat_sqr(GEN f)
returns f^2
.
GEN
famat_inv(GEN f)
returns f^{-1}
.
GEN
famat_inv_shallow(GEN f)
shallow version of famat_inv
.
GEN
famat_Z_gcd(GEN M, GEN n)
restrict the famat M
to the prime power dividing n
.
GEN
to_famat(GEN x, GEN k)
given an element x
and an exponent k
, returns the famat x^k
.
GEN
to_famat_shallow(GEN x, GEN k)
same, as a shallow function.
Note that it is trivial to break up a famat into its two constituent columns: gel(f,1)
and gel(f,2)
are the elements and exponents respectively. Conversely, mkmat2
builds a (shallow) famat from two t_COL
s of the same length.
The last two functions makes an assumption about the elements: they must be regular algebraic numbers (not famats) over a given number field:
GEN
famat_reduce(GEN f)
given a famat f
, returns a famat g
without repeated elements or 0 exponents, such that the expanded forms of f
and g
would be equal. Shallow function.
GEN
ZM_famat_limit(GEN f, GEN limit)
given a famat f
with t_INT
entries, returns a famat g
with all factors larger than limit
multiplied out as the last entry (with exponent 1
).
GEN
famat_to_nf(GEN nf, GEN f)
You normally never want to do this! This is a simplified form of nffactorback
, where we do not check the user input for consistency.
The description of famat_to_nf
says that you do not want to use this function. Then how do we recover genuine number field elements? Well, in most cases, we do not need to: most of the functions useful in this context accept famats as inputs, for instance nfsign
, nfsign_arch
, ideallog
and bnfisunit
. Otherwise, we can generally make good use of a quotient operation (modulo a fixed conductor, modulo ell-th powers); see the end of "Label se:Ideal_reduction".
@3Caveat. Receiving a famat input, bnfisunit
assumes that it is an algebraic integer, since this is expensive to check, and normally easy to ensure from the user's side; do not feed it ridiculous inputs.
GEN
famatsmall_reduce(GEN f)
as famat_reduce
, but for exponents given by a t_VECSMALL
.
@3Conversion to HNF.
GEN
idealhnf(GEN nf, GEN x)
returns the HNF of the ideal defined by x
: x
may be an algebraic number (defining a principal ideal), a maximal ideal (as given by idealprimedec
or idealfactor
), or a matrix whose columns give generators for the ideal. This last format is complicated, but useful to reduce general modules to the canonical form once in a while:
@3* if strictly less than N = [K:Q]
generators are given, x
is the Z_K
-module they generate,
@3* if N
or more are given, it is assumed that they form a Z-basis (that the matrix has maximal rank N
). This acts as mathnf
since the Z_K
-module structure is (taken for granted hence) not taken into account in this case.
Extended ideals are also accepted, their principal part being discarded.
GEN
idealhnf0(GEN nf, GEN x, GEN y)
returns the HNF of the ideal generated by the two algebraic numbers x
and y
.
The following low-level functions underlie the above two: they all assume that nf
is a true nf and perform no type checks:
GEN
idealhnf_principal(GEN nf, GEN x)
returns the ideal generated by the algebraic number x
.
GEN
idealhnf_shallow(GEN nf, GEN x)
is idealhnf
except that the result may not be suitable for gerepile
: if x
is already in HNF, we return x
, not a copy!
GEN
idealhnf_two(GEN nf, GEN v)
assuming a = v[1]
is a non-zero t_INT
and b = v[2]
is an algebraic integer, possibly given in regular representation by a t_MAT
(the multiplication table by b
, see zk_multable
), returns the HNF of a
Z_K+b
Z_K
.
@3Operations.
The basic ideal routines accept all nf
s (nf, bnf, bnr) and ideals in any form, including extended ideals, and return ideals in HNF, or an extended ideal when that makes sense:
GEN
idealadd(GEN nf, GEN x, GEN y)
returns x+y
.
GEN
idealdiv(GEN nf, GEN x, GEN y)
returns x/y
. Returns an extended ideal if x
or y
is an extended ideal.
GEN
idealmul(GEN nf, GEN x, GEN y)
returns xy
. Returns an extended ideal if x
or y
is an extended ideal.
GEN
idealsqr(GEN nf, GEN x)
returns x^2
. Returns an extended ideal if x
is an extended ideal.
GEN
idealinv(GEN nf, GEN x)
returns x^{-1}
. Returns an extended ideal if x
is an extended ideal.
GEN
idealpow(GEN nf, GEN x, GEN n)
returns x^n
. Returns an extended ideal if x
is an extended ideal.
GEN
idealpows(GEN nf, GEN ideal, long n)
returns x^n
. Returns an extended ideal if x
is an extended ideal.
GEN
idealmulred(GEN nf, GEN x, GEN y)
returns an extended ideal equal to xy
.
GEN
idealpowred(GEN nf, GEN x, GEN n)
returns an extended ideal equal to x^n
.
More specialized routines suffer from various restrictions:
GEN
idealdivexact(GEN nf, GEN x, GEN y)
returns x/y
, assuming that the quotient is an integral ideal. Much faster than idealdiv
when the norm of the quotient is small compared to Nx
. Strips the principal parts if either x
or y
is an extended ideal.
GEN
idealdivpowprime(GEN nf, GEN x, GEN pr, GEN n)
returns x
p^{-n}
, assuming x
is an ideal in HNF or a rational number, and pr
a prid attached to p. Not suitable for gerepileupto
since it returns x
when n = 0
.
GEN
idealmulpowprime(GEN nf, GEN x, GEN pr, GEN n)
returns x
p^{n}
, assuming x
is an ideal in HNF or a rational number, and pr
a prid attached to p. Not suitable for gerepileupto
since it retunrs x
when n = 0
.
GEN
idealprodprime(GEN nf, GEN v)
given a list v
of prime ideals in prid form, return their product. Assume that nf is a true nf structure.
GEN
idealprod(GEN nf, GEN v)
given a list v
of ideals, return their product. Assume that nf is a true nf structure.
GEN
idealHNF_mul(GEN nf, GEN x, GEN y)
returns xy
, assuming than nf
is a true nf, x
is an integral ideal in HNF and y
is an integral ideal in HNF or precompiled form (see below). For maximal speed, the second ideal y
may be given in precompiled form y = [a,b]
, where a
is a non-zero t_INT
and b
is an algebraic integer in regular representation (a t_MAT
giving the multiplication table by the fixed element): very useful when many ideals x
are going to be multiplied by the same ideal y
. This essentially reduces each ideal multiplication to an N x N
matrix multiplication followed by a N x 2N
modular HNF reduction (modulo xy
cap Z).
GEN
idealHNF_inv(GEN nf, GEN I)
returns I^{-1}
, assuming that nf
is a true nf and x
is a fractional ideal in HNF.
GEN
idealHNF_inv_Z(GEN nf, GEN I)
returns (I
cap Z).I^{-1}
, assuming that nf
is a true nf and x
is an integral fractional ideal in HNF. The result is an integral ideal in HNF.
@3Approximation.
GEN
idealaddtoone(GEN nf, GEN A, GEN B)
given to coprime integer ideals A
, B
, returns [a,b]
with a\in A
, b\in B
, such that a + b = 1
The result is reduced mod AB
, so a
, b
will be small.
GEN
idealaddtoone_i(GEN nf, GEN A, GEN B)
as idealaddtoone
except that nf
must be a true nf, and only a
is returned.
GEN
zkchineseinit(GEN nf, GEN A, GEN B, GEN AB)
given two coprime integral ideals A
and B
(in any form, preferably HNF) and their product AB
(in HNF form), initialize a solution to the Chinese remainder problem modulo AB
.
GEN
zkchinese(GEN zkc, GEN x, GEN y)
given zkc
from zkchineseinit
, and x
, y
two integral elements given as t_INT
or ZC
, return a z
modulo AB
such that z = x mod A
and z = y mod B
.
GEN
zkchinese1(GEN zkc, GEN x)
as zkchinese
for y = 1
; useful to lift elements in a nice way from (
Z_K/A_i)^*
to (
Z_K/
prod_i A_i)^*
.
GEN
hnfmerge_get_1(GEN A, GEN B)
given two square upper HNF integral matrices A
, B
of the same dimension n > 0
, return a
in the image of A
such that 1-a
is in the image of B
. (By abuse of notation we denote 1
the column vector [1,0,...,0]
.) If such an a
does not exist, return NULL
. This is the function underlying idealaddtoone
.
GEN
idealaddmultoone(GEN nf, GEN v)
given a list of n
(globally) coprime integer ideals (v[i])
returns an n
-dimensional vector a
such that a[i]\in v[i]
and sum a[i] = 1
. If [K:
Q] = N
, this routine computes the HNF reduction (with Gl_{nN}(
Z)
base change) of an N x nN
matrix; so it is well worth pruning "useless" ideals from the list (as long as the ideals remain globally coprime).
GEN
idealapprfact(GEN nf, GEN fx)
as idealappr
, except that x
must be given in factored form. (This is unchecked.)
GEN
idealcoprime(GEN nf, GEN x, GEN y)
. Given 2 integral ideals x
and y
, returns an algebraic number alpha such that alpha x
is an integral ideal coprime to y
.
GEN
idealcoprimefact(GEN nf, GEN x, GEN fy)
same as idealcoprime
, except that y
is given in factored form, as from idealfactor
.
GEN
idealchinese(GEN nf, GEN x, GEN y)
GEN
idealchineseinit(GEN nf, GEN x)
The PARI structure attached to maximal ideals is a prid (for prime ideal), usually produced by idealprimedec
and idealfactor
. In this section, we describe the format; other sections will deal with their daily use.
A prid attached to a maximal ideal p stores the following data: the underlying rational prime p
, the ramification degree e >= 1
, the residue field degree f >= 1
, a p
-uniformizer pi with valuation 1
at p and valuation 0
at all other primes dividing p
and a rescaled ``anti-uniformizer'' tau used to compute valuations. This tau is an algebraic integer such that tau/p
has valuation -1
at p and is integral at all other primes; in particular, the valuation of x\in
Z_K
is positive if and only if the algebraic integer x
tau is divisible by p
(easy to check for elements in t_COL
form).
GEN
pr_get_p(GEN pr)
returns p
. Shallow function.
GEN
pr_get_gen(GEN pr)
returns pi. Shallow function.
long
pr_get_e(GEN pr)
returns e
.
long
pr_get_f(GEN pr)
returns f
.
GEN
pr_get_tau(GEN pr)
returns zk_scalar_or_multable(
nf,
tau)
, which is the t_INT
1
iff p
is inert, and a ZM
otherwise. Shallow function.
int
pr_is_inert(GEN pr)
returns 1
if p
is inert, 0
otherwise.
GEN
pr_norm(GEN pr)
returns the norm p^f
of the maximal ideal.
ulong
upr_norm(GEN pr)
returns the norm p^f
of the maximal ideal, as an ulong
. Assume that the result does not overflow.
GEN
pr_inv(GEN pr)
return the fractional ideal p^{-1}
, in HNF.
GEN
pr_inv_p(GEN pr)
return the integral ideal p
p^{-1}
, in HNF.
GEN
idealprimedec(GEN nf, GEN p)
list of maximal ideals dividing the prime p
.
GEN
idealprimedec_limit_f(GEN nf, GEN p, long f)
as idealprimedec
, limiting the list to primes of residual degree <= f
if f
is non-zero.
GEN
idealprimedec_limit_norm(GEN nf, GEN p, GEN B)
as idealprimedec
, limiting the list to primes of norm <= B
, which must be a positive t_INT
.
GEN
idealprimedec_kummer(GEN nf, GEN Ti, long ei, GEN p)
let nf (true nf) correspond to K =
Q[X]/(T)
(T
monic ZX
). Let T =
prod_i T_i^{e_i} (mod p)
be the factorization of T
and let (f,g,h)
be as in Dedekind criterion for prime p
: f =
prod T_i
, g =
prod T_i^{e_i-1}
, h = (T - fg) / p
, and let D
be the gcd of (f,g,h)
in F_p[X]
. Let Ti
(FpX
) be one irreducible factor T_i
not dividing D
, with ei
= e_i
. This function returns the prime ideal attached to T_i
by Kummer / Dedekind criterion, namely p
Z_K + T_i(\bar{X})
Z_K
, which has ramification index e_i
over p
. Shallow function.
GEN
idealHNF_Z_factor(GEN x, GEN *pvN, GEN *pvZ)
given an integral (non-0
) ideal x
in HNF, compute both the factorization of Nx
and of x
capZ. This returns the vector of prime divisors of both and sets *pvN
and *pvZ
to the corresponding t_VECSMALL
vector of exponents for the factorization for the Norm and intersection with Z respetively.
GEN
nf_pV_to_prV(GEN)
{nf, GEN P} given a vector of rational primes P
, return the vector of all prime ideals above the P[i]
.
GEN
nf_deg1_prime(GEN nf)
let nf
be a true nf. This function returns a degree 1
(unramified) prime ideal not dividing nf.index
. In fact it returns an ideal above the smallest prime p >= [K:
Q]
satisfying those conditions.
GEN
prV_lcm_capZ(GEN L)
given a vector L
of prid (maximal ideals) return the squarefree positive integer generating their lcm intersected with Z. Not gerepile
-safe.
GEN
pr_uniformizer(GEN pr, GEN F)
given a prid attached to p / p
and F
in Z divisible exactly by p
, return an F
-uniformizer for pr
, i.e. a t
in Z_K
such that v_{
p}(t) = 1
and (t, F/
p) = 1
. Not gerepile
-safe.
GEN
idealramfrobenius(GEN nf, GEN gal, GEN pr, GEN ram)
Let K
be the number field defined by nf
and assume K/
Q be a Galois extension with Galois group given gal = galoisinit(nf)
, and that pr
is the prime ideal P in prid format, and that P is ramified, and ram
is its list of ramification groups as output by idealramgroups
. This function returns a permutation of gal.group
which defines an automorphism sigma in the decomposition group of P such that if p
is the unique prime number in P, then sigma(x) = x^p mod
P for all x\in
Z_K
.
GEN
idealfrobenius_aut(GEN nf, GEN gal, GEN pr, GEN aut)
faster version of idealfrobenius(nf, gal, pr
where aut
must be equal to nfgaloispermtobasis(nf, gal)
.
GEN
nfmodprinit(GEN nf, GEN pr)
returns an abstract modpr
structure, attached to reduction modulo the maximal ideal pr
, in idealprimedec
format. From this data we can quickly project any pr
-integral number field element to the residue field.
GEN
modpr_get_pr(GEN x)
return the pr
component from a modpr
structure.
GEN
modpr_get_p(GEN x)
return the p
component from a modpr
structure (underlying rational prime).
GEN
modpr_get_T(GEN x)
return the T
component from a modpr
structure: either NULL
(prime of degree 1
) or an irreducible FpX
defining the residue field over F_p
.
In library mode, it is often easier to use directly
GEN
nf_to_Fq_init(GEN nf, GEN *ppr, GEN *pT, GEN *pp)
concrete version of nfmodprinit
: nf
and *ppr
are the inputs, the return value is a modpr
and *ppr
, *pT
and *pp
are set as side effects.
The input *ppr
is either a maximal ideal or already a modpr
(in which case it is replaced by the underlying maximal ideal). The residue field is realized as F_p[X]/(T)
for some monic T\in
F_p[X]
, and we set *pT
to T
and *pp
to p
. Set T = NULL
if the prime has degree 1
and the residue field is F_p
.
In short, this receives (or initializes) a modpr
structure, and extracts from it T
, p
and p.
GEN
nf_to_Fq(GEN nf, GEN x, GEN modpr)
returns an Fq
congruent to x
modulo the maximal ideal attached to modpr
. The output is canonical: all elements in a given residue class are represented by the same Fq
.
GEN
Fq_to_nf(GEN x, GEN modpr)
returns an nf
element lifting the residue field element x
, either a t_INT
or an algebraic integer in algtobasis
format.
GEN
modpr_genFq(GEN modpr)
Returns an nf
element whose image by nf_to_Fq
is X \pmod T
, if deg T > 1
, else 1
.
GEN
zkmodprinit(GEN nf, GEN pr)
as nfmodprinit
, but we assume we will only reduce algebraic integers, hence do not initialize data allowing to remove denominators. More precisely, we can in fact still handle an x
whose rational denominator is not 0
in the residue field (i.e. if the valuation of x
is non-negative at all primes dividing p
).
GEN
zk_to_Fq_init(GEN nf, GEN *pr, GEN *T, GEN *p)
as nf_to_Fq_init
, able to reduce only p
-integral elements.
GEN
zk_to_Fq(GEN x, GEN modpr)
as nf_to_Fq
, for a p
-integral x
.
GEN
nfM_to_FqM(GEN M, GEN nf,GEN modpr)
reduces a matrix of nf
elements to the residue field; returns an FqM
.
GEN
FqM_to_nfM(GEN M, GEN modpr)
lifts an FqM
to a matrix of nf
elements.
GEN
nfV_to_FqV(GEN A, GEN nf,GEN modpr)
reduces a vector of nf
elements to the residue field; returns an FqV
with the same type as A
(t_VEC
or t_COL
).
GEN
FqV_to_nfV(GEN A, GEN modpr)
lifts an FqV
to a vector of nf
elements (same type as A
).
GEN
nfX_to_FqX(GEN Q, GEN nf,GEN modpr)
reduces a polynomial with nf
coefficients to the residue field; returns an FqX
.
GEN
FqX_to_nfX(GEN Q, GEN modpr)
lifts an FqX
to a polynomial with coefficients in nf
.
The following function is technical and may avoid computing a true nfmodpr
:
GEN
pr_basis_perm(GEN nf, GEN pr)
given a true nf structure and a prime ideal pr
above p
, return as a t_VECSMALL
the f(
p/p)
indices i
such that the nf.zk
[i]
mod p form an F_p
-basis of the residue field.
long
nfval(GEN nf, GEN x, GEN P)
return v_P(x)
@3Unsafe functions. assume nf is a genuine nf
structure, that P
, Q
are prid
.
long
ZC_nfval(GEN nf, GEN x, GEN P)
returns v_P(x)
, assuming x
is a ZC
, representing a non-zero algebraic integer.
long
ZC_nfvalrem(GEN nf, GEN x, GEN pr, GEN *newx)
returns v = v_P(x)
, assuming x
is a ZC
, representing a non-zero algebraic integer, and sets *newx
to x
tau^v
which is an algebraic integer coprime to p
.
int
ZC_prdvd(GEN nf, GEN x, GEN P)
returns 1
is P
divides x
and 0
otherwise. Assumes that x
is a ZC
, representing an algebraic integer. Faster than computing v_P(x)
.
int
pr_equal(GEN nf, GEN P, GEN Q)
returns 1 is P
and Q
represent the same maximal ideal: they must lie above the same p
and share the same e,f
invariants, but the p
-uniformizer and tau element may differ. Returns 0
otherwise.
``Signs'' of the real embeddings of number field element are represented in additive notation, using the standard identification (
Z/2
Z, +) \to ({-1,1}, x )
, s:--->(-1)^s
.
With respect to a fixed nf
structure, a selection of real places (a divisor at infinity) is normally given as a t_VECSMALL
of indices of the roots nf.roots
of the defining polynomial for the number field. For compatibility reasons, in particular under GP, the (obsolete) vec01
form is also accepted: a t_VEC
with gen_0
or gen_1
entries.
The following internal functions go back and forth between the two representations for the Archimedean part of divisors (GP: 0/1
vectors, library: list of indices):
GEN
vec01_to_indices(GEN v)
given a t_VEC
v
with t_INT
entries return as a t_VECSMALL
the list of indices i
such that v[i] != 0
. (Typically used with 0,1
-vectors but not necessarily so.) If v
is already a t_VECSMALL
, return it: not suitable for gerepile
in this case.
GEN
vecsmall01_to_indices(GEN v)
as
vec01_to_indices(zv_to_ZV(v));
GEN
indices_to_vec01(GEN p, long n)
return the 0/1
vector of length n
with ones exactly at the positions p[1], p[2],...
GEN
nfembed(GEN nf, GEN x, long k)
returns a floating point approximation of the k
-th embedding of x
(attached to the k
-th complex root in nf.roots
).
GEN
nfsign(GEN nf,GEN x)
x
being a number field element and nf
any form of number field, return the 0-1
-vector giving the signs of the r_1
real embeddings of x
, as a t_VECSMALL
. Linear algebra functions like Flv_add_inplace
then allow keeping track of signs in series of multiplications.
If x
is a t_VEC
of number field elements, return the matrix whose columns are the signs of the x[i]
.
GEN
nfsign_arch(GEN nf,GEN x,GEN arch)
arch
being a list of distinct real places, either in vec01
(t_VEC
with gen_0
or gen_1
entries) or indices
(t_VECSMALL
) form (see vec01_to_indices
), returns the signs of x
at the corresponding places. This is the low-level function underlying nfsign
.
int
nfchecksigns(GEN nf, GEN x, GEN pl)
pl is a t_VECSMALL
with r_1
components, all of which are in {-1,0,1}
. Return 1
if sigma_i(x)
pl[i] >= 0
for all i
, and 0
otherwise.
GEN
nfsign_units(GEN bnf, GEN archp, int add_tu)
archp
being a divisor at infinity in indices
form (or NULL
for the divisor including all real places), return the signs at archp
of a system of fundamental units for the field, in the same order as bnf.tufu
if add_tu
is set; and in the same order as bnf.fu
otherwise.
GEN
nfsign_from_logarch(GEN L, GEN invpi, GEN archp)
given L
the vector of the log sigma(x)
, where sigma runs through the (real or complex) embeddings of some number field, invpi
being a floating point approximation to 1/
pi, and archp
being a divisor at infinity in indices
form, return the signs of x
at the corresponding places. This is the low-level function underlying nfsign_units
; the latter is actually a trivial wrapper bnf
structures include the log sigma(x)
for a system of fundamental units of the field.
GEN
set_sign_mod_divisor(GEN nf, GEN x, GEN y, GEN sarch)
let f = f_0f_ oo
be a divisor, let sarch
be the output of nfarchstar(nf, f0, finf)
, and let x,y
be two number field elements. Returns yt
with t
integral, t = 1 mod^* f_0
such that x
and ty
have the same signs at f_ oo
; if x = NULL
, make ty
totally positive at f_ oo
.
GEN
nfarchstar(GEN nf, GEN f0, GEN finf)
for a divisor f = f_0f_ oo
represented by the integral ideal f0
in HNF and the finf
in indices
form, returns (
Z_K/f_ oo )^*
in a form suitable for computations mod f
. More precisely, returns [c, g, M,
finf]
, where c = [2,..., 2]
gives the cyclic structure of that group (#f_ oo
copies of Z/2
Z), g
a minimal system of independent generators, which are furthermore congruent to 1
mod f_0
(no condition if f_0 =
Z_K
), and M
is the matrix of signs of the g[i]
at f_ oo
, which is square and invertible over F_2
.
GEN
idealprincipalunits(GEN nf, GEN pr, long e)
returns the multiplicative group (1 +
pr) / (1 +
pr^e)
as an abelian group. Faster than idealstar
when the norm of pr is large, since it avoids (useless) work in the multiplicative group of the residue field.
nf
structureA number field K =
Q[X]/(T)
is defined by a monic T\in
Z[X]
. The low-level function computing a maximal order is
void
nfmaxord(nfmaxord_t *S, GEN T0, long flag)
, where the polynomial T_0
is squarefree with integer coefficients. Let K
be the étale algebra Q[X]/(T_0)
and let T = ZX_Q_normalize(T_0)
, i.e. T = C T_0(X/L)
is monic and integral for some C,Q\in
Q.
The structure nfmaxord_t
is initialized by the call; it has the following fields:
GEN T0, T, dT, dK; /* T0, T, discriminants of T and K */
GEN unscale; /* the integer L */
GEN index; /* index of power basis in maximal order */
GEN dTP, dTE; /* factorization of |dT|, primes / exponents */
GEN dKP, dKE; /* factorization of |dK|, primes / exponents */
GEN basis; /* Z-basis for maximal order of Q[X]/(T) */
@3The exponent vectors are t_VECSMALL
. The primes in dTP
and dKP
are pseudoprimes, not proven primes. We recommend restricting to T = T_0
, i.e. either to pass the input polynomial through ZX_Q_normalize
before the call, or to forget about T_0
and go on with the polynomial T
; otherwise unscale != 1
, all data is expressed in terms of T != T_0
, and needs to be converted to T_0
. For instance to convert the basis to Q[X]/(T_0)
:
RgXV_unscale(S.basis, S.unscale)
Instead of passing T
(monic ZX
), one can use the format [T,
listP]
as in nfbasis
or nfinit
, which computes an order which is maximal at a set of primes, but need not be the maximal order.
The flag
is an or-ed combination of the binary flags, both of them deprecated:
nf_PARTIALFACT
: do not try to fully factor dT
and only look for primes less than primelimit
. In that case, the elements in dTP
and dKP
need not all be primes. But the resulting dK
, index
and basis
are correct provided there exists no prime p > primelimit
such that p^2
divides the field discriminant dK
. This flag is deprecated: the [T,
listP]
format is safer and more flexible.
nf_ROUND2
: use the ROUND2 algorithm instead of the default ROUND4. This flag is deprecated: this algorithm is consistently slower.
void
nfinit_basic(nfmaxord_t *S, GEN T0)
a wrapper around nfmaxord
(without the deprecated flag
) that also accepts number field structures (nf, bnf,...) for T0
.
GEN
nfmaxord_to_nf(nfmaxord_t *S, GEN ro, long prec)
convert an nfmaxord_t
to an nf structure at precision prec
, where ro
is NULL
. The argument ro
may also be set to a vector with r_1 + r_2
components containing the roots of S- > T
suitably ordered, i.e. first r_1
t_REAL
roots, then r_2
t_COMPLEX
representing the conguate pairs, but this is strongly discouraged: the format is error-prone, and it is hard to compute the roots to the right accuracy in order to achieve prec
accuracy for the nf
. This function uses the integer basis S- > basis
as is, without performing LLL-reduction. Unless the basis is already known to be reduced, use rather the following higher-level function:
GEN
nfinit_complete(nfmaxord_t *S, long flag, long prec)
convert an nfmaxord_t
to an nf structure at precision prec
. The flag
has the same meaning as in nfinitall
. If S- > basis
is known to be reduced, it will be faster to use nfmaxord_to_nf
.
GEN
indexpartial(GEN T, GEN dT)
T
a monic separable ZX
, dT
is either NULL
(no information) or a multiple of the discriminant of T
. Let K =
Q[X]/(T)
and Z_K
its maximal order. Returns a multiple of the exponent of the quotient group Z_K/(
Z[X]/(T))
. In other word, a denominator d
such that d x\in
Z[X]/(T)
for all x\in
Z_K
.
We compute with arbitrary ideal representatives (in any of the various formats seen above), and call
GEN
bnfisprincipal0(GEN bnf, GEN x, long flag)
. The bnf
structure already contains information about the class group in the form \oplus_{i = 1}^n (
Z/d_i
Z) g_i
for canonical integers d_i
(with d_n | ... | d_1
all > 1
) and essentially random generators g_i
, which are ideals in HNF. We normally do not need the value of the g_i
, only that they are fixed once and for all and that any (non-zero) fractional ideal x
can be expressed uniquely as x = (t)
prod_{i = 1}^n g_i^{e_i}
, where 0 <= e_i < d_i
, and (t)
is some principal ideal. Computing e
is straightforward, but t
may be very expensive to obtain explicitly. The routine returns (possibly partial) information about the pair [e,t]
, depending on flag
, which is an or-ed combination of the following symbolic flags:
@3* nf_GEN
tries to compute t
. Returns [e,t]
, with t
an empty vector if the computation failed. This flag is normally useless in non-trivial situations since the next two serve analogous purposes in more efficient ways.
@3* nf_GENMAT
tries to compute t
in factored form, which is much more efficient than nf_GEN
if the class group is moderately large; imagine a small ideal x = (t)g^{10000}
: the norm of t
has 10000
as many digits as the norm of g
; do we want to see it as a vector of huge meaningless integers? The idea is to compute e
first, which is easy, then compute (t)
as x
prod g_i^{-e_i}
using successive idealmulred
, where the ideal reduction extracts small principal ideals along the way, eventually raised to large powers because of the binary exponentiation technique; the point is to keep this principal part in factored unexpanded form. Returns [e,t]
, with t
an empty vector if the computation failed; this should be exceedingly rare, unless the initial accuracy to which bnf
was computed was ridiculously low (and then bnfinit
should not have succeeded either). Setting/unsetting nf_GEN
has no effect when this flag is set.
@3* nf_GEN_IF_PRINCIPAL
tries to compute t
only if the ideal is principal (e = 0
). Returns gen_0
if the ideal is not principal. Setting/unsetting nf_GEN
has no effect when this flag is set, but setting/unsetting nf_GENMAT
is possible.
@3* nf_FORCE
in the above, insist on computing t
, even if it requires recomputing a bnf
from scratch. This is a last resort, and normally the accuracy of a bnf
can be increased without trouble, but it may be that some algebraic information simply cannot be recovered from what we have: see bnfnewprec
. It should be very rare, though.
In simple cases where you do not care about t
, you may use
GEN
isprincipal(GEN bnf, GEN x)
, which is a shortcut for bnfisprincipal0(bnf, x, 0)
.
The following low-level functions are often more useful:
GEN
isprincipalfact(GEN bnf, GEN C, GEN L, GEN f, long flag)
is about the same as bnfisprincipal0
applied to C
prod L[i]^{f[i]}
, where the L[i]
are ideals, the f[i]
integers and C
is either an ideal or NULL
(omitted). Make sure to include nf_GENMAT
in flag
!
GEN
isprincipalfact_or_fail(GEN bnf, GEN C, GEN L, GEN f)
is for delicate cases, where we must be more clever than nf_FORCE
(it is used when trying to increase the accuracy of a bnf, for instance). If performs
isprincipalfact(bnf,C, L, f, nf_GENMAT);
but if it fails to compute t
, it just returns a t_INT
, which is the estimated precision (in words, as usual) that would have been sufficient to complete the computation. The point is that nf_FORCE
does exactly this internally, but goes on increasing the accuracy of the bnf
, then discarding it, which is a major inefficiency if you intend to compute lots of discrete logs and have selected a precision which is just too low. (It is sometimes not so bad since most of the really expensive data is cached in bnf
anyway, if all goes well.) With this function, the caller may decide to increase the accuracy using bnfnewprec
(and keep the resulting bnf
!), or avoid the computation altogether. In any case the decision can be taken at the place where it is most likely to be correct.
void
bnftestprimes(GEN bnf, GEN B)
is an ingredient to certify unconditionnally a bnf
computed assuming GRH, cf. bnfcertify
. Running this function successfully proves that the classes of all prime ideals of norm <= B
belong to the subgroup of the class group generated by the factorbase used to compute the bnf
(equal to the class group under GRH). If the condition is not true, then (GRH is false and) the function will run forever.
If it is known that primes of norm less than B
generate the class group (through variants of Minkowski's convex body or Zimmert's twin classes theorems), then the true class group is proven to be a quotient of bnf.clgp
.
T_2
quadratic formWe assume the nf is a true nf
structure, attached to a number field K
of degree n
and signature (r_1,r_2)
. We saw that
GEN
nf_get_M(GEN nf)
returns the (r_1+r_2) x n
matrix M
giving the embeddings of K
, so that if v
is an n
-th dimensional t_COL
representing the element sum_{i = 1}^n v[i] w_i
of K
, then RgM_RgC_mul(M,v)
represents the embeddings of v
. Its first r_1
components are real numbers (t_INT
, t_FRAC
or t_REAL
, usually the latter), and the last r_2
are complex numbers (usually of t_COMPLEX
, but not necessarily for embeddings of rational numbers).
GEN
embed_T2(GEN x, long r1)
assuming x
is the vector of floating point embeddings of some algebraic number v
, i.e.
x = RgM_RgC_mul(nf_get_M(nf), algtobasis(nf,v));
@3returns T_2(v)
. If the floating point embeddings themselves are not needed, but only the values of T_2
, it is more efficient to restrict to real arithmetic and use
gnorml2( RgM_RgC_mul(nf_get_G(nf), algtobasis(nf,v)));
GEN
embednorm_T2(GEN x, long r1)
analogous to embed_T2
, applied to the gnorm
of the floating point embeddings. Assuming that
x = gnorm( RgM_RgC_mul(nf_get_M(nf), algtobasis(nf,v)) );
@3returns T_2(v)
.
GEN
embed_roots(GEN z, long r1)
given a vector z
of r_1+r_2
complex embeddings of the algebraic number v
, return the r_1+2r_2
roots of its characteristic polynomial. Shallow function.
GEN
embed_disc(GEN z, long r1, long prec)
given a vector z
of r_1+r_2
complex embeddings of the algebraic number v
, return a floating point approximation of the discriminant of its characteristic polynomial as a t_REAL
of precision prec
.
GEN
embed_norm(GEN x, long r1)
given a vector z
of r_1+r_2
complex embeddings of the algebraic number v
, return (a floating point approximation of) the norm of v
.
In the following routines nf is a true nf
, attached to a number field K
of degree n
:
GEN
nf_get_Gtwist(GEN nf, GEN v)
assuming v
is a t_VECSMALL
with r_1+r_2
entries, let
|| x ||_v^2 =
sum_{i = 1}^{r_1+r_2} 2^{v_i}
varepsilon_i|
sigma_i(x)|^2,
where as usual the sigma_i
are the (real and) complex embeddings and varepsilon_i = 1
, resp. 2
, for a real, resp. complex place. This is a twisted variant of the T_2
quadratic form, the standard Euclidean form on K\otimes
R. In applications, only the relative size of the v_i
will matter.
Let G_v\in M_n(
R)
be a square matrix such that if x\in K
is represented by the column vector X
in terms of the fixed Z-basis of Z_K
in nf, then
||x||_v^2 = ^t (G_v X).G_v X.
(This is a kind of Cholesky decomposition.) This function returns a rescaled copy of G_v
, rounded to nearest integers, specifically RM_round_maxrank
(G_v)
. Suitable for gerepileupto
, but does not collect garbage. For convenience, also allow v =
NULL
(nf_get_roundG
) and v
a t_MAT
as output from the function itself: in both these cases, shallow function.
GEN
nf_get_Gtwist1(GEN nf, long i)
. Simple special case. Returns the twisted G
matrix attached to the vector v
whose entries are all 0
except the i
-th one, which is equal to 10
.
GEN
idealpseudomin(GEN x, GEN G)
. Let x
, G
be two ZM
s, such that the product Gx
is well-defined. This returns a ``small'' integral linear combinations of the columns of x
, given by the LLL-algorithm applied to the lattice G x
. Suitable for gerepileupto
, but does not collect garbage.
In applications, x
is an integral ideal, G
approximates a Cholesky form for the T_2
quadratic form as returned by nf_get_Gtwist
, and we return a small element a
in the lattice (x,T_2)
. This is used to implement idealred
.
GEN
idealpseudomin_nonscalar(GEN x, GEN G)
. As idealpseudomin
, but we insist of returning a non-scalar a
(ZV_isscalar
is false), if the dimension of x
is > 1
.
In the interpretation where x
defines an integral ideal on a fixed Z_K
basis whose first element is 1
, this means that a
is not rational.
GEN
idealpseudored(GEN x, GEN G)
. As idealpseudomin
but we return the full reduced Z-basis of x
as a t_MAT
instead of a single vector.
GEN
idealred_elt(GEN nf, GEN x)
shortcut for
idealpseudomin(x, nf_get_roundG(nf))
Given an ideal x
this means finding a ``simpler'' ideal in the same ideal class. The public GP function is of course available
GEN
idealred0(GEN nf, GEN x, GEN v)
finds an a\in K^*
such that (a) x
is integral of small norm and returns it, as an ideal in HNF. What ``small'' means depends on the parameter v
, see the GP description. More precisely, a
is returned by idealpseudomin
((x_
Z) x^(-1),G)
divided by x_
Z, where x_
Z = (x
cap Z)
and where G
is nf_get_Gtwist
(
nf, v)
for v != NULL
and nf_get_roundG
(
nf)
otherwise.
@3Usually one sets v = NULL
to obtain an element of small T_2
norm in x
:
GEN
idealred(GEN nf, GEN x)
is a shortcut for idealred0(nf,x,NULL)
.
The function idealred
remains complicated to use: in order not to lose information x
must be an extended ideal, otherwise the value of a
is lost. There is a subtlety here: the principal ideal (a)
is easy to recover, but a
itself is an instance of the principal ideal problem which is very difficult given only an nf (once a bnf structure is available, bnfisprincipal0
will recover it).
GEN
idealmoddivisor(GEN bnr, GEN x)
A proof-of-concept implementation, useless in practice. If bnr
is attached to some modulus f
, returns a ``small'' ideal in the same class as x
in the ray class group modulo f
. The reason why this is useless is that using extended ideals with principal part in a computation, there is a simple way to reduce them: simply reduce the generator of the principal part in (
Z_K/f)^*
.
GEN
famat_to_nf_moddivisor(GEN nf, GEN g, GEN e, GEN bid)
given a true nf attached to a number field K
, a bid structure attached to a modulus f
, and an algebraic number in factored form prod g[i]^{e[i]}
, such that (g[i],f) = 1
for all i
, returns a small element in Z_K
congruent to it mod f
. Note that if f
contains places at infinity, this includes sign conditions at the specified places.
A simpler case when the conductor has no place at infinity:
GEN
famat_to_nf_modideal_coprime(GEN nf, GEN g, GEN e, GEN f, GEN expo)
as above except that the ideal f
is now integral in HNF (no need for a full bid), and we pass the exponent of the group (
Z_K/f)^*
as expo
; any multiple will also do, at the expense of efficiency. Of course if a bid for f
is available, if is easy to extract f
and the exact value of expo
from it (the latter is the first elementary divisor in the group structure). A useful trick: if you set expo
to any positive integer, the result is correct up to expo
-th powers, hence exact if expo
is a multiple of the exponent; this is useful when trying to decide whether an element is a square in a residue field for instance! (take expo
= 2
).
GEN
nf_to_Fp_coprime(GEN nf, GEN x, GEN modpr)
a low-level simple variant of famat_to_nf_modideal_coprime
: nf is a true nf structure, modpr
is from zkmodprinit
attached to a prime of degree 1
above the prime number p
, and x
is either a number field element or a famat
factorization matrix. We finally assume that no component of x
has a denominator p
.
What to do when the g[i]
are not coprime to f
, but only prod g[i]^{e[i]}
is? Then the situation is more complicated, and we advise to solve it one prime divisor of f
at a time. Let v
the valuation attached to a maximal ideal pr
and assume v(f) = k > 0
:
GEN
famat_makecoprime(GEN nf, GEN g, GEN e, GEN pr, GEN prk, GEN expo)
returns an element in (
Z_K/pr^k)^*
congruent to the product prod g[i]^{e[i]}
, assumed to be globally coprime to f
. As above, expo
is any positive multiple of the exponent of (
Z_K/pr^k)^*
, for instance (Nv-1)p^{k-1}
, if p
is the underlying rational prime. You may use other values of expo
(see the useful trick in famat_to_nf_modideal_coprime
).
GEN
Idealstarprk(GEN nf, GEN pr, long k, long flag)
same as Idealstar
for I = pr^k
Under GP, a class-field theoretic description of a number field is given by a triple A
, B
, C
, where the defining set [A,B,C]
can have any of the following forms: [
bnr]
, [
bnr,
subgroup]
, [
bnf,
modulus]
, [
bnf,
modulus,
subgroup]
. You can still use directly all of (libpari
's routines implementing) GP's functions as described in Chapter 3, but they are often awkward in the context of libpari
programming. In particular, it does not make much sense to always input a triple A,B,C
because of the fringe [
bnf,
modulus,
subgroup]
. The first routine to call, is thus
GEN
Buchray(GEN bnf, GEN mod, long flag)
initializes a bnr structure from bnf
and modulus mod
. flag
is an or-ed combination of nf_GEN
(include generators) and nf_INIT
(if omitted, do not return a bnr, only the ray class group as an abelian group). In fact, a single value of flag
actually makes sense: nf_GEN | nf_INIT
to initialize a proper bnr: removing nf_GEN
saves very little time, but the corresponding crippled bnr structure will raise errors in most class field theoretic functions. Possibly also 0 to quickly compute the ray class group structure; bnrclassno
is faster if we only need the order of the ray class group.
Now we have a proper bnr encoding a bnf
and a modulus, we no longer need the [
bnf,
modulus]
and [
bnf,
modulus,
subgroup]
forms, which would internally call Buchray
anyway. Recall that a subgroup H
is given by a matrix in HNF, whose column express generators of H
on the fixed generators of the ray class group that stored in our bnr. You may also code the trivial subgroup by NULL
.
GEN
bnrconductor(GEN bnr, GEN H, long flag)
see the documentation of the GP function.
GEN
bnrconductor_i(GEN bnr, GEN H, long flag)
shallow variant of bnrconductor
. Useful when flag = 2
and the conductor is the bnr
modulus: avoids copying the bnr
(wasteful).
long
bnrisconductor(GEN bnr, GEN H)
returns 1 is the class field defined by the subgroup H
(of the ray class group mod f
coded in bnr
) has conductor f
. Returns 0 otherwise.
GEN
bnrchar_primitive(GEN bnr, GEN chi, GEN bnrc)
Given a normalized character chi = [d,c]
on bnr.clgp
(see char_normalize
) of conductor bnrc.mod
, compute the primitive character chic
on bnrc.clgp
equivalent to chi
, given as a normalized character [D,C]
: chic(bnrc.gen[i])
is zeta_D^{C[i]}
, where D
is minimal. It is easier to use bnrconductor_i(bnr,chi,2)
, but the latter recomputes bnrc
for each new character.
GEN
bnrdisc(GEN bnr, GEN H, long flag)
returns the discriminant and signature of the class field defined by bnr
and H
. See the description of the GP function for details. flag is an or-ed combination of the flags rnf_REL
(output relative data) and rnf_COND
(return 0 unless the modulus is the conductor).
GEN
bnrsurjection(GEN BNR, GEN bnr)
BNR
and bnr
defined over the same field K
, for moduli F
and f
with F | f
, returns the matrix of the canonical surjection Cl_K(F)\to Cl_K(f)
(giving the image of the fixed ray class group generators of BNR
in terms of the ones in bnr
). BNR
must include the ray class group generators.
GEN
ABC_to_bnr(GEN A, GEN B, GEN C, GEN *H, int addgen)
This is a quick conversion function designed to go from the too general (inefficient) A
, B
, C
form to the preferred bnr, H
form for class fields. Given A
, B
, C
as explained above (omitted entries coded by NULL
), return the attached bnr, and set H
to the attached subgroup. If addgen
is 1
, make sure that if the bnr needed to be computed, then it contains generators.
GEN
nfissquarefree(GEN nf, GEN P)
given P
a polynomial with coefficients in nf, return 1
is P
is squarefree, and 0
otherwise. If is allowed (though less efficient) to replace nf by a monic ZX
defining the field.
GEN
rnfequationall(GEN A, GEN B, long *pk, GEN *pLPRS)
A
is either an nf type (corresponding to a number field K
) or an irreducible ZX
defining a number field K
. B
is an irreducible polynomial in K[X]
. Returns an absolute equation C
(over Q) for the number field K[X]/(B)
. C
is the characteristic polynomial of b + k a
for some roots a
of A
and b
of B
, and k
is a small rational integer. Set *pk
to k
.
If pLPRS
is not NULL
set it to [h_0, h_1]
, h_i\in
Q[X]
, where h_0+h_1 Y
is the last non-constant polynomial in the pseudo-Euclidean remainder sequence attached to A(Y)
and B(X-kY)
, leading to C = Res_Y(A(Y), B(Y-kX))
. In particular a := -h_0/h_1
is a root of A
in Q[X]/(C)
, and X - ka
is a root of B
.
GEN
nf_rnfeq(GEN A, GEN B)
wrapper around rnfequationall
to allow mapping K\to L
(eltup
) and converting elements of L
between absolute and relative form (reltoabs
, abstorel
), without computing a full rnf structure, which is useful if the relative integral basis is not required. In fact, since A
may be a t_POL
or an nf, the integral basis of the base field is not needed either. The return value is the same as rnf_get_map
. Shallow function.
GEN
nf_rnfeqsimple(GEN nf, GEN relpol)
as nf_rnfeq
except some fields are omitted, so that only the abstorel
operation is supported. Shallow function.
GEN
eltabstorel(GEN rnfeq, GEN x)
rnfeq
is as given by rnf_get_map
(but in this case rnfeltabstorel
is more robust), nf_rnfeq
or nf_rnfeqsimple
, return x
as an element of L/K
, i.e. as a t_POLMOD
with t_POLMOD
coefficients. Shallow function.
GEN
eltabstorel_lift(GEN rnfeq, GEN x)
same as eltabstorel
, except that x
is returned in partially lifted form, i.e. as a t_POL
with t_POLMOD
coefficients.
GEN
eltreltoabs(GEN rnfeq, GEN x)
rnfeq
is as given by rnf_get_map
(but in this case rnfeltreltoabs
is more robust) or nf_rnfeq
, return x
in absolute form.
void
nf_nfzk(GEN nf, GEN rnfeq, GEN *zknf, GEN *czknf)
rnfeq
as given by nf_rnfeq
, nf
a true nf structure, set *zknf
and *czknf
to a suitable representation of nf.zk
allowing quick computation of the map K\to L
by the function nfeltup
, without computing a full rnf structure, which is useful if the relative integral basis is not required. The computed values are the same as in rnf_get_nfzk
. Shallow function.
GEN
nfeltup(GEN nf, GEN x, GEN zknf, GEN czknf)
zknf
and czknf
are initialized by nf_nfzk
or rnf_get_nfzk
(but in this case rnfeltup
is more robust); nf
is a true nf structure for K
, returns x \in K
as a (lifted) element of L
, in absolute form.
GEN
Rg_nffix(const char *f, GEN T, GEN c, int lift)
given a ZX
T
and a ``coefficient'' c
supposedly belonging to Q[y]/(T)
, check whether this is a the case and return a cleaned up version of c
. The string f
is the calling function name, used to report errors.
This means that c
must be one of t_INT
, t_FRAC
, t_POL
in the variable y
with rational coefficients, or t_POLMOD
modulo T
which lift to a rational t_POL
as above. The cleanup consists in the following improvements:
@3* t_POL
coefficients are reduced modulo T
.
@3* t_POL
and t_POLMOD
belonging to Q are converted to rationals, t_INT
or t_FRAC
.
@3* if lift
is non-zero, convert t_POLMOD
to t_POL
, and otherwise convert t_POL
to t_POLMOD
s modulo T
.
GEN
RgX_nffix(const char *f, GEN T, GEN P, int lift)
check whether P
is a polynomials with coefficients in the number field defined by the absolute equation T(y) = 0
, where T
is a ZX
and returns a cleaned up version of P
. This checks whether P
is indeed a t_POL
with variable compatible with coefficients in Q[y]/(T)
, i.e.
varncmp(varn(P), varn(T)) < 0
@3and applies Rg_nffix
to each coefficient.
GEN
RgV_nffix(const char *f, GEN T, GEN P, int lift)
as RgX_nffix
for a vector of coefficients.
GEN
polmod_nffix(const char *f, GEN rnf, GEN x, int lift)
given a t_POLMOD
x
supposedly defining an element of rnf, check this and perform Rg_nffix
cleanups.
GEN
polmod_nffix2(const char *f, GEN T, GEN P, GEN x, int lift)
as in polmod_nffix
, where the relative extension is explicitly defined as L = (
Q[y]/(T))[x]/(P)
, instead of by an rnf
structure.
long
numberofconjugates(GEN T, long pinit)
returns a quick multiple for the number of Q-automorphism of the (integral, monic) t_POL
T
, from modular factorizations, starting from prime pinit
(you can set it to 2
). This upper bounds often coincides with the actual number of conjugates. Of course, you should use nfgaloisconj
to be sure.
GEN
nfroots_if_split(GEN *pt, GEN T)
let *pt
point either to a number field structure or an irreducible ZX
, defining a number field K
. Given T
a monic squarefree polynomial with coefficients in Z_K
, return the list of roots of pol
in K
if the polynomial splits completely, and NULL
otherwise. In other words, this checks whether K[X]/(T)
is normal over K
(hence Galois since T
is separable by assumption).
In the case where *pT
is a ZX
, the function has to compute internally a conditional nf
attached to K
, whose nf.zk
may not define the maximal order Z_K
(see nfroots
); *pT
is then replaced by the conditional nf
to avoid losing that information.
GEN
nfrootsof1(GEN nf)
returns a two-component vector [w,z]
where w
is the number of roots of unity in the number field nf, and z
is a primitive w
-th root of unity.
GEN
nfcyclotomicunits(GEN nf, GEN zu)
where zu
is as output by nfrootsof1(nf)
, return the vector of the cyclotomic units in nf
expressed over the integral basis.
Still provided for backward compatibility, but should not be used in new programs. They will eventually disappear.
GEN
zidealstar(GEN nf, GEN x)
short for Idealstar(nf,x,nf_GEN)
GEN
zidealstarinit(GEN nf, GEN x)
short for Idealstar(nf,x,nf_INIT)
GEN
zidealstarinitgen(GEN nf, GEN x)
short for Idealstar(nf,x,nf_GEN|nf_INIT)
GEN
buchimag(GEN D, GEN c1, GEN c2, GEN gCO)
short for
Buchquad(D,gtodouble(c1),gtodouble(c2), /*ignored*/0)
GEN
buchreal(GEN D, GEN gsens, GEN c1, GEN c2, GEN RELSUP, long prec)
short for
Buchquad(D,gtodouble(c1),gtodouble(c2), prec)
The following use a naming scheme which is error-prone and not easily extensible; besides, they compute generators as per nf_GEN
and not nf_GENMAT
. Don't use them:
GEN
isprincipalforce(GEN bnf,GEN x)
GEN
isprincipalgen(GEN bnf, GEN x)
GEN
isprincipalgenforce(GEN bnf, GEN x)
GEN
isprincipalraygen(GEN bnr, GEN x)
, use bnrisprincipal
.
@3Variants on polred
: use polredbest
.
GEN
factoredpolred(GEN x, GEN fa)
GEN
factoredpolred2(GEN x, GEN fa)
GEN
smallpolred(GEN x)
GEN
smallpolred2(GEN x)
, use Polred
.
GEN
polred0(GEN x, long flag, GEN p)
GEN
polredabs(GEN x)
GEN
polredabs2(GEN x)
GEN
polredabsall(GEN x, long flun)
@3Superseded by bnrdisclist0
:
GEN
discrayabslist(GEN bnf,GEN listes)
GEN
discrayabslistarch(GEN bnf, GEN arch, long bound)
@3Superseded by idealappr
(flag is ignored)
GEN
idealappr0(GEN nf, GEN x, long flag)
This section describes the data structure output by the function galoisinit
. This will be called a gal
structure in the following.
gal
structureThe functions below expect a gal
structure and are shallow. See the documentation of galoisinit
for the meaning of the member functions.
GEN
gal_get_pol(GEN gal)
returns gal.pol
GEN
gal_get_p(GEN gal)
returns gal.p
GEN
gal_get_e(GEN gal)
returns the integer e
such that gal.mod == gal.p^e
.
GEN
gal_get_mod(GEN gal)
returns gal.mod
.
GEN
gal_get_roots(GEN gal)
returns gal.roots
.
GEN
gal_get_invvdm(GEN gal)
gal[4]
.
GEN
gal_get_den(GEN gal)
return gal[5]
.
GEN
gal_get_group(GEN gal)
returns gal.group
.
GEN
gal_get_gen(GEN gal)
returns gal.gen
.
GEN
gal_get_orders(GEN gal)
returns gal.orders
.
GEN
nfgaloispermtobasis(GEN nf, GEN gal)
return the images of the field generator by the automorphisms gal.orders
expressed on the integral basis nf.zk
.
GEN
nfgaloismatrix(GEN nf, GEN s)
returns the ZM
attached to the automorphism s
, seen as a linear operator expressend on the number field integer basis. This allows to use
M = nfgaloismatrix(nf, s);
sx = ZM_ZC_mul(M, x); /* or RgM_RgC_mul(M, x) if x is not integral */
instead of
sx = nfgaloisapply(nf, s, x);
for an algebraic integer x
.
void
check_quaddisc(GEN x, long *s, long *mod4, const char *f)
checks whether the GEN
x
is a quadratic discriminant (t_INT
, not a square, congruent to 0,1
modulo 4
), and raise an exception otherwise. Set *s
to the sign of x
and *mod4
to x
modulo 4
(0 or 1).
void
check_quaddisc_real(GEN x, long *mod4, const char *f)
as check_quaddisc
; check that signe(x)
is positive.
void
check_quaddisc_imag(GEN x, long *mod4, const char *f)
as check_quaddisc
; check that signe(x)
is negative.
t_QFI
, t_QFR
GEN
qfi(GEN x, GEN y, GEN z)
creates the t_QFI
(x,y,z)
.
GEN
qfr(GEN x, GEN y, GEN z, GEN d)
creates the t_QFR
(x,y,z)
with distance component d
.
GEN
qfr_1(GEN q)
given a t_QFR
q
, return the unit form q^0
.
GEN
qfi_1(GEN q)
given a t_QFI
q
, return the unit form q^0
.
int
qfb_equal1(GEN q)
returns 1 if the t_QFI
or t_QFR
q
is the unit form.
GEN
qficomp(GEN x, GEN y)
compose the two t_QFI
x
and y
, then reduce the result. This is the same as gmul(x,y)
.
GEN
qfrcomp(GEN x, GEN y)
compose the two t_QFR
x
and y
, then reduce the result. This is the same as gmul(x,y)
.
GEN
qfisqr(GEN x)
as qficomp(x,y)
.
GEN
qfrsqr(GEN x)
as qfrcomp(x,y)
.
@3Same as above, without reducing the result:
GEN
qficompraw(GEN x, GEN y)
GEN
qfrcompraw(GEN x, GEN y)
GEN
qfisqrraw(GEN x)
GEN
qfrsqrraw(GEN x)
GEN
qfbcompraw(GEN x, GEN y)
compose two t_QFI
s or two t_QFR
s, without reduce the result.
GEN
powgi(GEN x, GEN n)
computes x^n
(will work for many more types than t_QFI
and t_QFR
, of course). Reduce the result.
GEN
qfrpow(GEN x, GEN n)
computes x^n
for a t_QFR
x
, reducing along the way. If the distance component is initially 0
, leave it alone; otherwise update it.
GEN
qfbpowraw(GEN x, long n)
compute x^n
(pure composition, no reduction), for a t_QFI
or t_QFR
x
.
GEN
qfipowraw(GEN x, long n)
as qfbpowraw
, for a t_QFI
x
.
GEN
qfrpowraw(GEN x, long n)
as qfbpowraw
, for a t_QFR
x
.
GEN
qfi_order(GEN q, GEN o)
assuming that the t_QFI
q
has order dividing o
, compute its order in the class group. The order can be given in all formats allowed by generic discrete log functions, the preferred format being [ord, fa]
(t_INT
and its factorization).
GEN
qfi_log(GEN a, GEN g, GEN o)
given a t_QFI
a
and assuming that the t_QFI
g
has order o
, compute an integer k
such that a^k = g
. Return cgetg(1, t_VEC)
if there are no solutions. Uses a generic Pollig-Hellman algorithm, then either Shanks (small o
) or Pollard rho (large o
) method. The order can be given in all formats allowed by generic discrete log functions, the preferred format being [ord, fa]
(t_INT
and its factorization).
GEN
qfi_Shanks(GEN a, GEN g, long n)
given a t_QFI
a
and assuming that the t_QFI
g
has (small) order n
, compute an integer k
such that a^k = g
. Return cgetg(1, t_VEC)
if there are no solutions. Directly uses Shanks algorithm, which is inefficient when n
is composite.
The following functions underly qfbsolve
; p
denotes a prime number.
GEN
qfisolvep(GEN Q, GEN p)
solves Q(x,y) = p
over the integers, for a t_QFI
Q
. Return gen_0
if there are no solutions.
GEN
qfrsolvep(GEN Q, GEN p)
solves Q(x,y) = p
over the integers, for a t_QFR
Q
. Return gen_0
if there are no solutions.
long
cornacchia(GEN d, GEN p, GEN *px, GEN *py)
solves x^2+ dy^2 = p
over the integers, where d > 0
. Return 1
if there is a solution (and store it in *x
and *y
), 0
otherwise.
long
cornacchia2(GEN d, GEN p, GEN *px, GEN *py)
as cornacchia
, for the equation x^2 + dy^2 = 4p
.
GEN
primeform_u(GEN x, ulong p)
t_QFI
whose first coefficient is the prime p
.
GEN
primeform(GEN x, GEN p, long prec)
Unfortunately, t_QFR
s are very inefficient, and are only provided for backward compatibility.
@3* they do not contain needed quantities, which are thus constantly recomputed (the discriminant D
, sqrt {D}
and its integer part),
@3* the distance component is stored in logarithmic form, which involves computing one extra logarithm per operation. It is much more efficient to store its exponential, computed from ordinary multiplications and divisions (taking exponent overflow into account), and compute its logarithm at the very end.
Internally, we have two representations for real quadratic forms:
@3* qfr3
, a container [a,b,c]
with at least 3 entries: the three coefficients; the idea is to ignore the distance component.
@3* qfr5
, a container with at least 5 entries [a,b,c,e,d]
: the three coefficients a t_REAL
d
and a t_INT
e
coding the distance component 2^{Ne} d
, in exponential form, for some large fixed N
.
It is a feature that qfr3
and qfr5
have no specified length or type. It implies that a qfr5
or t_QFR
will do whenever a qfr3
is expected. Routines using these objects all require a global context, provided by a struct qfr_data *
:
struct qfr_data {
GEN D; /* discriminant, t_INT */
GEN sqrtD; /* sqrt(D), t_REAL */
GEN isqrtD; /* floor(sqrt(D)), t_INT */
};
void
qfr_data_init(GEN D, long prec, struct qfr_data *S)
given a discriminant D > 0
, initialize S
for computations at precision prec
( sqrt {D}
is computed to that initial accuracy).
@3All functions below are shallow, and not stack clean.
GEN
qfr3_comp(GEN x, GEN y, struct qfr_data *S)
compose two qfr3
, reducing the result.
GEN
qfr3_pow(GEN x, GEN n, struct qfr_data *S)
compute x^n
, reducing along the way.
GEN
qfr3_red(GEN x, struct qfr_data *S)
reduce x
.
GEN
qfr3_rho(GEN x, struct qfr_data *S)
perform one reduction step; qfr3_red
just performs reduction steps until we hit a reduced form.
GEN
qfr3_to_qfr(GEN x, GEN d)
recover an ordinary t_QFR
from the qfr3
x
, adding distance component d
.
Before we explain qfr5
, recall that it corresponds to an ideal, that reduction corresponds to multiplying by a principal ideal, and that the distance component is a clever way to keep track of these principal ideals. More precisely, reduction consists in a number of reduction steps, going from the form (a,b,c)
to rho(a,b,c) = (c, -b mod 2c, *)
; the distance component is multiplied by (a floating point approximation to) (b +
sqrt {D}) / (b -
sqrt {D})
.
GEN
qfr5_comp(GEN x, GEN y, struct qfr_data *S)
compose two qfr5
, reducing the result, and updating the distance component.
GEN
qfr5_pow(GEN x, GEN n, struct qfr_data *S)
compute x^n
, reducing along the way.
GEN
qfr5_red(GEN x, struct qfr_data *S)
reduce x
.
GEN
qfr5_rho(GEN x, struct qfr_data *S)
perform one reduction step.
GEN
qfr5_dist(GEN e, GEN d, long prec)
decode the distance component from exponential (qfr5
-specific) to logarithmic form (as in a t_QFR
).
GEN
qfr_to_qfr5(GEN x, long prec)
convert a t_QFR
to a qfr5
with initial trivial distance component ( = 1
).
GEN
qfr5_to_qfr(GEN x, GEN d)
, assume x
is a qfr5
and d
was the original distance component of some t_QFR
that we converted using qfr_to_qfr5
to perform efficiently a number of operations. Convert x
to a t_QFR
with the correct (logarithmic) distance component.
GEN
ZM_hnf(GEN x)
returns the upper triangular Hermite Normal Form of the ZM
x
(removing 0
columns), using the ZM_hnfall
algorithm. If you want the true HNF, use ZM_hnfall(x, NULL, 0)
.
GEN
ZM_hnfmod(GEN x, GEN d)
returns the HNF of the ZM
x
(removing 0
columns), assuming the t_INT
d
is a multiple of the determinant of x
. This is usually faster than ZM_hnf
(and uses less memory) if the dimension is large, > 50
say.
GEN
ZM_hnfmodid(GEN x, GEN d)
returns the HNF of the matrix (x | d Id)
(removing 0
columns), for a ZM
x
and a t_INT
d
.
GEN
ZM_hnfmodall(GEN x, GEN d, long flag)
low-level function underlying the ZM_hnfmod
variants. If flag
is 0
, calls ZM_hnfmod(x,d)
; flag
is an or-ed combination of:
@3* hnf_MODID
call ZM_hnfmodid
instead of ZM_hnfmod
,
@3* hnf_PART
return as soon as we obtain an upper triangular matrix, saving time. The pivots are non-negative and give the diagonal of the true HNF, but the entries to the right of the pivots need not be reduced, i.e. they may be large or negative.
@3* hnf_CENTER
returns the centered HNF, where the entries to the right of a pivot p
are centered residues in [-p/2, p/2[
, hence smallest possible in absolute value, but possibly negative.
GEN
ZM_hnfmodall_i(GEN x, GEN d, long flag)
as ZM_hnfmodall
without final garbage collection. Not gerepile
-safe.
GEN
ZM_hnfall(GEN x, GEN *U, long remove)
returns the upper triangular HNF H
of the ZM
x
; if U
is not NULL
, set if to the matrix U
such that x U = H
. If remove = 0
, H
is the true HNF, including 0
columns; if remove = 1
, delete the 0
columns from H
but do not update U
accordingly (so that the integer kernel may still be recovered): we no longer have x U = H
; if remove = 2
, remove 0
columns from H
and update U
so that x U = H
. The matrix U
is square and invertible unless remove = 2
.
This routine uses a naive algorithm which is potentially exponential in the dimension (due to coefficient explosion) but is fast in practice, although it may require lots of memory. The base change matrix U
may be very large, when the kernel is large.
GEN
ZM_hnfall_i(GEN x, GEN *U, long remove)
as ZM_hnfall
without final garbage collection. Not gerepile
-safe.
GEN
ZM_hnfperm(GEN A, GEN *ptU, GEN *ptperm)
returns the hnf H = P A U
of the matrix P A
, where P
is a suitable permutation matrix, and U\in Gl_n(
Z)
. P
is chosen so as to (heuristically) minimize the size of U
; in this respect it is less efficient than ZM_hnflll
but usually faster. Set *ptU
to U
and *pterm
to a t_VECSMALL
representing the row permutation attached to P = (
delta_{i,perm[i]}
. If ptU
is set to NULL
, U
is not computed, saving some time; although useless, setting ptperm
to NULL
is also allowed.
GEN
ZM_hnf_knapsack(GEN x)
given a ZM
x
, compute its HNF h
. Return h
if it has the knapsack property: every column contains only zeroes and ones and each row contains a single 1
; return NULL
otherwise. Not suitable for gerepile.
GEN
ZM_hnflll(GEN x, GEN *U, int remove)
returns the HNF H
of the ZM
x
; if U
is not NULL
, set if to the matrix U
such that x U = H
. The meaning of remove
is the same as in ZM_hnfall
.
This routine uses the LLL variant of Havas, Majewski and Mathews, which is polynomial time, but rather slow in practice because it uses an exact LLL over the integers instead of a floating point variant; it uses polynomial space but lots of memory is needed for large dimensions, say larger than 300. On the other hand, the base change matrix U
is essentially optimally small with respect to the L_2
norm.
GEN
ZM_hnfcenter(GEN M)
. Given a ZM
in HNF M
, update it in place so that non-diagonal entries belong to a system of centered residues. Not suitable for gerepile.
Some direct applications: the following routines apply to upper triangular integral matrices; in practice, these come from HNF algorithms.
GEN
hnf_divscale(GEN A, GEN B,GEN t)
A
an upper triangular ZM
, B
a ZM
, t
an integer, such that C := tA^{-1}B
is integral. Return C
.
GEN
hnf_invscale(GEN A, GEN t)
A
an upper triangular ZM
, t
an integer such that C := tA^{-1}
is integral. Return C
. Special case of hnf_divscale
when B
is the identity matrix.
GEN
hnf_solve(GEN A, GEN B)
A
a ZM
in upper HNF (not necessarily square), B
a ZM
or ZC
. Return A^{-1}B
if it is integral, and NULL
if it is not.
GEN
hnf_invimage(GEN A, GEN b)
A
a ZM
in upper HNF (not necessarily square), b
a ZC
. Return A^{-1}B
if it is integral, and NULL
if it is not.
int
hnfdivide(GEN A, GEN B)
A
and B
are two upper triangular ZM
. Return 1
if A^{-1} B
is integral, and 0
otherwise.
@3Smith Normal Form.
GEN
ZM_snf(GEN x)
returns the Smith Normal Form (vector of elementary divisors) of the ZM
x
.
GEN
ZM_snfall(GEN x, GEN *U, GEN *V)
returns ZM_snf(x)
and sets U
and V
to unimodular matrices such that U x V = D
(diagonal matrix of elementary divisors). Either (or both) U
or V
may be NULL
in which case the corresponding matrix is not computed.
GEN
ZV_snfall(GEN d, GEN *U, GEN *V)
here d
is a ZV
; same as ZM_snfall
applied to diagonal(d)
, but faster.
GEN
ZM_snfall_i(GEN x, GEN *U, GEN *V, int returnvec)
same as ZM_snfall
, except that, depending on the value of returnvec
, we either return a diagonal matrix (as in ZM_snfall
, returnvec
is 0) or a vector of elementary divisors (as in ZM_snf
, returnvec
is 1) .
void
ZM_snfclean(GEN d, GEN U, GEN V)
assuming d
, U
, V
come from d = ZM_snfall(x, &U, &V)
, where U
or V
may be NULL
, cleans up the output in place. This means that elementary divisors equal to 1 are deleted and U
, V
are updated. The output is not suitable for gerepileupto
.
void
ZV_snf_trunc(GEN D)
given a vector D
of elementary divisors (i.e. a ZV
such that d_i | d_{i+1}
), truncate it in place to leave out the trivial divisors (equal to 1
).
GEN
ZM_snf_group(GEN H, GEN *U, GEN *Uinv)
this function computes data to go back and forth between an abelian group (of finite type) given by generators and relations, and its canonical SNF form. Given an abstract abelian group with generators g = (g_1,...,g_n)
and a vector X = (x_i)\in
Z^n
, we write g X
for the group element sum_i x_i g_i
; analogously if M
is an n x r
integer matrix g M
is a vector containing r
group elements. The group neutral element is 0
; by abuse of notation, we still write 0
for a vector of group elements all equal to the neutral element. The input is a full relation matrix H
among the generators, i.e. a ZM
(not necessarily square) such that gX = 0
for some X\in
Z^n
if and only if X
is in the integer image of H
, so that the abelian group is isomorphic to Z^n/Im H
. The routine assumes that H
is in HNF; replace it by its HNF if it is not the case. (Of course this defines the same group.)
Let G
a minimal system of generators in SNF for our abstract group: if the d_i
are the elementary divisors (... | d_2 | d_1
), each G_i
has either infinite order (d_i = 0
) or order d_i > 1
. Let D
the matrix with diagonal (d_i)
, then
G D = 0, G = g U_{inv}, g = G U,
for some integer matrices U
and U_{inv}
. Note that these are not even square in general; even if square, there is no guarantee that these are unimodular: they are chosen to have minimal entries given the known relations in the group and only satisfy D | (U U_{inv} - Id)
and H | (U_{inv}U - Id)
.
The function returns the vector of elementary divisors (d_i)
; if U
is not NULL
, it is set to U
; if Uinv
is not NULL
it is set to U_{inv}
. The function is not memory clean.
GEN
ZV_snf_group(GEN d, GEN *newU, GEN *newUi)
, here d
is a ZV
; same as ZM_snf_group
applied to diagonal(d)
, but faster.
The following 3 routines underly the various matrixqz
variants. In all case the m x n
t_MAT
x
is assumed to have rational (t_INT
and t_FRAC
) coefficients
GEN
QM_ImQ_hnf(GEN x)
returns an HNF basis for Im_
Q x
cap Z^n
.
GEN
QM_ImZ_hnf(GEN x)
returns an HNF basis for Im_
Z x
cap Z^n
.
GEN
QM_minors_coprime(GEN x, GEN D)
, assumes m >= n
, and returns a matrix in M_{m,n}(
Z)
with the same Q-image as x
, such that the GCD of all n x n
minors is coprime to D
; if D
is NULL
, we want the GCD to be 1
.
The following routines are simple wrappers around the above ones and are normally useless in library mode:
GEN
hnf(GEN x)
checks whether x
is a ZM
, then calls ZM_hnf
. Normally useless in library mode.
GEN
hnfmod(GEN x, GEN d)
checks whether x
is a ZM
, then calls ZM_hnfmod
. Normally useless in library mode.
GEN
hnfmodid(GEN x,GEN d)
checks whether x
is a ZM
, then calls ZM_hnfmodid
. Normally useless in library mode.
GEN
hnfall(GEN x)
calls ZM_hnfall(x, &U, 1)
and returns [H, U]
. Normally useless in library mode.
GEN
hnflll(GEN x)
calls ZM_hnflll(x, &U, 1)
and returns [H, U]
. Normally useless in library mode.
GEN
hnfperm(GEN x)
calls ZM_hnfperm(x, &U, &P)
and returns [H, U, P]
. Normally useless in library mode.
GEN
smith(GEN x)
checks whether x
is a ZM
, then calls ZM_snf
. Normally useless in library mode.
GEN
smithall(GEN x)
checks whether x
is a ZM
, then calls ZM_snfall(x, &U, &V)
and returns [U,V,D]
. Normally useless in library mode.
@3Some related functions over K[X]
, K
a field:
GEN
gsmith(GEN A)
the input matrix must be square, returns the elementary divisors.
GEN
gsmithall(GEN A)
the input matrix must be square, returns the [U,V,D]
, D
diagonal, such that UAV = D
.
GEN
RgM_hnfall(GEN A, GEN *pB, long remove)
analogous to ZM_hnfall
.
GEN
smithclean(GEN z)
cleanup the output of smithall
or gsmithall
(delete elementary divisors equal to 1
, updating base change matrices).
The basic GP functions and their immediate variants are normally not very useful in library mode. We briefly list them here for completeness, see the documentation of qflll
and qflllgram
for details:
@3* GEN
qflll0(GEN x, long flag)
GEN
lll(GEN x)
flag = 0
GEN
lllint(GEN x)
flag = 1
GEN
lllkerim(GEN x)
flag = 4
GEN
lllkerimgen(GEN x)
flag = 5
GEN
lllgen(GEN x)
flag = 8
@3* GEN
qflllgram0(GEN x, long flag)
GEN
lllgram(GEN x)
flag = 0
GEN
lllgramint(GEN x)
flag = 1
GEN
lllgramkerim(GEN x)
flag = 4
GEN
lllgramkerimgen(GEN x)
flag = 5
GEN
lllgramgen(GEN x)
flag = 8
The basic workhorse underlying all integral and floating point LLLs is
GEN
ZM_lll(GEN x, double D, long flag)
, where x
is a ZM
; D \in ]1/4,1[
is the Lovász constant determining the frequency of swaps during the algorithm: a larger values means better guarantees for the basis (in principle smaller basis vectors) but longer running times (suggested value: D = 0.99
).
@3Important. This function does not collect garbage and its output is not suitable for either gerepile
or gerepileupto
. We expect the caller to do something simple with the output (e.g. matrix multiplication), then collect garbage immediately.
@3flag
is an or-ed combination of the following flags:
@3* LLL_GRAM
. If set, the input matrix x
is the Gram matrix ^t v v
of some lattice vectors v
.
@3* LLL_INPLACE
. If unset, we return the base change matrix U
, otherwise the transformed matrix x U
or ^t U x U
(LLL_GRAM
). Implies LLL_IM
(see below).
@3* LLL_KEEP_FIRST
. The first vector in the output basis is the same one as was originally input. Provided this is a shortest non-zero vector of the lattice, the output basis is still LLL-reduced. This is used to reduce maximal orders of number fields with respect to the T_2
quadratic form, to ensure that the first vector in the output basis corresponds to 1
(which is a shortest vector).
@3* LLL_COMPATIBLE
. This is a no-op on 64-bit kernels; on 32-bit kernels, restrict to 64-bit-compatible accuracies in the course of LLL algorithms. This is very likely to produce identical results on all kernels, but this is not guaranteed.
The last three flags are mutually exclusive, either 0 or a single one must be set:
@3* LLL_KER
If set, only return a kernel basis K
(not LLL-reduced).
@3* LLL_IM
If set, only return an LLL-reduced lattice basis T
. (This is implied by LLL_INPLACE
).
@3* LLL_ALL
If set, returns a 2-component vector [K, T]
corresponding to both kernel and image.
GEN
lllfp(GEN x, double D, long flag)
is a variant for matrices with inexact entries: x
is a matrix with real coefficients (types t_INT
, t_FRAC
and t_REAL
), D
and flag are as in ZM_lll
. The matrix is rescaled, rounded to nearest integers, then fed to ZM_lll
. The flag LLL_INPLACE
is still accepted but less useful (it returns an LLL-reduced basis attached to rounded input, instead of an exact base change matrix).
GEN
ZM_lll_norms(GEN x, double D, long flag, GEN *ptB)
slightly more general version of ZM_lll
, setting *ptB
to a vector containing the squared norms of the Gram-Schmidt vectors (b_i^*)
attached to the output basis (b_i)
, b_i^ *= b_i +
sum_{j < i}
mu_{i,j} b_j^*
.
GEN
lllintpartial_inplace(GEN x)
given a ZM
x
of maximal rank, returns a partially reduced basis (b_i)
for the space spanned by the columns of x
: |b_i +- b_j| >= |b_i|
for any two distinct basis vectors b_i
, b_j
. This is faster than the LLL algorithm, but produces much larger bases.
GEN
lllintpartial(GEN x)
as lllintpartial_inplace
, but returns the base change matrix U
from the canonical basis to the b_i
, i.e. x U
is the output of lllintpartial_inplace
.
GEN
RM_round_maxrank(GEN G)
given a matrix G
with real floating point entries and independent columns, let G_e
be the rescaled matrix 2^e G
rounded to nearest integers, for e >= 0
. Finds a small e
such that the rank of G_e
is equal to the rank of G
(its number of columns) and return G_e
. This is useful as a preconditioning step to speed up LLL reductions, see nf_get_Gtwist
. Suitable for gerepileupto
, but does not collect garbage.
GEN
ZC_hnfremdiv(GEN x, GEN y, GEN *Q)
assuming y
is an invertible ZM
in HNF and x
is a ZC
, returns the ZC
R
equal to x
mod y
(whose i
-th entry belongs to [-y_{i,i}/2, y_{i,i}/2[
). Stack clean unless x
is already reduced (in which case, returns x
itself, not a copy). If Q
is not NULL
, set it to the ZC
such that x = yQ + R
.
GEN
ZM_hnfdivrem(GEN x, GEN y, GEN *Q)
reduce each column of the ZM
x
using ZC_hnfremdiv
. If Q
is not NULL
, set it to the ZM
such that x = yQ + R
.
GEN
ZC_hnfrem(GEN x, GEN y)
alias for ZC_hnfremdiv(x,y,NULL)
.
GEN
ZM_hnfrem(GEN x, GEN y)
alias for ZM_hnfremdiv(x,y,NULL)
.
GEN
ZC_reducemodmatrix(GEN v, GEN y)
Let y
be a ZM, not necessarily square, which is assumed to be LLL-reduced (otherwise, very poor reduction is expected). Size-reduces the ZC v
modulo the Z-module Y
spanned by y
: if the columns of y
are denoted by (y_1,..., y_{n-1})
, we return y_n = v
modulo Y
, such that the Gram-Schmidt coefficients mu_{n,j}
are less than 1/2
in absolute value for all j < n
. In short, y_n
is almost orthogonal to Y
.
GEN
ZM_reducemodmatrix(GEN v, GEN y)
Let y
be as in ZC_reducemodmatrix
, and v
be a ZM. This returns a matrix v
which is congruent to v
modulo the Z-module spanned by y
, whose columns are size-reduced. This is faster than repeatedly calling ZC_reducemodmatrix
on the columns since most of the Gram-Schmidt coefficients can be reused.
GEN
ZC_reducemodlll(GEN v, GEN y)
Let y
be an arbitrary ZM, LLL-reduce it then call ZC_reducemodmatrix
.
GEN
ZM_reducemodlll(GEN v, GEN y)
Let y
be an arbitrary ZM, LLL-reduce it then call ZM_reducemodmatrix
.
Besides the above functions, which were specific to integral input, we also have:
GEN
reducemodinvertible(GEN x, GEN y)
y
is an invertible matrix and x
a t_COL
or t_MAT
of compatible dimension. Returns x - y\lfloor y^{-1}x \rceil
, which has small entries and differs from x
by an integral linear combination of the columns of y
. Suitable for gerepileupto
, but does not collect garbage.
GEN
closemodinvertible(GEN x, GEN y)
returns x - reducemodinvertible(x,y)
, i.e. an integral linear combination of the columns of y
, which is close to x
.
GEN
reducemodlll(GEN x,GEN y)
LLL-reduce the non-singular ZM
y
and call reducemodinvertible
to find a small representative of x
mod y
Z^n
. Suitable for gerepileupto
, but does not collect garbage.
A finite abelian group G
in GP format is given by its Smith Normal Form as a pair [h,d]
or triple [h,d,g]
. Here h
is the cardinality of G
, (d_i)
is the vector of elementary divisors, and (g_i)
is a vector of generators. In short, G = \oplus_{i <= n} (
Z/d_i
Z) g_i
, with d_n | ... | d_2 | d_1
and prod d_i = h
.
Let e(x) :=
exp (2i
pi x)
. For ease of exposition, we restrict to complex-valued characters, but everything applies to more general fields K
where e
denotes a morphism (
Q,+) \to (K^*, x )
such that e(a/b)
denotes a b
-th root of unity.
A character on the abelian group \oplus (
Z/d_j
Z) g_j
is given by a row vector chi = [a_1,...,a_n]
such that chi(
prod g_j^{n_j}) = e(
sum a_j n_j / d_j)
.
GEN
cyc_normalize(GEN d)
shallow function. Given a vector (d_i)_{i <= n}
of elementary divisors for a finite group (no d_i
vanish), returns the vector D = [1]
if n = 0
(trivial group) and [d_1, d_1/d_2,..., d_1/d_n]
otherwise. This will allow to define characters as chi(
prod g_j^{x_j}) = e(
sum_j x_j a_j D_j / D_1)
, see char_normalize
.
GEN
char_normalize(GEN chi, GEN ncyc)
shallow function. Given a character chi
= (a_j)
and ncyc from cyc_normalize
above, returns the normalized representation [d, (n_j)]
, such that chi(
prod g_j^{x_j}) =
zeta_d^{
sum_j n_j x_j}
, where zeta_d = e(1/d)
and d
is minimal. In particular, d
is the order of chi
. Shallow function.
GEN
char_simplify(GEN D, GEN N)
given a quasi-normalized character [D, (N_j)]
such that chi(
prod g_j^{x_j}) =
zeta_D^{
sum_j N_j x_j}
, but where we only assume that D
is a multiple of the character order, return a normalized character [d, (n_j)]
with d
minimal. Shallow function.
GEN
char_denormalize(GEN cyc, GEN d, GEN n)
given a normalized representation [d, n]
(where d
need not be minimal) of a character on the abelian group with abelian divisors cyc
, return the attached character (where the image of each generator g_i
is given in terms of roots of unity of different orders cyc[i]
).
GEN
charconj(GEN cyc, GEN chi)
return the complex conjugate of chi
.
GEN
charmul(GEN cyc, GEN a, GEN b)
return the product character a x b
.
GEN
chardiv(GEN cyc, GEN a, GEN b)
returns the character a / b = a x \overline{b}
.
int
char_check(GEN cyc, GEN chi)
return 1
if chi
is a character compatible with cyclic factors cyc
, and 0
otherwise.
GEN
cyc2elts(GEN d)
given a t_VEC
d = (d_1,...,d_n)
of non-negative integers, return the vector of all t_VECSMALL
s of length n
whose i
-th entry lies in [0,d_i]
. Assumes that the product of the d_i
fits in a long
.
The functions in this section are specific to characters on (
Z/N
Z)^*
. The argument G
is a special bid
structure as returned by znstar0(N, nf_INIT)
. In this case, there are additional ways to input character via Conrey's representation. The character chi
is either a t_INT
(Conrey label), a t_COL
(a Conrey logarithm) or a t_VEC
(generic character on bid.gen
as explained in the previous subsection). The following low-level functions are called by GP's generic character functions.
int
zncharcheck(GEN G, GEN chi)
return 1
if chi
is a valid character and 0
otherwise.
GEN
zncharconj(GEN G, GEN chi)
as charconj
.
GEN
znchardiv(GEN G, GEN a, GEN b)
as chardiv
.
GEN
zncharker(GEN G, GEN chi)
as charker
.
GEN
znchareval(GEN G, GEN chi, GEN n, GEN z)
as chareval
.
GEN
zncharmul(GEN G, GEN a, GEN b)
as charmul
.
GEN
zncharorder(GEN G, GEN chi)
as charorder
.
The following functions handle characters in Conrey notation (attached to Conrey generators, not G.gen
):
int
znconrey_check(GEN cyc, GEN chi)
return 1
if chi
is a valid Conrey logarithm and 0
otherwise.
GEN
znconrey_normalized(GEN G, GEN chi)
return normalized character attached to chi
, as in char_normalize
but on Conrey generators.
GEN
znconreyfromchar(GEN G, GEN chi)
return Conrey logarithm attached to the generic (t_VEC
, on G.gen
)
GEN
znconreyfromchar_normalized(GEN G, GEN chi)
return normalized Conrey character attached to the generic (t_VEC
, on G.gen
) character chi
.
GEN
znconreylog_normalize(GEN G, GEN m)
given a Conrey logarithm m
(t_COL
), return the attached normalized Conrey character, as in char_normalize
but on Conrey generators.
GEN
Zideallog(GEN G, GEN x)
return the znconreylog
of x
expressed on G.gen
, i.e. the ordinary discrete logarithm from ideallog
.
Low-level routines underlying alginit
.
GEN
alg_csa_table(GEN nf, GEN mt, long v, long flag)
algebra defined by a multiplication table.
GEN
alg_cyclic(GEN rnf, GEN aut, GEN b, long flag)
cyclic algebra.
GEN
alg_hasse(GEN nf, long d, GEN hi, GEN hf, long v, long flag)
algebra defined by local Hasse invariants.
GEN
alg_hilbert(GEN nf, GEN a, GEN b, long v, long flag)
quaternion algebra.
GEN
alg_matrix(GEN nf, long n, long v, GEN L, long flag)
matrix algebra.
void
checkalg(GEN a)
raise an exception if a
was not initialized by alginit
.
long
alg_type(GEN al)
internal function called by algtype
: assume al
was created by alginit
(thereby saving a call to checkalg
). Return values are symbolic rather than numeric:
@3* al_NULL
: not a valid algebra.
@3* al_TABLE
: table algebra output by algtableinit
.
@3* al_CSA
: central simple algebra output by alginit
and represented by a multiplication table over its center.
@3* al_CYCLIC
: central simple algebra output by alginit
and represented by a cyclic algebra.
long
alg_model(GEN al, GEN x)
given an element x
in algebra al, check for inconsistencies (raise a type error) and return the representation model used for x
:
@3* al_ALGEBRAIC
: basistoalg
form, algebraic representation.
@3* al_BASIS
: algtobasis
form, column vector on the integral basis.
@3* al_MATRIX
: left multiplication table.
@3* al_TRIVIAL
: trivial algebra of degree 1
; can be understood as both basis or algebraic form (since e_1 = 1
).
All these routines assume their argument was initialized by alginit
and provide minor speedups compared to the GP equivalent. The routines returning a GEN
are shallow.
long
alg_get_absdim(GEN al)
low-level version of algabsdim
.
long
alg_get_dim(GEN al)
low-level version of algdim
.
long
alg_get_degree(GEN al)
low-level version of algdegree
.
GEN
alg_get_aut(GEN al)
low-level version of algaut
.
GEN
alg_get_auts(GEN al)
, given a cyclic algebra al = (L/K,
sigma,b)
of degree n
, returns the vector of sigma^i
, 1 <= i < n
.
GEN
alg_get_b(GEN al)
low-level version of algb
.
GEN
alg_get_basis(GEN al)
low-level version of albasis
.
GEN
alg_get_center(GEN al)
low-level version of algcenter
.
GEN
alg_get_char(GEN al)
low-level version of algchar
.
GEN
alg_get_hasse_f(GEN al)
low-level version of alghassef
.
GEN
alg_get_hasse_i(GEN al)
low-level version of alghassei
.
GEN
alg_get_invbasis(GEN al)
low-level version of alginvbasis
.
GEN
alg_get_multable(GEN al)
low-level version of algmultable
.
GEN
alg_get_relmultable(GEN al)
low-level version of algrelmultable
.
GEN
alg_get_splittingfield(GEN al)
low-level version of algsplittingfield
.
GEN
alg_get_abssplitting(GEN al)
returns the absolute nf structure attached to the rnf returned by algsplittingfield
.
GEN
alg_get_splitpol(GEN al)
returns the relative polynomial defining the rnf returned by algsplittingfield
.
GEN
alg_get_splittingdata(GEN al)
low-level version of algsplittingdata
.
GEN
alg_get_splittingbasis(GEN al)
the matrix Lbas from algsplittingdata
GEN
alg_get_splittingbasisinv(GEN al)
the matrix Lbasinv from algsplittingdata
.
GEN
alg_get_tracebasis(GEN al)
returns the traces of the basis elements; used by algtrace
.
GEN
alg_change_overorder_shallow(GEN al, GEN ord)
GEN
alg_complete(GEN rnf, GEN aut, GEN hi, GEN hf, int maxord)
GEN
alg_ordermodp(GEN al, GEN p)
GEN
algbasischarpoly(GEN al, GEN x, long v)
GEN
algbasismul(GEN al, GEN x, GEN y)
GEN
algbasismultable(GEN al, GEN x)
GEN
algbasismultable_Flm(GEN mt, GEN x, ulong m)
GEN
algleftordermodp(GEN al, GEN Ip, GEN p)
GEN
bnfgwgeneric(GEN bnf, GEN Lpr, GEN Ld, GEN pl, long var)
GEN
bnrgwsearch(GEN bnr, GEN Lpr, GEN Ld, GEN pl)
void
checkhasse(GEN nf, GEN hi, GEN hf, long n)
long
cyclicrelfrob(GEN rnf, GEN auts, GEN pr)
GEN
hassecoprime(GEN hi, GEN hf, long n)
GEN
hassedown(GEN nf, long n, GEN hi, GEN hf)
GEN
hassewedderburn(GEN hi, GEN hf, long n)
long
localhasse(GEN rnf, GEN cnd, GEN pl, GEN auts, GEN b, long k)
GEN
nfgwkummer(GEN nf, GEN Lpr, GEN Ld, GEN pl, long var)
\newpage