SageMath Interface

Real Embedded Number Fields for SageMath

This wraps e-antic for SageMath providing number fields with less of a focus on number theory but fast exact ball arithmetic of the kind that is usually required for classical geometry.

class pyeantic.real_embedded_number_field.CoercionNumberFieldRenf(domain)

A coercion from RealEmbeddedNumberField to a SageMath NumberField.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: KK = RealEmbeddedNumberField(K)
sage: coercion = K.convert_map_from(KK)

TESTS:

sage: from pyeantic.real_embedded_number_field import CoercionNumberFieldRenf
sage: isinstance(coercion, CoercionNumberFieldRenf)
True
section()

Return the inverse of this coercion.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: KK = RealEmbeddedNumberField(K)
sage: K.coerce_map_from(KK).section()
Conversion map:
  From: Number Field in a with defining polynomial x^2 - 2 with a = 1.414213562373095?
  To:   Real Embedded Number Field in a with defining polynomial x^2 - 2 with a = 1.414213562373095?
class pyeantic.real_embedded_number_field.RealEmbeddedNumberField(embedded, category=None)

See RealEmbeddedNumberField in __init__.py for details.

Element

alias of RealEmbeddedNumberFieldElement

an_element()

Return a typical element in this number field.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K.an_element()
(a ~ 1.4142136)
characteristic()

Return zero, the characteristic of this number field.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: RealEmbeddedNumberField(K).characteristic()
0
degree()

Return the absolute degree of this number field.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K.degree()
2
gen()

Return the generator of this number field, i.e., a root of its defining polynomial.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K.gen()
(a ~ 1.4142136)
is_field(*args, **kwargs)

Return whether this number field is a field, i.e., return True.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K.is_field()
True
random_element()

Return a randomly generated element in this number field.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K.random_element().parent() is K
True
class pyeantic.real_embedded_number_field.RealEmbeddedNumberFieldElement(parent, value)

An element of a RealEmbeddedNumberField, i.e., a wrapper of e-antic’s renf_elem_class.

..NOTES:

This class wraps a renf_elem_class (which, at no additional runtime cost wraps a renf_elem.) At the moment it’s not possible to use a renf_elem_class directly in SageMath. Changing this might lead to a slight improvement in performance.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField  # random output due to deprecation warnings
sage: K = NumberField(x^2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: a = K.gen()

TESTS:

sage: TestSuite(a).run()

Verify that #192 has been resolved:

sage: R.<x> = QQ[]
sage: K.<b> = NumberField(x^2 - 2, embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K(b)
(b ~ 1.4142136)
sage: R.<y> = QQ[]
sage: K.<b> = NumberField(y^2 - 2, embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K(b)
(b ~ 1.4142136)

Verify that #226 has been resolved:

sage: K.<a> = NumberField(x^70 - 71*x^68 + 2414*x^66 - 52327*x^64 + 812240*x^62 - 9613968*x^60 + 90223392*x^58 \
....: - 689161713*x^56 + 4364690849*x^54 - 23231419035*x^52 + 104960312886*x^50 - 405528481605*x^48 + 1347179362620*x^46 \
....: - 3862867084860*x^44 + 9584557428600*x^42 - 20606798471490*x^40 + 38403578969595*x^38 - 61997934676405*x^36 \
....: + 86563154076490*x^34 - 104261288816825*x^32 + 107941099010360*x^30 - 95604973409176*x^28 + 72014135814704*x^26 \
....: - 45791597230002*x^24 + 24357232569150*x^22 - 10717182330426*x^20 + 3847193657076*x^18 - 1107525446734*x^16 \
....: + 250205084312*x^14 - 43138807640*x^12 + 5471263408*x^10 - 485354012*x^8 + 28001193*x^6 - 937839*x^4 + 14910*x^2 - 71, \
....: embedding=1.999510553370486)
sage: from pyeantic import RealEmbeddedNumberField
sage: RealEmbeddedNumberField(K)
Real Embedded Number Field in a with defining polynomial x^70 - 71*x^68 + 2414*x^66 - 52327*x^64 + 812240*x^62 - 9613968*x^60 + 90223392*x^58
- 689161713*x^56 + 4364690849*x^54 - 23231419035*x^52 + 104960312886*x^50 - 405528481605*x^48 + 1347179362620*x^46
- 3862867084860*x^44 + 9584557428600*x^42 - 20606798471490*x^40 + 38403578969595*x^38 - 61997934676405*x^36
+ 86563154076490*x^34 - 104261288816825*x^32 + 107941099010360*x^30 - 95604973409176*x^28 + 72014135814704*x^26
- 45791597230002*x^24 + 24357232569150*x^22 - 10717182330426*x^20 + 3847193657076*x^18 - 1107525446734*x^16
+ 250205084312*x^14 - 43138807640*x^12 + 5471263408*x^10 - 485354012*x^8 + 28001193*x^6 - 937839*x^4 + 14910*x^2 - 71
with a = 1.999510553370486?
minpoly(var='x')

Return the minimal polynomial of this element over the rationals.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K.gen().minpoly()
x^2 - 2
vector()

Return a vector representation of this element in terms of the basis of the number field.

EXAMPLES:

sage: from pyeantic import RealEmbeddedNumberField
sage: K = NumberField(x**2 - 2, 'a', embedding=sqrt(AA(2)))
sage: K = RealEmbeddedNumberField(K)
sage: K.gen().vector()
(0, 1)