| Home | Trees | Indices | Help |
|---|
|
|
Field of Algebraic Numbers
AUTHOR:
-- Carl Witty (2007-01-27): initial version
-- Carl Witty (2007-10-29): massive rewrite to support complex
as well as real numbers
This is an implementation of the algebraic numbers (the complex
numbers which are the zero of a polynomial in ZZ[x]; in other words,
the algebraic closure of QQ, with an embedding into CC). All
computations are exact. We also include an implementation of the
algebraic reals (the intersection of the algebraic numbers with RR).
The field of algebraic numbers is available with abbreviation QQbar;
the field of algebraic reals has abbreviation AA.
As with many other implementations of the algebraic numbers, we try
hard to avoid computing a number field and working in the number
field; instead, we use floating-point interval arithmetic whenever
possible (basically whenever we need to prove non-equalities), and
resort to symbolic computation only as needed (basically to prove
equalities).
Algebraic numbers exist in one of the following forms:
* a rational number
* the product of a rational number and an n'th root of unity
* the sum, difference, product, or quotient of algebraic numbers
* the negation, inverse, absolute value, norm, real part,
imaginary part, or complex conjugate of an algebraic number
* a particular root of a polynomial, given as a polynomial with
algebraic coefficients, an isolating interval (given as a
RealIntervalFieldElement) which encloses exactly one root, and
the multiplicity of the root
* a polynomial in one generator, where the generator is an algebraic
number given as the root of an irreducible polynomial with integral
coefficients and the polynomial is given as a NumberFieldElement
The multiplicative subgroup of the algebraic numbers generated
by the rational numbers and the roots of unity is handled particularly
efficiently, as long as these roots of unity come from the QQbar.zeta()
method. Cyclotomic fields in general are fairly efficient, again
as long as they are derived from QQbar.zeta().
An algebraic number can be coerced into ComplexIntervalField (or
RealIntervalField, for algebraic reals); every algebraic number has a
cached interval of the highest precision yet calculated.
Everything is done with intervals except for comparisons. By default,
comparisons compute the two algebraic numbers with 128-bit precision
intervals; if this does not suffice to prove that the numbers are different,
then we fall back on exact computation.
Note that division involves an implicit comparison of the divisor against
zero, and may thus trigger exact computation.
Also, using an algebraic number in the leading coefficient of
a polynomial also involves an implicit comparison against zero, which
again may trigger exact computation.
Note that we work fairly hard to avoid computing new number fields;
to help, we keep a lattice of already-computed number fields and
their inclusions.
EXAMPLES:
sage: sqrt(AA(2)) > 0
True
sage: (sqrt(5 + 2*sqrt(QQbar(6))) - sqrt(QQbar(3)))^2 == 2
True
sage: AA((sqrt(5 + 2*sqrt(6)) - sqrt(3))^2) == 2
True
For a monic cubic polynomial x^3 + b*x^2 + c*x + d with roots
s1,s2,s3, the discriminant is defined as (s1-s2)^2(s1-s3)^2(s2-s3)^2
and can be computed as b^2c^2 - 4*b^3d - 4*c^3 + 18*bcd - 27*d^2.
We can test that these definitions do give the same result.
sage: def disc1(b, c, d):
... return b^2*c^2 - 4*b^3*d - 4*c^3 + 18*b*c*d - 27*d^2
sage: def disc2(s1, s2, s3):
... return ((s1-s2)*(s1-s3)*(s2-s3))^2
sage: x = polygen(AA)
sage: p = x*(x-2)*(x-4)
sage: cp = AA.common_polynomial(p)
sage: d, c, b, _ = p.list()
sage: s1 = AA.polynomial_root(cp, RIF(-1, 1))
sage: s2 = AA.polynomial_root(cp, RIF(1, 3))
sage: s3 = AA.polynomial_root(cp, RIF(3, 5))
sage: disc1(b, c, d) == disc2(s1, s2, s3)
True
sage: p = p + 1
sage: cp = AA.common_polynomial(p)
sage: d, c, b, _ = p.list()
sage: s1 = AA.polynomial_root(cp, RIF(-1, 1))
sage: s2 = AA.polynomial_root(cp, RIF(1, 3))
sage: s3 = AA.polynomial_root(cp, RIF(3, 5))
sage: disc1(b, c, d) == disc2(s1, s2, s3)
True
sage: p = (x-sqrt(AA(2)))*(x-AA(2).nth_root(3))*(x-sqrt(AA(3)))
sage: cp = AA.common_polynomial(p)
sage: d, c, b, _ = p.list()
sage: s1 = AA.polynomial_root(cp, RIF(1.4, 1.5))
sage: s2 = AA.polynomial_root(cp, RIF(1.7, 1.8))
sage: s3 = AA.polynomial_root(cp, RIF(1.2, 1.3))
sage: disc1(b, c, d) == disc2(s1, s2, s3)
True
We can coerce from symbolic expressions:
sage: QQbar(sqrt(-5))
[2.2360679774997893 .. 2.2360679774997899]*I
sage: AA(sqrt(2) + sqrt(3))
[3.1462643699419721 .. 3.1462643699419726]
sage: AA(sqrt(2)) + sqrt(3)
[3.1462643699419721 .. 3.1462643699419726]
sage: sqrt(2) + QQbar(sqrt(3))
[3.1462643699419721 .. 3.1462643699419726]
sage: QQbar(I)
1*I
sage: AA(I)
Traceback (most recent call last):
...
TypeError: Cannot coerce algebraic number with non-zero imaginary part to algebraic real
sage: QQbar(I * golden_ratio)
[1.6180339887498946 .. 1.6180339887498950]*I
sage: AA(golden_ratio)^2 - AA(golden_ratio)
1
sage: QQbar((-8)^(1/3))
[0.99999999999999988 .. 1.0000000000000003] + [1.7320508075688771 .. 1.7320508075688775]*I
sage: AA((-8)^(1/3))
-2
sage: QQbar((-4)^(1/4))
[1.0000000000000000 .. 1.0000000000000000] + [1.0000000000000000 .. 1.0000000000000000]*I
sage: AA((-4)^(1/4))
Traceback (most recent call last):
...
TypeError: Cannot coerce algebraic number with non-zero imaginary part to algebraic real
Note the different behavior in taking roots: for AA we prefer real roots if they exist, but for QQbar we take the principal root:
sage: AA(-1)^(1/3)
-1
sage: QQbar(-1)^(1/3)
[0.49999999999999994 .. 0.50000000000000012] + [0.86602540378443859 .. 0.86602540378443871]*I
We can explicitly coerce from QQ[I]. (Technically, this is not quite
kosher, since QQ[I] doesn't come with an embedding; we don't know
whether the field generator is supposed to map to +I or -I. We assume
that for any quadratic field with polynomial x^2+1, the generator maps
to +I.)
sage: K.<im> = QQ[I]
sage: pythag = QQbar(3/5 + 4*im/5); pythag
4/5*I + 3/5
sage: pythag.abs() == 1
True
However, implicit coercion from QQ[I] is not allowed.
sage: QQbar(1) + im
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for '+': 'Algebraic Field' and 'Number Field in I with defining polynomial x^2 + 1'
We can implicitly coerce from algebraic reals to algebraic numbers:
sage: a = QQbar(1); print a, a.parent()
1 Algebraic Field
sage: b = AA(1); print b, b.parent()
1 Algebraic Real Field
sage: c = a + b; print c, c.parent()
2 Algebraic Field
Some computation with radicals:
sage: phi = (1 + sqrt(AA(5))) / 2
sage: phi^2 == phi + 1
True
sage: tau = (1 - sqrt(AA(5))) / 2
sage: tau^2 == tau + 1
True
sage: phi + tau == 1
True
sage: tau < 0
True
sage: rt23 = sqrt(AA(2/3))
sage: rt35 = sqrt(AA(3/5))
sage: rt25 = sqrt(AA(2/5))
sage: rt23 * rt35 == rt25
True
The Sage rings AA and QQbar can decide equalities between radical
expressions (over the reals and complex numbers respectively):
sage: a = AA((2/(3*sqrt(3)) + 10/27)^(1/3) - 2/(9*(2/(3*sqrt(3)) + 10/27)^(1/3)) + 1/3)
sage: a
[0.99999999999999988 .. 1.0000000000000003]
sage: a == 1
True
Algebraic numbers which are known to be rational print as rationals;
otherwise they print as intervals (with 53-bit precision).
sage: AA(2)/3
2/3
sage: QQbar(5/7)
5/7
sage: QQbar(1/3 - 1/4*I)
-1/4*I + 1/3
sage: two = QQbar(4).nth_root(4)^2; two
[1.9999999999999997 .. 2.0000000000000005]
sage: two == 2; two
True
2
sage: phi
[1.6180339887498946 .. 1.6180339887498950]
We can find the real and imaginary parts of an algebraic number (exactly).
sage: r = QQbar.polynomial_root(x^5 - x - 1, CIF(RIF(0.1, 0.2), RIF(1.0, 1.1))); r
[0.18123244446987538 .. 0.18123244446987541] + [1.0839541013177105 .. 1.0839541013177108]*I
sage: r.real()
[0.18123244446987538 .. 0.18123244446987541]
sage: r.imag()
[1.0839541013177105 .. 1.0839541013177108]
sage: r.minpoly()
x^5 - x - 1
sage: r.real().minpoly()
x^10 + 3/16*x^6 + 11/32*x^5 - 1/64*x^2 + 1/128*x - 1/1024
sage: r.imag().minpoly() # this takes a long time (143s on my laptop)
x^20 - 5/8*x^16 - 95/256*x^12 - 625/1024*x^10 - 5/512*x^8 - 1875/8192*x^6 + 25/4096*x^4 - 625/32768*x^2 + 2869/1048576
We can find the absolute value and norm of an algebraic number exactly.
(Note that we define the norm as the product of a number and its
complex conjugate; this is the algebraic definition of norm, if we
view QQbar as AA[I].)
sage: r = QQbar((-8)^(1/3)); r
[0.99999999999999988 .. 1.0000000000000003] + [1.7320508075688771 .. 1.7320508075688775]*I
sage: r.abs() == 2
True
sage: r.norm() == 4
True
sage: (r+I).norm().minpoly()
x^2 - 10*x + 13
sage: r = AA.polynomial_root(x^2 - x - 1, RIF(-1, 0)); r
[-0.61803398874989491 .. -0.61803398874989479]
sage: r.abs().minpoly()
x^2 + x - 1
We can compute the multiplicative order of an algebraic number.
sage: QQbar(-1/2 + I*sqrt(3)/2).multiplicative_order()
3
sage: QQbar(-sqrt(3)/2 + I/2).multiplicative_order()
12
sage: QQbar.zeta(12345).multiplicative_order()
12345
Cyclotomic fields are very fast as long as we only multiply and divide:
sage: z3_3 = QQbar.zeta(3) * 3
sage: z4_4 = QQbar.zeta(4) * 4
sage: z5_5 = QQbar.zeta(5) * 5
sage: z6_6 = QQbar.zeta(6) * 6
sage: z20_20 = QQbar.zeta(20) * 20
sage: z3_3 * z4_4 * z5_5 * z6_6 * z20_20
7200
And they're still pretty fast even if you add and subtract, and trigger
exact computation:
sage: (z3_3 + z4_4 + z5_5 + z6_6 + z20_20)._exact_value()
4*zeta60^15 + 5*zeta60^12 + 9*zeta60^10 + 20*zeta60^3 - 3 where a^16 + a^14 - a^10 - a^8 - a^6 + a^2 + 1 = 0 and a in [0.99452189536827328 .. 0.99452189536827341] + [0.10452846326765345 .. 0.10452846326765352]*I
The paper _ARPREC: An Arbitrary Precision Computation Package_ discusses
this result. Evidently it is difficult to find, but we can easily
verify it.
sage: alpha = QQbar.polynomial_root(x^10 + x^9 - x^7 - x^6 - x^5 - x^4 - x^3 + x + 1, RIF(1, 1.2))
sage: lhs = alpha^630 - 1
sage: rhs_num = (alpha^315 - 1) * (alpha^210 - 1) * (alpha^126 - 1)^2 * (alpha^90 - 1) * (alpha^3 - 1)^3 * (alpha^2 - 1)^5 * (alpha - 1)^3
sage: rhs_den = (alpha^35 - 1) * (alpha^15 - 1)^2 * (alpha^14 - 1)^2 * (alpha^5 - 1)^6 * alpha^68
sage: rhs = rhs_num / rhs_den
sage: lhs
[2.6420403358193507e44 .. 2.6420403358193520e44]
sage: rhs
[2.6420403358193503e44 .. 2.6420403358193520e44]
sage: lhs - rhs
[-7.9344219392947342e28 .. 8.1800756658404269e28]
sage: lhs == rhs
True
sage: lhs - rhs
0
sage: lhs._exact_value()
-242494609856316402264822833062350847769474540*a^9 + 862295472068289472491654837785947906234680703*a^8 - 829559238431038252116584538075753012193290520*a^7 - 125882239615006638366472766103700441555126185*a^6 + 1399067970863104691667276008776398309383579345*a^5 - 1561176687069361567616835847286958553574223422*a^4 + 761706318888840943058230840550737823821027895*a^3 + 580740464974951394762758666210754821723780266*a^2 - 954587496403409756503464154898858512440951323*a + 546081123623099782018260884934770383777092602 where a^10 - 4*a^9 + 5*a^8 - a^7 - 6*a^6 + 9*a^5 - 6*a^4 - a^3 + 5*a^2 - 4*a + 1 = 0 and a in [0.44406334400909258 .. 0.44406334400909265]
We can pickle and unpickle algebraic fields (and they are globally unique):
sage: loads(dumps(AlgebraicField())) is AlgebraicField()
True
sage: loads(dumps(AlgebraicRealField())) is AlgebraicRealField()
True
We can pickle and unpickle algebraic numbers:
sage: loads(dumps(QQbar(10))) == QQbar(10)
True
sage: loads(dumps(QQbar(5/2))) == QQbar(5/2)
True
sage: loads(dumps(QQbar.zeta(5))) == QQbar.zeta(5)
True
sage: t = QQbar(sqrt(2)); type(t._descr)
<class 'sage.rings.qqbar.ANRoot'>
sage: loads(dumps(t)) == QQbar(sqrt(2))
True
sage: t.exactify(); type(t._descr)
<class 'sage.rings.qqbar.ANExtensionElement'>
sage: loads(dumps(t)) == QQbar(sqrt(2))
True
sage: t = ~QQbar(sqrt(2)); type(t._descr)
<class 'sage.rings.qqbar.ANUnaryExpr'>
sage: loads(dumps(t)) == 1/QQbar(sqrt(2))
True
sage: t = QQbar(sqrt(2)) + QQbar(sqrt(3)); type(t._descr)
<class 'sage.rings.qqbar.ANBinaryExpr'>
sage: loads(dumps(t)) == QQbar(sqrt(2)) + QQbar(sqrt(3))
True
|
|||
| _uniq_alg | |||
| _uniq_alg_r | |||
| AlgebraicField_common | |||
|
AlgebraicRealField The field of algebraic reals. |
|||
|
AlgebraicField The field of algebraic numbers. |
|||
|
AlgebraicGeneratorRelation A simple class for maintaining relations in the lattice of algebraic extensions. |
|||
|
AlgebraicGenerator An AlgebraicGenerator represents both an algebraic number alpha and the number field QQ[alpha]. |
|||
|
ANDescr An AlgebraicNumber or AlgebraicReal is a wrapper around an ANDescr object. |
|||
|
AlgebraicNumber_base This is the common base class for algebraic numbers (complex numbers which are the zero of a polynomial in ZZ[x]) and algebraic reals (algebraic numbers which happen to be real). |
|||
|
AlgebraicNumber The class for algebraic numbers (complex numbers which are the roots of a polynomial with integer coefficients). |
|||
| AlgebraicReal | |||
|
ANRational The subclass of ANDescr that represents an arbitrary rational. |
|||
|
ANRootOfUnity The subclass of ANDescr that represents a rational multiplied by a root of unity. |
|||
|
AlgebraicPolynomialTracker Keeps track of a polynomial used for algebraic numbers. |
|||
|
ANRoot The subclass of ANDescr that represents a particular root of a polynomial with algebraic coefficients. |
|||
|
ANExtensionElement The subclass of ANDescr that represents a number field element in terms of a specific generator. |
|||
| ANUnaryExpr | |||
| ANBinaryExpr | |||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
|
|||
CC = Complex Field with 53 bits of precision
|
|||
CIF = Complex Interval Field with 53 bits of precision
|
|||
_obj =
|
|||
_obj_r =
|
|||
is_SymbolicExpressionRing =
|
|||
AA =
|
|||
QQbar =
|
|||
_short_prec_seq =
|
|||
QQx = Univariate Polynomial Ring in x over Rational Field
|
|||
QQx_x = x
|
|||
QQy = Univariate Polynomial Ring in y over Rational Field
|
|||
QQy_y = y
|
|||
QQxy = Multivariate Polynomial Ring in x, y over Rational Field
|
|||
QQxy_x = x
|
|||
QQxy_y = y
|
|||
algebraic_generator_counter = 2
|
|||
_mul_algo =
|
|||
_add_algo =
|
|||
_descriptors =
|
|||
QQbarPoly = Univariate Polynomial Ring in x over Algebraic Field
|
|||
AAPoly = Univariate Polynomial Ring in x over Algebraic Real F
|
|||
qq_generator = Trivial generator
|
|||
_cyclotomic_gen_cache =
|
|||
ax = x
|
|||
RR_1_10 = 0.100000000000000
|
|||
QQ_0 = 0
|
|||
QQ_1 = 1
|
|||
QQ_1_2 = 1/2
|
|||
QQ_1_4 = 1/4
|
|||
AA_0 = 0
|
|||
QQbar_I_nf = Number Field in I with defining polynomial x^2 + 1
|
|||
QQbar_I_generator = Number Field in I with defining polynomial
|
|||
QQbar_I = 1*I
|
|||
AA_hash_offset = 1/123456789
|
|||
QQbar_hash_offset = 1/987654321*I + 1/123456789
|
|||
ZZX_x = x
|
|||
AA_golden_ratio =
|
|||
a =
|
|||
a1 =
|
|||
b =
|
|||
b1 =
|
|||
key =
|
|||
|
|||
Checks whether the rational r is an exact d'th power. If so, returns
the d'th root of r; otherwise, returns None.
EXAMPLES:
sage: from sage.rings.qqbar import rational_exact_root
sage: rational_exact_root(16/81, 4)
2/3
sage: rational_exact_root(8/81, 3) is None
True
|
Takes a monic polynomial and rescales the variable to get a monic
polynomial with "integral" coefficients. Works on any univariate
polynomial whose base ring has a denominator() method that returns
integers; for example, the base ring might be QQ or a number
field.
Returns the scale factor and the new polynomial.
(Inspired by Pari's primitive_pol_to_monic().)
We assume that coefficient denominators are "small"; the algorithm
factors the denominators, to give the smallest possible scale factor.
EXAMPLES:
sage: from sage.rings.qqbar import clear_denominators
sage: _.<x> = QQ['x']
sage: clear_denominators(x + 3/2)
(2, x + 3)
sage: clear_denominators(x^2 + x/2 + 1/4)
(2, x^2 + x + 1)
|
Find the polynomial of lowest discriminant that generates
the same field as poly, out of those returned by the Pari
polred routine. Returns a triple: (elt_fwd, elt_back, new_poly),
where new_poly is the new polynomial, elt_fwd is a polynomial
expression for a root of the new polynomial in terms of a root
of the original polynomial, and elt_back is a polynomial expression
for a root of the original polynomial in terms of a root of the
new polynomial.
EXAMPLES:
sage: from sage.rings.qqbar import do_polred
sage: _.<x> = QQ['x']
sage: do_polred(x^2-5)
(-1/2*x + 1/2, -2*x + 1, x^2 - x - 1)
sage: do_polred(x^2-x-11)
(-1/3*x + 2/3, -3*x + 2, x^2 - x - 1)
sage: do_polred(x^3 + 123456)
(-1/4*x, -4*x, x^3 - 1929)
|
intv_fn is a function that takes a precision and returns an
interval of that precision containing some particular root of pol.
(It must return better approximations as the precision increases.)
pol is an irreducible polynomial with rational coefficients.
Returns an interval containing at most one root of pol.
EXAMPLES:
sage: from sage.rings.qqbar import isolating_interval
sage: _.<x> = QQ['x']
sage: isolating_interval(lambda prec: sqrt(RealIntervalField(prec)(2)), x^2 - 2)
[1.41421356237309504876 .. 1.41421356237309504888]
And an example that requires more precision:
sage: delta = 10^(-70)
sage: p = (x - 1) * (x - 1 - delta) * (x - 1 + delta)
sage: isolating_interval(lambda prec: RealIntervalField(prec)(1 + delta), p)
[1.00000000000000000000000000000000000000000000000000000000000000000000009999999999999999999999999999999999999999999999999999999999999999999999999999999999998 .. 1.00000000000000000000000000000000000000000000000000000000000000000000010000000000000000000000000000000000000000000000000000000000000000000000000000000000014]
The function also works with complex intervals and complex roots:
sage: p = x^2 - x + 13/36
sage: isolating_interval(lambda prec: ComplexIntervalField(prec)(1/2, 1/3), p)
[0.500000000000000000000 .. 0.500000000000000000000] + [0.333333333333333333315 .. 0.333333333333333333343]*I
|
l is a list of some sort. fn is a function which maps an element
of l and a precision into an interval (either real or complex) of
that precision, such that for sufficient precision, exactly one
element of l results in an interval containing 0. Returns that
one element of l.
EXAMPLES:
sage: from sage.rings.qqbar import find_zero_result
sage: _.<x> = QQ['x']
sage: delta = 10^(-70)
sage: p1 = x - 1
sage: p2 = x - 1 - delta
sage: p3 = x - 1 + delta
sage: p2 == find_zero_result(lambda p, prec: p(RealIntervalField(prec)(1 + delta)), [p1, p2, p3])
True
|
If the interval v (which may be real or complex) includes some
purely real numbers, return v' containing v such that
v' == v'.conjugate(). Otherwise return v unchanged. (Note that if
v' == v'.conjugate(), and v' includes one non-real root of a real
polynomial, then v' also includes the conjugate of that root.
Also note that the diameter of the return value is at most twice
the diameter of the input.)
EXAMPLES:
sage: from sage.rings.qqbar import conjugate_expand
sage: conjugate_expand(CIF(RIF(0, 1), RIF(1, 2)))
[0.00000000000000000 .. 1.0000000000000000] + [1.0000000000000000 .. 2.0000000000000000]*I
sage: conjugate_expand(CIF(RIF(0, 1), RIF(0, 1)))
[0.00000000000000000 .. 1.0000000000000000] + [-1.0000000000000000 .. 1.0000000000000000]*I
sage: conjugate_expand(CIF(RIF(0, 1), RIF(-2, 1)))
[0.00000000000000000 .. 1.0000000000000000] + [-2.0000000000000000 .. 2.0000000000000000]*I
sage: conjugate_expand(RIF(1, 2))
[1.0000000000000000 .. 2.0000000000000000]
|
If the interval v includes some purely real numbers, return
a real interval containing only those real numbers. Otherwise
return v unchanged.
If v includes exactly one root of a real polynomial, and v was
returned by conjugate_expand(), then conjugate_shrink(v) still
includes that root, and is a RealIntervalFieldElement iff the root
in question is real.
EXAMPLES:
sage: from sage.rings.qqbar import conjugate_shrink
sage: conjugate_shrink(RIF(3, 4))
[3.0000000000000000 .. 4.0000000000000000]
sage: conjugate_shrink(CIF(RIF(1, 2), RIF(1, 2)))
[1.0000000000000000 .. 2.0000000000000000] + [1.0000000000000000 .. 2.0000000000000000]*I
sage: conjugate_shrink(CIF(RIF(1, 2), RIF(0, 1)))
[1.0000000000000000 .. 2.0000000000000000]
sage: conjugate_shrink(CIF(RIF(1, 2), RIF(-1, 2)))
[1.0000000000000000 .. 2.0000000000000000]
|
Given a sequence of elements of either \code{AA} or \code{QQbar}
(or a mixture), computes a number field containing all of these
elements, these elements as members of that number field, and a
homomorphism from the number field back to \code{AA} or
\code{QQbar}.
This may not return the smallest such number field, unless
\var{minimal}=\code{True} is specified.
Also, a single number can be passed, rather than a sequence; and
any values which are not elements of \code{AA} or \code{QQbar}
will automatically be coerced to \code{QQbar}.
This function may be useful for efficiency reasons: doing exact
computations in the corresponding number field will be faster
than doing exact computations directly in \code{AA} or \code{QQbar}.
EXAMPLES:
We can use this to compute the splitting field of a polynomial.
(Unfortunately this takes an unreasonably long time for non-toy
examples.)
sage: x = polygen(QQ)
sage: p = x^3 + x^2 + x + 17
sage: rts = p.roots(ring=QQbar, multiplicities=False)
sage: splitting = number_field_elements_from_algebraics(rts)[0]; splitting
Number Field in a with defining polynomial y^6 + 169*y^4 + 7968*y^2 + 121088
sage: p.roots(ring=splitting)
[(-9/2176*a^4 - 1121/2176*a^2 - 1625/136, 1), (9/17408*a^5 + 9/4352*a^4 + 1121/17408*a^3 + 1121/4352*a^2 + 1489/1088*a + 1489/272, 1), (-9/17408*a^5 + 9/4352*a^4 - 1121/17408*a^3 + 1121/4352*a^2 - 1489/1088*a + 1489/272, 1)]
sage: rt2 = AA(sqrt(2)); rt2
[1.4142135623730949 .. 1.4142135623730952]
sage: rt3 = AA(sqrt(3)); rt3
[1.7320508075688771 .. 1.7320508075688775]
sage: qqI = QQbar.zeta(4); qqI
1*I
sage: z3 = QQbar.zeta(3); z3
[-0.50000000000000012 .. -0.49999999999999994] + [0.86602540378443859 .. 0.86602540378443871]*I
sage: rt2b = rt3 + rt2 - rt3; rt2b
[1.4142135623730949 .. 1.4142135623730952]
sage: rt2c = z3 + rt2 - z3; rt2c
[1.4142135623730949 .. 1.4142135623730952] + [-2.7105054312137611e-19 .. 2.7105054312137611e-19]*I
sage: number_field_elements_from_algebraics(rt2)
(Number Field in a with defining polynomial y^2 - 2, a, Ring morphism:
From: Number Field in a with defining polynomial y^2 - 2
To: Algebraic Real Field
Defn: a |--> [1.4142135623730949 .. 1.4142135623730952])
sage: number_field_elements_from_algebraics((rt2,rt3))
(Number Field in a with defining polynomial y^4 - 4*y^2 + 1, [-a^3 + 3*a, -a^2 + 2], Ring morphism:
From: Number Field in a with defining polynomial y^4 - 4*y^2 + 1
To: Algebraic Real Field
Defn: a |--> [0.51763809020504147 .. 0.51763809020504159])
We've created \code{rt2b} in such a way that \sage doesn't initially know
that it's in a degree-2 extension of $\QQ$.
sage: number_field_elements_from_algebraics(rt2b)
(Number Field in a with defining polynomial y^4 - 4*y^2 + 1, -a^3 + 3*a, Ring morphism:
From: Number Field in a with defining polynomial y^4 - 4*y^2 + 1
To: Algebraic Real Field
Defn: a |--> [0.51763809020504147 .. 0.51763809020504159])
We can specify \code{minimal=True} if we want the smallest number field.
sage: number_field_elements_from_algebraics(rt2b, minimal=True)
(Number Field in a with defining polynomial y^2 - 2, a, Ring morphism:
From: Number Field in a with defining polynomial y^2 - 2
To: Algebraic Real Field
Defn: a |--> [1.4142135623730949 .. 1.4142135623730952])
Things work fine with rational numbers, too.
sage: number_field_elements_from_algebraics((QQbar(1/2), AA(17)))
(Rational Field, [1/2, 17], Ring morphism:
From: Rational Field
To: Algebraic Real Field
Defn: 1 |--> 1)
Or we can just pass in symbolic expressions, as long as they can be
coerced into \code{QQbar}.
sage: number_field_elements_from_algebraics((sqrt(7), sqrt(9), sqrt(11)))
(Number Field in a with defining polynomial y^4 - 9*y^2 + 1, [-a^3 + 8*a, 3, -a^3 + 10*a], Ring morphism:
From: Number Field in a with defining polynomial y^4 - 9*y^2 + 1
To: Algebraic Real Field
Defn: a |--> [0.33543673964540460 .. 0.33543673964540466])
Here we see an example of doing some computations with number field
elements, and then mapping them back into \code{QQbar}.
sage: (fld,nums,hom) = number_field_elements_from_algebraics((rt2, rt3, qqI, z3))
sage: fld,nums,hom
(Number Field in a with defining polynomial y^8 - y^4 + 1, [-a^5 + a^3 + a, a^6 - 2*a^2, a^6, -a^4], Ring morphism:
From: Number Field in a with defining polynomial y^8 - y^4 + 1
To: Algebraic Field
Defn: a |--> [-0.25881904510252080 .. -0.25881904510252073] - [0.96592582628906820 .. 0.96592582628906832]*I)
sage: (nfrt2, nfrt3, nfI, nfz3) = nums
sage: hom(nfrt2)
[1.4142135623730949 .. 1.4142135623730952] + [-5.4210108624275222e-19 .. 7.8604657505199072e-19]*I
sage: nfrt2^2
2
sage: nfrt3^2
3
sage: nfz3 + nfz3^2
-1
sage: nfI^2
-1
sage: sum = nfrt2 + nfrt3 + nfI + nfz3; sum
2*a^6 - a^5 - a^4 + a^3 - 2*a^2 + a
sage: hom(sum)
[2.6462643699419721 .. 2.6462643699419726] + [1.8660254037844385 .. 1.8660254037844389]*I
sage: hom(sum) == rt2 + rt3 + qqI + z3
True
sage: [hom(n) for n in nums] == [rt2, rt3, qqI, z3]
True
TESTS:
sage: number_field_elements_from_algebraics(rt3)
(Number Field in a with defining polynomial y^2 - 3, a, Ring morphism:
From: Number Field in a with defining polynomial y^2 - 3
To: Algebraic Real Field
Defn: a |--> [1.7320508075688771 .. 1.7320508075688775])
sage: number_field_elements_from_algebraics((rt2,qqI))
(Number Field in a with defining polynomial y^4 + 1, [a^3 - a, -a^2], Ring morphism:
From: Number Field in a with defining polynomial y^4 + 1
To: Algebraic Field
Defn: a |--> [-0.70710678118654758 .. -0.70710678118654746] + [0.70710678118654746 .. 0.70710678118654758]*I)
Note that for the first example, where \sage doesn't realize that
the number is real, we get a homomorphism to \code{QQbar}; but with
\code{minimal=True}, we get a homomorphism to \code{AA}. Also note
that the exact answer depends on a Pari function that gives
different answers for 32-bit and 64-bit machines.
sage: number_field_elements_from_algebraics(rt2c)
(Number Field in a with defining polynomial y^4 + 2*y^2 + 4, 1/2*a^3, Ring morphism:
From: Number Field in a with defining polynomial y^4 + 2*y^2 + 4
To: Algebraic Field
Defn: a |--> [-0.70710678118654758 .. -0.70710678118654746] + [1.2247448713915889 .. 1.2247448713915892]*I) # 32-bit
Defn: a |--> [-0.70710678118654758 .. -0.70710678118654746] - [1.2247448713915889 .. 1.2247448713915892]*I) # 64-bit
sage: number_field_elements_from_algebraics(rt2c, minimal=True)
(Number Field in a with defining polynomial y^2 - 2, a, Ring morphism:
From: Number Field in a with defining polynomial y^2 - 2
To: Algebraic Real Field
Defn: a |--> [1.4142135623730949 .. 1.4142135623730952])
|
|
|||
_obj
|
_obj_r
|
is_SymbolicExpressionRing
|
AA
|
QQbar
|
_mul_algo
|
_add_algo
|
_descriptors
|
AAPoly
|
_cyclotomic_gen_cache
|
QQbar_I_generator
|
AA_golden_ratio
|
| Home | Trees | Indices | Help |
|---|
| Generated by Epydoc 3.0beta1 on Thu Jul 17 04:23:29 2008 | http://epydoc.sourceforge.net |