Contents:
Hypercircles
This is a set of different algorithms related to the reparametrization problem and adds the class Hypercircle.
The git-aware user may use my github branch: https://github.com/lftabera/sage/tree/hypercircles
You may also download directly the module from http://personales.unican.es/taberalf/Documentos/Hypercircle.zip, unzip it and load from your Sage session:
sage: load('hypercircle.py')
sage: u=random_linear_fraction(QQ[I]['t'])
sage: H=Hypercircle([u])
In this case, ignore the import sentence of the examples and tests.
In any case, it is advisable to run at lest Sage 6.1.1 and apply the trac patch:
- patch #8558 fast gcd for polynomials with number field coefficients.
TESTS:
The following is a subtle error that can happen if, for a parameter t, Phi(t) is attained by Phibeta on two parameters, tb and infinity. This test shows that the method woks in this case:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle sage: N = NumberField(x^2-2, 'a') sage: QQx=QQ['x'] sage: D = QQx.random_element(2) sage: x = QQx.gen() sage: C = [(x**2-2)*QQx.random_element(2)/D for i in range(2)] sage: a = N.gen() sage: C = [f((a*x-a)/(x+1)) for f in C] sage: H = Hypercircle(C) sage: H.parametrization() [(1/2*x^2 + 1/2)/x, (1/4*a*x^2 - 1/4*a)/x]
This is a class representing a hypercircle for a extension QQ in QQ(alpha)
Accesing the elements of the Hypercircle are interpreting as accesing elements of the standard parametrization.
It is initialized by a proper parametrization of a curve C in QQ(alpha) represented by a list of rational functions. The hypercircle will be associated to the parametrization. In particular, if one wants to compute the hypercircle generated by a unit u, one can call Hypercircle in the parametrization [inverse_unit(u)] defined by the inverse of u.
While some features work if the ground field is different from QQ this is not assured to work.
INPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, random_linear_fraction
sage: N.<I> = QQ[I]
sage: K.<t> = N[]
sage: u = t
sage: H = Hypercircle([u])
sage: H
Hypercircle over Number Field in I with defining polynomial x^2 + 1
sage: H.parametrization()
[t, 0]
sage: H[0] # Extract a term form the parametrization
t
sage: H(2/3) # compute the point corresponding to parameter 2/3
[2/3, 0]
sage: H.compute_associated_unit(0) #See the documentation
t
sage: v = random_linear_fraction(NumberField(x^6-2,'a')['x']);
sage: H2 = Hypercircle([v])
sage: alpha = v.base_ring().gen()
sage: sum([H2[i]*alpha**i for i in range(6)])
x
Return the base field K
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = 1/(t-a)
sage: H=Hypercircle([u])
sage: H.K() is QQ
True
Return the number field K(alpha)
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = 1/(t-a)
sage: H=Hypercircle([u])
sage: H.K_alpha() is N
True
Return the primitive element of the extension.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = 1/(t-a)
sage: H=Hypercircle([u])
sage: H.alpha()
a
Return the ambient dimension of the hypercircles
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = 1/(t-a)
sage: H=Hypercircle([u])
sage: H.ambient_dimension() == N.degree()
True
Return the conic hypercircle associated to small_place_unit.
We have to be quite carful since since the ground field is QQ[alpha] although the standard parametrization is over QQ[beta].
If the hypercircle is of degree 1 or the small place is of degree 1, then returns a line.
If beta is not in QQ[alpha] then it computes the hypercircle from small_place_unit.
If beta is in QQ[alpha] it projects using relativize and then computes the conic from small_place_unit in the projection.
The result is the conic hypercircle for the extension QQ in Q[beta].
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, inverse_unit
sage: N.<a> = NumberField(x^5 + 15*x^4 + 20*x^3 + 15*x^2 + 6*x - 1)
sage: K.<t> = N['t']
sage: u = (a*t+1)/(t-2*a)
sage: H1=Hypercircle([u])
sage: H1.degree()
5
sage: H1.small_place_degree()
2
sage: conic = H1.birational_conic(); conic
Hypercircle over Number Field in beta with defining polynomial
x^2 - 2883
sage: conic.parametrization()
[(1/2*t^2 - 1776455701/1302*beta*t + 4529581007743677/2)/(t -
1776455701/1302*beta - 28998415), (1/5766*beta*t^2 -
28998415/2883*beta*t - 1509860335914559/1922*beta)/(t -
1776455701/1302*beta - 28998415)]
sage: conic.ideal('R')
Ideal (R0^2 - 2883*R1^2 - 57996830*R0*R2 + 55070126731/7*R1*R2 -
4529581007743677*R2^2) of Multivariate Polynomial Ring in R0, R1, R2
over Rational Field
Note that, since N is of odd degree, we can easily define an odd divisor in the conic.:
sage: odd = inverse_unit(H1.small_place_unit())(0);odd
(-1256591910/186889*beta + 644442070590/1308223)*a^4 +
(-19357049580/186889*beta + 9863976192090/1308223)*a^3 +
(-32933742355/186889*beta + 15916658306895/1308223)*a^2 +
(-31764862315/186889*beta + 14651617450070/1308223)*a +
15728077502199/11587118*beta + 78287276434345/2616446
sage: oddpoint = conic(odd); oddpoint
[644442070590/1308223*a^4 + 9863976192090/1308223*a^3 +
15916658306895/1308223*a^2 + 14651617450070/1308223*a +
78287276434345/2616446, -1256591910/186889*a^4 -
19357049580/186889*a^3 - 32933742355/186889*a^2 -
31764862315/186889*a + 15728077502199/11587118]
sage: pol0 = oddpoint[0].absolute_minpoly()
sage: pol0.degree()
5
sage: pol1 = oddpoint[1].absolute_minpoly()
sage: pol1.degree()
5
sage: NumberField(pol0, 'g').is_isomorphic(NumberField(pol1, 'g'))
True
An example where we have to relativize the hypercircle before computing the conic:
sage: var('x')
x
sage: N.<alpha>=NumberField(x^6-2*x^3-17,'alpha')
sage: K.<t>=N[]
sage: u = (-alpha^3*t + 1/3*alpha^3 - 1/3)/(-alpha*t + 1)
sage: H = Hypercircle([u])
sage: v = H.small_place_unit(); v
((-36*beta + 114)*alpha^2 + (-51*beta + 306)*alpha + 867)/(t +
(-51*beta + 306)*alpha^2 + 867*alpha)
sage: len(H.small_place_beta_minpoly().roots(H.K_alpha()))
2
beta is in K_alpha, so we cannot use v to compute the conic. Internally we use relativize, but this is transparent to the user, although probably slower:
sage: C = H.birational_conic(); C
Hypercircle over Number Field in beta with defining
polynomial x^2 - 2
sage: par = inverse_unit(u(v))(0); par
867/2*beta + 2601
sage: C(par)
[2601, 867/2]
Note that par is in QQ[beta][alpha]:
sage: w = C.compute_associated_unit(par[0]); w
((306*beta - 102)*t - 153*beta - 1683)/(t - beta + 8)
sage: u(v(v.parent(w)))
t + 8
Another example where beta is of degree one:
sage: N = NumberField(x^3+x+4, 'a')
sage: K = N['t']
sage: a = N.gen()
sage: t = K.gen()
sage: u = (a*t+a)/(t-1)
sage: H = Hypercircle([u])
sage: H.degree()
3
sage: C = H.birational_conic(); C
Hypercircle over Number Field in c with defining polynomial x
sage: C.degree()
1
sage: C.ambient_dimension()
1
Compute an associated unit of the hypercircle from t0 where t0 is either a parameter that gives a rational point or a list or coordinates of a rational point.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u= (t+a)/(t-a)
sage: H1 = Hypercircle([u])
sage: H2 = Hypercircle([u])
sage: H1.degree()
1
sage: H1(a)
[0, 1, 0]
sage: H1.compute_associated_unit(a)
a*t
sage: H2.compute_associated_unit([0, 1, 0])
a*t
A non linear, non primitive case:
sage: N.<a> = NumberField(x^4-2)
sage: K.<t> = N[]
False
Computes an associated unit from an odd divisor.
The hypercircle must be a plane conic or a line.
INPUT:
OUTPUT:
TESTS:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: N = NumberField(x^3+2*x+4, 'a')
sage: K = N['t']
sage: u = random_linear_fraction(K)
sage: H = Hypercircle([u])
sage: v = H.small_place_unit()
sage: v1 = inverse_unit(v)
sage: C = H.birational_conic()
sage: divisor = v1(0).minpoly()
sage: i = 1
sage: while divisor.degree() != 3:
... divisor = v1(i).minpoly()
... i = i + 1
sage: w = C.compute_associated_unit_from_odd_divisor(divisor)
sage: u1 = v(v.parent(w))
sage: drop_beta(simsim(u(u1))) in QQ['t'].fraction_field()
True
Another example in which the conic is a line because beta is in QQ:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N = NumberField(x^3+x+4, 'a')
sage: K = N['t']
sage: a = N.gen()
sage: t = K.gen()
sage: u = (t+a)/(t-a)
sage: H = Hypercircle([u])
sage: C = H.birational_conic(); C
Hypercircle over Number Field in c with defining polynomial x
An explicit one:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: u = (2*a*t-a^2)/(t+a+2-a^2)
sage: H = Hypercircle([u])
sage: v = H.small_place_unit(); v
(((1/36*beta - 5/36)*a^2 + (1/9*beta - 1/18)*a - 1/36*beta -
31/36)*t + (3190/9*beta - 30586/9)*a^2 + (11381/9*beta + 22357/9)*a
- 540*beta - 24964/3)/(t + (790/3*beta + 566/3)*a^2 + (164/3*beta -
6812/3)*a + 1580/9*beta + 209356/9)
sage: C = H.birational_conic()
sage: par = inverse_unit(v)(0); par
(-10442/81*beta + 99926/81)*a^2 + (29872/81*beta - 172624/81)*a -
44836/81*beta - 579140/81
sage: C(par)
[99926/81*a^2 - 172624/81*a - 579140/81, -10442/81*a^2 + 29872/81*a
- 44836/81]
sage: w = C.compute_associated_unit_from_odd_divisor(par.minpoly())
sage: w
((239552/99*beta + 657400/99)*t - 523552/99*beta - 504560/99)/(t -
1/11*beta - 41/11)
sage: u1 = simsim(v(v.parent(w))); u1
((41402330111/282972916566*a^2 + 251048485799/282972916566*a -
200797966037/282972916566)*t - 47061795353/141486458283*a^2 -
344736433073/141486458283*a + 263070188339/141486458283)/(t -
1248255360/15720717587*a^2 + 566134272/15720717587*a -
48682972618/15720717587)
sage: u1 = simplify_unit(u1); u1
(1/2*a*t - 2*a^2 - 8*a + 2)/(t + a^2 - 7)
sage: H(u1(0))
[-22/27, 25/27, 4/27]
sage: simsim(u(u1))
8/(t - 8)
Compute a place of degree 1 or 2 of the hypercircle.
Note that, since the method uses LLL to basis reduction, it only works for absolute number fields. In particular, the hypercircle must be defined over QQ.
The method tries to help with the hell of number fields. If we start with a number field QQ[alpha]. It will compute a field QQ[beta] such that the hypercircle has points in QQ[beta]. These structures are incompatible, so it also computes a new field QQ[gamma] = QQ[alpha,beta] as well as relative representations QQ[alpha][beta], QQ[beta][alpha] and isomorphisms with QQ[gamma]
NOTE: With improvements in coercion, there may be some morphisms that are not needed.
INPUT:
OUPUT:
If the place is of degree 1, a dictionary with the following keys:
If the place is of degree 2, a dictionary with the following keys:
EXAMPLES:
This used to fail:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, random_linear_fraction
sage: var('x')
x
sage: N.<alpha> = NumberField(x^3-2,'alpha')[x]
sage: u = random_linear_fraction(N, reduced = True)
sage: w = Hypercircle([u])
sage: w.compute_small_place()
sage: P = w._small_place_structure()
sage: P['degree_place']
1
sage: [foo in QQ for foo in P['rational_point_witness']]
[True, True, True]
sage: [p(P['parameter_to_QQ']) for p in w] == P['rational_point_witness']
True
sage: p = P['parameter_to_QQ']
sage: u = random_linear_fraction(N)
sage: w = Hypercircle([u])
sage: w.compute_small_place()
sage: P = w._small_place_structure()
sage: sorted(P.keys())
['K_alpha_beta', 'K_alpha_beta_to_NewK', 'K_alpha_to_NewK', 'K_beta',
'K_beta_alpha', 'K_beta_alpha_to_NewK', 'K_beta_to_NewK', 'NewK',
'NewK_to_K_alpha_beta', 'NewK_to_K_beta_alpha', 'W_D', 'beta',
'degree_place', 'gamma', 'parameter_to_beta_in_K_beta_alpha',
'parameter_to_beta_in_NewK', 'point_beta_in_K_beta_alpha',
'point_beta_in_NewK', 'witness_in_K_beta_alpha', 'witness_in_NewK']
sage: len(P['point_beta_in_K_beta_alpha'])
3
sage: P['beta'].minpoly().degree()
2
sage: Phi_ba = P['witness_in_K_beta_alpha']
sage: Kba = P['K_beta_alpha']
sage: parameter = P['parameter_to_beta_in_K_beta_alpha']
sage: point_beta = P['point_beta_in_K_beta_alpha']
sage: beta = P['beta']
sage: x = Phi_ba[0].parent().gen()
sage: W = P['witness_in_K_beta_alpha']
sage: [foo(parameter) for foo in W] == point_beta
True
sage: str(alpha) in str(point_beta)
False
sage: str(beta) in str(point_beta)
True
Rational points can be found during the process:
sage: var('x')
x
sage: N.<a> = NumberField(x^3-2)
sage: K.<t> = N[]
sage: u = (a*t+1)/t
sage: H = Hypercircle([u])
sage: H.compute_small_place()
sage: P = H._small_place_structure()
sage: P['degree_place']
1
sage: P['rational_point_witness']
[0, 0, 0]
A non-reduced case with a rational point found:
sage: var('x')
x
sage: N.<a> = NumberField(x^3-2)
sage: K.<t> = N[]
sage: u = (a*t - a^2 + 1)/(t - a)
sage: H = Hypercircle([u])
sage: H(0)
[70/433*a^2 - 145/433*a + 22/433, 13/433*a^2 + 4/433*a + 301/433,
-78/433*a^2 - 24/433*a - 74/433]
sage: H.compute_small_place()
sage: P = H._small_place_structure()
sage: P['degree_place']
1
sage: P['rational_point_witness']
[0, 1, 0]
TEST:
This used to fail:
sage: var('x')
x
sage: N.<a> = NumberField(x^4-5, 'a')
sage: K.<x> = N[]
sage: u = (a*x-(a**2+1))/((4-a**2)*x+(a**2+3))
sage: S = Hypercircle([u])
sage: S.compute_small_place()
EXAMPLES:
sage: var('x')
x
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = ((a-1)*t+a+3)/(a**2*t-a)
sage: H = Hypercircle([u])
sage: H.compute_small_place()
sage: H.small_place_degree()
2
sage: H.small_place_coordinates()
[13/8*beta - 121/8, -5/8*beta + 41/8, 1/2*beta - 7/2]
sage: H.small_place_beta_minpoly()
x^2 - 73
sage: t0 = H.small_place_parameter(); t0
(1/2*beta - 7/2)*a^2 + (-5/8*beta + 41/8)*a + 13/8*beta - 121/8
sage: H(t0)
[13/8*beta - 121/8, -5/8*beta + 41/8, 1/2*beta - 7/2]
sage: v = random_linear_fraction(K, reduced='true')
sage: H1 = Hypercircle([v,v])
sage: H1.compute_small_place(beta_name = 'jj', gamma_name = 'gg')
sage: H1.small_place_beta()
1
Return the degree of the hypercircle
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = (a*t+a)/(t-a)
sage: H=Hypercircle([u])
sage: H.degree()
3
sage: u = (t+a)/(t-a)
sage: H = Hypercircle([u])
sage: H.degree()
1
Return the implicit ideal of the hypercircle
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = (a*t+a)/(t-a)
sage: H=Hypercircle([u])
sage: I = H.ideal(); I
Ideal (Y0^2 - 2*Y1*Y2 + Y0*Y3 + 2*Y2*Y3, Y0*Y1 - 2*Y2^2, -Y1^2 +
Y0*Y2 + Y1*Y3 + Y2*Y3) of Multivariate Polynomial Ring in Y0, Y1,
Y2, Y3 over Rational Field
sage: [foo(Y0=H[0],Y1=H[1],Y2=H[2],Y3=1) for foo in I.gens()]
[0, 0, 0]
sage: H.ideal(name='x')
Ideal (x0^2 - 2*x1*x2 + x0*x3 + 2*x2*x3, x0*x1 - 2*x2^2, -x1^2 +
x0*x2 + x1*x3 + x2*x3) of Multivariate Polynomial Ring in x0, x1,
x2, x3 over Rational Field
Check if the hypercircle is a conic.
TEST:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, random_linear_fraction, inverse_unit
sage: N = QQ[I]
sage: K = N['t']
sage: u1 = random_linear_fraction(K)
sage: H1 = Hypercircle([u1])
sage: H1.is_conic()
True
sage: u2 = random_linear_fraction(K, reduced=True)
sage: H2 = Hypercircle([inverse_unit(u2)])
sage: H2.is_conic()
False
Return wether the hypercircle is a line or not.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = (a*t+a)/(t-a)
sage: H = Hypercircle([u])
sage: H.is_line()
False
sage: H = Hypercircle([t])
sage: H.is_line()
True
Check if the hypercircle is primitive
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = (a*t+a)/(t-a)
sage: H=Hypercircle([u])
sage: H.is_primitive()
True
sage: u = 1/(t-a)
sage: H=Hypercircle([u])
sage: H.is_primitive()
False
Return the standar parametrization of the hypercircle
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = 1/(t-a)
sage: H=Hypercircle([u])
sage: H.parametrization()
[t - a, 1, 0]
Return the minimal polynomial of alpha
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = 1/(t-a)
sage: H=Hypercircle([u])
sage: H.polmin()
x^3 - 2
Return a rational parametrization if precomputed or if we already have a unit. See also compute_rational_parametrization. The result is cached
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, inverse_unit
sage: N.<a> = NumberField(x^4+3)
sage: K.<t> = N[]
sage: u = ((2*a+a**3)*t+1)/((1-a)*t-2*a)
sage: H = Hypercircle([u])
sage: H.set_unit(inverse_unit(u))
sage: H.rational_parametrization()
[(-3/2*t^4 - 41/4*t^3 - 21*t^2 - 18*t)/(t^4 + 9*t^3 + 57/2*t^2 +
42*t + 147/4), (1/2*t^4 + 13/4*t^3 + 29/4*t^2 + 3*t + 21/4)/(t^4 +
9*t^3 + 57/2*t^2 + 42*t + 147/4), (1/2*t^4 + 11/4*t^3 + 7*t^2 +
43/4*t)/(t^4 + 9*t^3 + 57/2*t^2 + 42*t + 147/4), (1/2*t^4 + 9/4*t^3
+ 9/4*t^2 + 15/4*t + 7/2)/(t^4 + 9*t^3 + 57/2*t^2 + 42*t + 147/4)]
If beta is an algebraic integer of K_alpha, compute the hypercircle associated to K(beta) in K(alpha).
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^6 + 6*x^5 + 15*x^4 + 20*x^3 + 15*x^2 + 6*x - 1)
sage: K.<t> = N['t']
sage: u = (a*t+1)/(t-a)
sage: H=Hypercircle([u])
sage: H.degree()
6
sage: H.small_place_degree()
2
sage: beta = a^3 + 3*a^2 + 3*a
sage: beta.minpoly()
x^2 + 2*x - 1
If the hypercircle is of degree 6 and beta is of degree 2, then the hypercircle associated to u with base field K[beta] must be of degree 3.:
sage: H2, phi, mat = H.relativize(beta)
sage: H2.degree()
3
sage: H2.ambient_dimension()
3
sage: phi
Ring morphism:
From: Number Field in a with defining polynomial x^6 + 6*x^5 +
15*x^4 + 20*x^3 + 15*x^2 + 6*x - 1
To: Number Field in a with defining polynomial x^3 + 3*x^2 + 3*x -
beta over its base field
Defn: a |--> a
sage: mat
[ 1 0 0 beta -3*beta
6*beta]
[ 0 1 0 -3 beta + 9
-3*beta - 18]
[ 0 0 1 -3 6
beta - 9]
sage: N2 = phi.codomain(); N2
Number Field in a with defining polynomial x^3 + 3*x^2 + 3*x - beta
over its base field
sage: N.is_isomorphic(N2)
True
Compute directly the hypercircle over the field N2:
sage: uu = N2[t].fraction_field()(u); uu
(a*t + 1)/(t - a)
sage: H3 = Hypercircle([uu])
sage: H2.parametrization() == H3.parametrization()
True
Check that the matrix mat is the isomorphism searched:
sage: t0 = N.random_element()
sage: mat * vector(H(t0)) == vector(H2(phi(t0)))
True
Declare an associated unit to the hypercircle. If check is True, compute also an associated rational parametrization and cache it. If simplify is True, an equivalent unit with possibly smaller coefficients will be used instead.
As a side effect, it will change the internal unit and, as a side effect if check = True, change the internar rational parametrization.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, inverse_unit
sage: N.<a> = NumberField(x^3-2)
sage: K.<t> = N[]
sage: u = (t+a)/(t-a)
sage: v = inverse_unit(u)
sage: H = Hypercircle([v])
sage: H.set_unit(t) # not an associated unit
Traceback (most recent call last):
...
TypeError: not a constant polynomial
sage: H.set_unit(u)
sage: H.rational_parametrization()
[(t^3 + 2)/(t^3 - 2), 2*t^2/(t^3 - 2), 2*t/(t^3 - 2)]
sage: H.set_unit((a*t + 1)/(-a*t + 1)) # another unit
sage: H.rational_parametrization()
[(-t^3 - 1/2)/(t^3 - 1/2), -t/(t^3 - 1/2), -t^2/(t^3 - 1/2)]
Return the primitive element of QQ[beta] used to define this field.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = ((a-1)*t+a+3)/(a**2*t-a)
sage: H = Hypercircle([u])
sage: H.small_place_beta()
beta
sage: HH = Hypercircle([u])
sage: HH.compute_small_place(beta_name = 'jj')
sage: HH.small_place_beta()
jj
sage: v = [((59*a^2 + 21*a + 84)*t + 47*a^2 + 54*a + 54)/t]
sage: H1 = Hypercircle(v)
sage: H1.small_place_beta()
1
Retuns the minimal polynomia of beta that is the defining polynomial of QQ[beta].
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = ((a-1)*t+a+3)/(a**2*t-a)
sage: H = Hypercircle([u])
sage: H.small_place_beta_minpoly()
x^2 - 73
sage: v = [((59*a^2 + 21*a + 84)*t + 47*a^2 + 54*a + 54)/t]
sage: H1 = Hypercircle(v)
sage: H1.small_place_beta_minpoly()
x - 1
Return a list representing the coordinates of a point in an extension of degree at most 2 over QQ
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = ((a-1)*t+a+3)/(a**2*t-a)
sage: H = Hypercircle([u])
sage: t0=H.small_place_parameter(); t0
(1/2*beta - 7/2)*a^2 + (-5/8*beta + 41/8)*a + 13/8*beta - 121/8
sage: H(t0)
[13/8*beta - 121/8, -5/8*beta + 41/8, 1/2*beta - 7/2]
sage: H.small_place_coordinates()
[13/8*beta - 121/8, -5/8*beta + 41/8, 1/2*beta - 7/2]
sage: v = [((59*a^2 + 21*a + 84)*t + 47*a^2 + 54*a + 54)/t]
sage: H1 = Hypercircle(v)
sage: t1=H1.small_place_parameter(); t1
2255/26446*a^2 + 580/13223*a + 1735/26446
sage: H1(t1)
[1735/26446, 580/13223, 2255/26446]
sage: H1.small_place_coordinates()
[1735/26446, 580/13223, 2255/26446]
Return the degree of a small place computed. The degree is garanteed to be 1 or 2.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = ((a-1)*t+a+3)/(a**2*t-a)
sage: H = Hypercircle([u])
sage: H.compute_small_place()
sage: H.small_place_degree()
2
sage: v = [((59*a^2 + 21*a + 84)*t + 47*a^2 + 54*a + 54)/t]
sage: H1 = Hypercircle(v)
sage: H1.compute_small_place()
sage: H1.small_place_degree()
1
Returns true is self.small_place_beta() defines a subfield of K[alpha]
Examples:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, random_linear_fraction, inverse_unit
sage: N.<alpha>=NumberField(x^6-2*x^3-17,'alpha')
sage: K.<t>=N[]
sage: v = (-alpha^3*t + 1/3*alpha^3 - 1/3)/(-alpha*t + 1)
sage: H = Hypercircle([v])
sage: H.small_place_in_subfield()
True
sage: v = ((alpha+1)*t-alpha**4)/(t+alpha**2-alpha)
sage: H = Hypercircle([v])
sage: H.small_place_in_subfield()
False
sage: v = random_linear_fraction(K, reduced=True)
sage: H = Hypercircle([ inverse_unit(v)])
sage: H.small_place_in_subfield()
True
Return a parameter in QQ(beta, alpha) such that the image, by the parametrization, is in QQ(beta)
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = ((a-1)*t+a+3)/(a**2*t-a)
sage: H = Hypercircle([u])
sage: t0=H.small_place_parameter(); t0
(1/2*beta - 7/2)*a^2 + (-5/8*beta + 41/8)*a + 13/8*beta - 121/8
sage: H(t0)
[13/8*beta - 121/8, -5/8*beta + 41/8, 1/2*beta - 7/2]
sage: v = [((59*a^2 + 21*a + 84)*t + 47*a^2 + 54*a + 54)/t]
sage: H1 = Hypercircle(v)
sage: t1=H1.small_place_parameter(); t1
2255/26446*a^2 + 580/13223*a + 1735/26446
sage: H1(t1)
[1735/26446, 580/13223, 2255/26446]
Return a unit that reparametrizes the hypercircle over QQ[beta].
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, simsim
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = ((a-1)*t+a+3)/(a**2*t-a)
sage: H = Hypercircle([u])
sage: ubeta = H.small_place_unit(); ubeta
(((1/2*beta - 7/2)*a^2 + (-5/8*beta + 41/8)*a + 13/8*beta - 121/8)*t
+ (255/2*beta - 1259)*a^2 + (384*beta + 12803/2)*a + 16945/8*beta +
74901/8)/(t + (4461/8*beta + 44765/8)*a^2 + (17635/8*beta +
141935/8)*a + 1487/2*beta + 9812)
sage: map(simsim, H(ubeta))
[((13/8*beta - 121/8)*t^3 + (183391/16*beta + 1134121/16)*t^2 +
(350112143/16*beta + 2855716015/16)*t + 434921034785763/128*beta +
3715851728706871/128)/(t^3 + 14107/2*t^2 + (15409306305/32*beta +
131798889317/32)*t - 225946444725621/32*beta - 1930463829601393/32),
((-5/8*beta + 41/8)*t^3 + (32285/4*beta + 160159/2)*t^2 +
(-6579417895/32*beta - 56095033211/32)*t - 93993070662403/128*beta -
803026632222679/128)/(t^3 + 14107/2*t^2 + (15409306305/32*beta +
131798889317/32)*t - 225946444725621/32*beta - 1930463829601393/32),
((1/2*beta - 7/2)*t^3 + (34803/16*beta + 174853/16)*t^2 +
(-1326724227/32*beta - 11418602135/32)*t + 200770410597589/128*beta
+ 1715344616464353/128)/(t^3 + 14107/2*t^2 + (15409306305/32*beta +
131798889317/32)*t - 225946444725621/32*beta - 1930463829601393/32)]
sage: v = [((59*a^2 + 21*a + 84)*t + 47*a^2 + 54*a + 54)/t]
sage: H1 = Hypercircle(v)
sage: vbeta = simsim(H1.small_place_unit()); vbeta
((2255/26446*a^2 + 580/13223*a + 1735/26446)*t + 79827/26446*a^2 +
20532/13223*a + 61419/26446)/(t - 45347/13223*a^2 + 17819/26446*a -
3841/26446)
sage: par_base_field = map(simsim, H1(vbeta)); par_base_field
[(1735/26446*t^3 + 144843/52892*t^2 + 974961/26446*t +
41451453/52892)/(t^3 + 50179/3778*t^2 + 495569/52892*t -
28845451/26446), (580/13223*t^3 + 23773/52892*t^2 - 574608/13223*t -
8236341/52892)/(t^3 + 50179/3778*t^2 + 495569/52892*t -
28845451/26446), (2255/26446*t^3 + 99807/26446*t^2 + 1303349/52892*t
- 3937719/52892)/(t^3 + 50179/3778*t^2 + 495569/52892*t -
28845451/26446)]
sage: sum([par_base_field[i] * a**i for i in range(3)]) == vbeta
True
sage: simsim(v[0](vbeta))
(2867/5*t + 13706/5)/(t + 177/5)
sage: H2 = Hypercircle([t])
sage: w = H2.small_place_unit(); w
t
The following example used to fail, It is a primitive hypercircle such that we find a point on a subfield of degree 2. Hence, it needs relativize:
sage: var('x')
x
sage: N.<alpha>=NumberField(x^6-2*x^3-17,'alpha')
sage: K.<t>=N[]
sage: v = (-alpha^3*t + 1/3*alpha^3 - 1/3)/(-alpha*t + 1)
sage: H = Hypercircle([v])
sage: H.small_place_unit()
((-36*beta + 114)*alpha^2 + (-51*beta + 306)*alpha + 867)/(t +
(-51*beta + 306)*alpha^2 + 867*alpha)
sage: parbeta = map(simsim,H(H.small_place_unit()))
sage: S = str(parbeta)
sage: 'beta' in S
True
sage: 'alpha' in S
False
sage: parbeta[0]
((136/699*beta + 323/699)*t^6 + (-169932/233*beta - 100572/233)*t^5
+ (-202469643/233*beta + 76760712/233)*t^4 + (224927147781/233*beta
+ 680118435243/233)*t^3 + (-310588970748192/233*beta -
321170435162694/233)*t^2 + (81307917817026876/233*beta +
107073950888577564/233)*t - 5484690082401920928/233*beta -
14990067970691614488/233)/(t^6 + (106488/233*beta + 359856/233)*t^5
+ (-1583088246/233*beta - 1056775896/233)*t^4 +
(677541778002/233*beta + 1127117771766/233)*t^3 +
(3177467265953376/233*beta + 5570316058860162/233)*t^2 +
(-6201603752918368986/233*beta - 9890189581198153140/233)*t +
3701691290153002346406/233*beta + 5448920353023639800187/233)
This used to fail:
sage: N.<alpha> = NumberField(x^4 - 26*x^2 + 49, 'alpha')
sage: K.<t> = N[]
sage: Phi = [((-5/7*alpha^3 + 165/7*alpha)*t^3 + (45*alpha^3 + 15*alpha^2 + 315*alpha + 135)*t^2 + (31266/7*alpha^3 + 3150*alpha^2 - 44358/7*alpha - 4410)*t + 110952*alpha^3 + 116454*alpha^2 - 219996*alpha - 231556)/(25*t^4 + (-5/7*alpha^3 + 300*alpha^2 + 865/7*alpha)*t^3 + (945*alpha^3 + 35265*alpha^2 + 315*alpha - 65955)*t^2 + (523382/7*alpha^3 + 1719810*alpha^2 - 970146/7*alpha - 3488310)*t + 1811964*alpha^3 + 31409338*alpha^2 - 3674832*alpha - 64192860), (50*t^3 + (-15/14*alpha^3 + 450*alpha^2 + 2595/14*alpha)*t^2 + (945*alpha^3 + 35265*alpha^2 + 315*alpha - 66255)*t + 261706/7*alpha^3 + 859005*alpha^2 - 487668/7*alpha - 1744155)/(25*t^4 + (-5/7*alpha^3 + 300*alpha^2 + 865/7*alpha)*t^3 + (945*alpha^3 + 35265*alpha^2 + 315*alpha - 65955)*t^2 + (523382/7*alpha^3 + 1719810*alpha^2 - 970146/7*alpha - 3488310)*t + 1811964*alpha^3 + 31409338*alpha^2 - 3674832*alpha - 64192860)]
sage: H = Hypercircle(Phi)
sage: H.degree() # It is not primitive
2
sage: H.small_place_beta_minpoly()
x^2 + 2
sage: H.small_place_unit()
((3/2870*alpha^3 - 3*alpha^2 - 2969/2870*alpha - 61/41*beta)*t +
(45/328*beta - 57/2870)*alpha^3 + (15/328*beta + 99/82)*alpha^2 +
(315/328*beta + 1518/1435)*alpha + 357/328*beta)/(t +
15/2296*beta*alpha^3 - 495/2296*beta*alpha - 33/82)
sage: map(simsim, H(H.small_place_unit()))
[(-61/41*beta*t^2 + 111/82*beta*t - 9/164*beta)/(t^2 - 33/41*t +
27/82), (-2969/2870*t^2 + 3036/1435*t - 10611/11480)/(t^2 - 33/41*t
+ 27/82), -3, (3/2870*t^2 - 57/1435*t + 207/11480)/(t^2 - 33/41*t +
27/82)]
Returns the standard parametrization of the hypercircle.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3+2*x+5)
sage: K.<t>=N[]
sage: u = 1/(t-a**3)
sage: H=Hypercircle([u])
sage: H.standard_parametrization()
[t + 2*a, -2, 0]
Return the variable of the parametrization
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle
sage: N.<a> = NumberField(x^3-2)
sage: K.<s>=N[]
sage: u = 1/(s-a)
sage: H = Hypercircle([u])
sage: H.t()
s
sage: H2 = Hypercircle([u], name='var9')
sage: H2.t()
var9
Return an associated unit
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import Hypercircle, simsim
sage: N.<a> = NumberField(x^3-2)
sage: K.<t>=N[]
sage: u = (a*t+a)/(t-a)
sage: H = Hypercircle([u])
sage: H.unit()
Traceback (most recent call last):
...
ValueError: Associated unit not yet discovered
sage: H.compute_associated_unit(0)
(-a^2 + a + 2)/(t + a^2 - a)
sage: [simsim(P(H.unit())) for P in H]
[(2*t^2 - 4*t + 2)/(t^3 + 6*t + 2), (t^2 + 4*t + 4)/(t^3 + 6*t + 2),
(-t^2 - t + 2)/(t^3 + 6*t + 2)]
Write a multivariate polynomial P in K(alpha)[t] as a vector in K[t]. This is the inverse of sums_alpha.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import alpha_components, sums_alpha
sage: N.<a> = QQ[2^(1/3)]
sage: K.<x,y> = N['x,y']
sage: f = K.random_element(5)
sage: P = alpha_components(f)
sage: P[0] in QQ['x,y']
True
sage: f - sums_alpha(P, 1, a)
0
TESTS:
sage: K1=PolynomialRing(QQ[sqrt(2)], ('t','x','y'), order=TermOrder('degrevlex', 1) + TermOrder('degrevlex',2))
sage: t,x,y=K1.gens()
sage: a=K1.base_ring().gen()
sage: a**2
2
sage: f=a*x-y
sage: alpha_components(f)
(-y, x)
Compute the composition of the parametrization represented by Num and Den with unit.
This is optimized for dealing with number fields as coefficients.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import composition_by_unit, random_linear_fraction
sage: N = NumberField(x^3-2,'a')
sage: K.<t> = N[]
sage: Num = [K.random_element(3) for i in range(7)]
sage: Den = 0
sage: while Den == 0: Den = K.random_element(5)
...
sage: unit = random_linear_fraction(K)
sage: NN, DD = composition_by_unit(Num, Den, unit)
sage: [NN[i]/DD - (Num[i]/Den)(unit) for i in range(7)]
[0, 0, 0, 0, 0, 0, 0]
Compute the conjugate of a polynomial.
INPUT:
OUTPUT:
K[t] -> L[t] or K(t) -> L(t)
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import conjugate_pol
sage: N.<I> = QQ[I]
sage: conj = N.automorphisms()[1]
sage: conj
Ring endomorphism of Number Field in I with defining polynomial x^2 + 1
Defn: I |--> -I
sage: K.<x> = N[]
sage: f = (1-I)*x-(2+3*I)
sage: conjugate_pol(f, conj, K)
(I + 1)*x + 3*I - 2
sage: f = K.random_element(20)
sage: f + conjugate_pol(f, conj, K) in QQ[x]
True
sage: I*(f - conjugate_pol(f, conj, K)) in QQ[x]
True
A more interesting example:
sage: var('y')
y
sage: N.<a> = NumberField(y^5-2, 'a')
sage: K.<x> =N[]
sage: f = x^4 + a*x^3 + a^2*x^2 + a^3*x + a^4
sage: M.<b> = NumberField(f)
sage: conj1 = N.hom([M(a)])
sage: conj1
Ring morphism:
From: Number Field in a with defining polynomial y^5 - 2
To: Number Field in b with defining polynomial x^4 + a*x^3 + a^2*x^2 +
a^3*x + a^4 over its base field
Defn: a |--> a
sage: conj2 = N.hom([1/2*a^4*b^2])
sage: conj2
Ring morphism:
From: Number Field in a with defining polynomial y^5 - 2
To: Number Field in b with defining polynomial x^4 + a*x^3 + a^2*x^2 +
a^3*x + a^4 over its base field
Defn: a |--> 1/2*a^4*b^2
sage: p = a + x + a^3*x^2
sage: p1 = conjugate_pol(p, conj1, M[x]); p1
a^3*x^2 + x + a
sage: p2 = conjugate_pol(p, conj2, M[x]); p2
a^2*b*x^2 + x + 1/2*a^4*b^2
sage: sage.rings.polynomial.all.is_Polynomial(p1)
True
sage: sage.rings.polynomial.all.is_Polynomial(p2)
True
sage: q = p / p.parent().one(); q
a^3*x^2 + x + a
sage: q1 = conjugate_pol(q, conj1, M[x])
sage: sage.rings.polynomial.all.is_Polynomial(q1)
False
Take a fraction in N(beta)(t) but whose coefficients are in N. Output de same rational fraction in N(t).
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import drop_beta
sage: N = QQ[sqrt(2), sqrt(3)]
sage: s2, s3 = N.gens()
sage: K.<x> = N[]
sage: u = (s3*x^2 + (1 - s3)*x +1) / (s3 * x^2 +1 -3*s3)
sage: v = drop_beta(u)
sage: u == v
True
sage: u.parent() is K.fraction_field()
True
sage: v.parent() is K.fraction_field()
False
sage: v.parent() is N.base_ring()[x].fraction_field()
True
Compute the implicit ideal of an hypercircle using the algorithm in Fast computation of the implicit ideal of a hypercircle.
INPUT:
OUTPUT:
One a set of generators is computed, perform a reduction to return less generators. echelon uses linear algebra, grobner use groebner bases, dummy does not perform any reduction and return all generators.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: N.<I> = QQ[I]
sage: K.<x> = N[]
sage: u = (I*x-2)/(x+1+3*I)
sage: w = witness([u])
sage: implicite_hc(w)
Ideal (Y0^2 + Y1^2 + Y0*Y2 + 5*Y1*Y2 + 6*Y2^2) of Multivariate
Polynomial Ring in Y0, Y1, Y2 over Rational Field
sage: var('x')
x
sage: N.<a> = NumberField(x^3-2)
sage: K = N['p']
sage: u = random_linear_fraction(K)
sage: w = witness([u])
sage: I = implicite_hc(w)
sage: L = I.gens()
sage: len(L)
3
sage: L[0](w[0], w[1], w[2], 1)
0
sage: L[1](w[0], w[1], w[2], 1)
0
sage: L[2](w[0], w[1], w[2], 1)
0
A non-primitive case:
sage: var('x')
x
sage: N.<alpha> = NumberField(x^4-2)
sage: K.<t> = N[]
sage: u = ((1-alpha**2)*t+(1+alpha**2*3) ) /((1+3*alpha**2)*t+5*alpha**2-7)
sage: H = witness([u])
sage: I = implicite_hc(H)
sage: I
Ideal (Y0^2 - 2*Y2^2 - 1/2*Y0*Y4 - 5*Y2*Y4 + 13/2*Y4^2, Y1, Y3) of
Multivariate Polynomial Ring in Y0, Y1, Y2, Y3, Y4 over Rational Field
Another examples, where the linear part is not constant:
sage: var('x')
x
sage: N.<a> = NumberField(x^6 + 6*x^5 + 15*x^4 + 20*x^3 + 15*x^2 + 6*x - 1)
sage: K.<t> = N['t']
sage: u = ((-a^3 - 3*a^2 - 3*a + 1)*t - 2*a^3 - 6*a^2 - 6*a - 1)/((13*a^3 + 39*a^2 + 39*a + 14)*t - 7*a^3 - 21*a^2 - 21*a - 11)
sage: H = witness([u])
sage: I = implicite_hc(H)
sage: I = I.gens()
sage: [I[i](H[0], H[1], H[2], H[3], H[4], H[5], 1) for i in range(len(I))]
[0, 0, 0, 0, 0]
sage: I[1]
Y1 - 3*Y3
REFERENCE:
AGGM algorithm
Compute the inverse of a linear fraction.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import inverse_unit, random_linear_fraction
sage: K.<a,b,c,d> = QQ['a,b,c,d']
sage: L.<x> = K.fraction_field()['tr']
sage: u = (a*x+b)/(c*x+d)
sage: inverse_unit(u)
(-d*tr + b)/(c*tr - a)
Note that the variable may be different from the variable of u:
sage: K1.<xx> = QQ['xx']
sage: K2.<t> = QQ['t']
sage: u = random_linear_fraction(K1)
sage: v = inverse_unit(u, t)
sage: u(v)
t
sage: v(u)
xx
Check if u is a linear fraction over a univariate polynomial rings.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import is_com_unit
sage: K.<x> = QQ[]
sage: u = (2*x+1)/(4*x+2)
sage: is_com_unit(u)
False
sage: u = (2*x+1)/(1*x)
sage: is_com_unit(u)
True
sage: u = K.fraction_field()(x)
sage: is_com_unit(u)
True
sage: u = K.fraction_field()(0)
sage: is_com_unit(u)
False
sage: u = K.fraction_field()(1+x^3)
sage: is_com_unit(u)
False
Check if a parametrization is a parametrization of a hypercircle.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: N = NumberField(x^3-2, 'a')
sage: K.<t> = N[]
sage: u = random_linear_fraction(K)
sage: Phi = witness([u], name=t)
sage: v = is_hypercircle(Phi, inverse_unit(u,t)(0))
sage: simsim(u(v)) in QQ[t].fraction_field()
True
TODO:
Support for non-primitive hypercircles, more tests
REFERENCE:
ACA algorithm
Auxiliary function to is_hypercircle. It accepts the following data:
INPUT:
Such that the parametrization Num / Den is a unit parametrization passing through the origin in standard form.
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import is_hypercircle_unit_reduced_standar
sage: N.<a> = QQ[I]
sage: K.<t> = N[]
sage: d = 2
sage: Num = 29*t^2 - 26*a*t, -29*a*t^2
sage: Den = 58*t - 26*a
sage: is_hypercircle_unit_reduced_standar(Num, Den, K, t, a, d)
1/(t - 29/26*I)
TODO:
Support for non-primitive hypercircles, more tests
REFERENCE:
ACA algorithm
Check if a parametrization is a parametrization of a reduced hypercircle.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: N = NumberField(x^3-2, 'a')
sage: K.<t> = N[]
sage: u = random_linear_fraction(K, reduced = True)
sage: Phi = witness([u], name=t)
sage: v = is_reduced_hypercircle(Phi)
sage: simsim(u(v)) in QQ[t].fraction_field()
True
TODO:
Support for non-primitive hypercircles, more tests
REFERENCE:
ACA algorithm
Compute the gcd of a list of polys.
The difference with the general gcd is that, once a gcd is computed, before computed next, check if the gcd already divides next element. In most cases, for long lists, one only needs to check a few (expect one) gcd to obtain the gcd of the whole list, and division is cheaper, at least in the conjugation process of the hypercircle model.
INPUT:
OUTPUT:
EXAPLES:
sage: from sage.geometry.hypercircles.hypercircle import my_gcd
sage: K.<x> = QQ[x]
sage: f = [K.random_element() for i in range(5)]
sage: gcd(f) - my_gcd(f)
0
sage: h = K.random_element()
sage: f = [ h*m for m in f]
sage: gcd(f) - my_gcd(f)
0
Algorithm to compute the inverse of an element over an absolute number field. This is faster than current sage code for big number fields with big defining coefficients and big s, maybe due to flint faster xgcd.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import my_inverse_big_absnfield
sage: N = NumberField(x^7+1213451*x^4 -135156164614,'a')
sage: c = 0
sage: while c.is_zero(): c = N.random_element()
...
sage: my_inverse_big_absnfield(c) * c
1
Compute the monic lcm of a list of polynomials, avoiding calling to Singular.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import my_lcm
sage: K.<x> = NumberField(x^3-2, 'a')[]
sage: f = [K.random_element(3) for i in range(3)]
sage: l = my_lcm(f)
sage: [ l % m for m in f]
[0, 0, 0]
sage: f = [x-1, x+1, x^2-1]
sage: my_lcm(f)
x^2 - 1
TEST:
sage: from sage.geometry.hypercircles.hypercircle import my_quo_rem sage: N = QQ[I] sage: K.<t> = N[] sage: f = K.random_element(10) sage: g = K.random_element(10) sage: f.quo_rem(g) == my_quo_rem(f,g) True
Compute Newton sums of the generator of a NumberField
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import newton_sums
sage: N=NumberField(x^3+x^2-1,'b')
sage: newton_sums(N)
[3, -1, 1]
sage: f = x^5-3
sage: N.<b> = NumberField(f)
sage: vector(newton_sums(N)) - vector((b**i).trace() for i in range(5))
(0, 0, 0, 0, 0)
Take a list of rational functions representing a parametrization and reduce it to common denominator.
INPUT:
OUTPUT:
The length of Num equals the length of Phi and Phi[i] = Num[i]/Den
gcd of Num[0], ..., Num[-1], Den is 1
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import parametrization_to_common_denominator
sage: K = QQ[x].fraction_field()
sage: Phi = [K.random_element(3) for i in range(5)]
sage: Num, Den = parametrization_to_common_denominator(Phi)
sage: [Num[i] / Den - Phi[i] for i in range(5)]
[0, 0, 0, 0, 0]
sage: from sage.rings.polynomial.polynomial_element import is_Polynomial
sage: is_Polynomial(Den)
True
sage: [is_Polynomial(i) for i in Num]
[True, True, True, True, True]
Return a linear fraction on the fraction field of K
INPUT:
OUPUT:
A linear fraction u with monic denominator:
Warning
Note that if the base_ring is of relative degree 1, reduced is ignored, since this does not make sense.
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: L.<t> = NumberField(x^3-2, 'a')[]
sage: u = random_linear_fraction(L)
sage: is_com_unit(u)
True
sage: Phi = witness([u])
sage: l = is_hypercircle(Phi, inverse_unit(u,t)(0))
sage: is_com_unit(l)
True
sage: is_reduced_hypercircle(Phi)
Traceback (most recent call last):
...
ValueError: Not reduced hypercircle
sage: v = random_linear_fraction(L, reduced=True)
sage: is_com_unit(v)
True
sage: Psi = witness([v])
sage: l = is_reduced_hypercircle(Psi)
sage: is_com_unit(l)
True
Compute a parameter t0 such that Phi(t0) is a rational point:
INPUT:
OUPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import rational_point_conic
sage: N.<beta, alpha> = NumberField([x^2+23, x^3+2], 'beta, alpha')
sage: K.<t> = N[]
sage: conic = [(1/2*t^2 - 10/69*beta*t - 49/3)/(t - 10/69*beta - 20/3), (-1/46*beta*t^2 + 20/69*beta*t - 49/69*beta)/(t - 10/69*beta - 20/3)]
sage: t_odd = (-1/108*alpha^2 + 8/27*alpha + 19/54)*beta - 187/108*alpha^2 + 38/27*alpha + 421/54
sage: conic[0](t_odd), conic[1](t_odd) # An odd point
(-187/108*alpha^2 + 38/27*alpha + 421/54, -1/108*alpha^2 + 8/27*alpha +
19/54)
sage: t0 = rational_point_conic(conic, t_odd)
sage: t0
5/6*beta + 11/2
sage: conic[0](t0), conic[1](t0)
(11/2, 5/6)
TODO:
Do not ask the parameter to give an element in QQ[alpha]. Needed if the original curve is of odd degree but alpha is of even degree.
Compute the trace of an element on a relative number field over its base field. Uses newton_sums.
This is intended to be used with the method newton_sums.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import rel_trace, newton_sums
sage: N = NumberField([x^3-2,x^2+4], 'a,b')
sage: nn = newton_sums(N)
sage: r = N.random_element()
sage: rel_trace(nn, r) - r.trace(N.base_ring())
0
Take a linear fraction U with coefficients on a number field and return a linear fraction V such that V = U(S), where S is a linear fraction in QQ.
It is expected that V has smaller coefficients that U.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: K.<x> = QQ[I][]
sage: I = QQ[I].gen()
sage: u = (125*x+I) / (125*x-I)
sage: simplify_unit(u)
(x + I)/(x - I)
sage: u = random_linear_fraction(K)
sage: v = simplify_unit(u)
sage: w = inverse_unit(u)(v)
sage: w in K.fraction_field()
True
Given a rational fraction a / b, return a canonical representative with monic denominator.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import simsim
sage: K.<x> = QQ[x]
sage: f = (3*x^2+1) /(3*x^3)
sage: f
(3*x^2 + 1)/(3*x^3)
sage: simsim(f)
(x^2 + 1/3)/x^3
sage: f = K.fraction_field().random_element(20,10**6,10)
sage: f - simsim(f)
0
sage: f = NumberField(x^4+2,'a')[x].fraction_field().random_element(10, 10,10)
sage: g = simsim(f)
sage: f - g
0
sage: g.denominator().leading_coefficient()
1
Compute the sum
sum( [( Nu[i] * alpha**i ) for i in range(len(Nu)) ] ) / De
Horner is not significant here alpha**i is immediate. It is expected to be used to deal with unit-parametrizations in terms of alpha.
INPUT:
Expected to be used if [X/De for X in Nu] is a unit parametrization.
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: N.<i> = NumberField(x^2+1)
sage: K.<t> = N[]
sage: Nu = [t^2 - i*t, -i*t^2 + 2*i*t]
sage: De = 2*t - i - 2
sage: sums_alpha(Nu, De, i)
t
A random example:
sage: N.<a> = NumberField(x^3-2, 'a')
sage: K.<t> = N[]
sage: u = random_linear_fraction(K)
sage: r = witness([u],name = 't')
sage: D = my_lcm([f.denominator() for f in r])
sage: Nu = [K(f*D) for f in r]
sage: sums_alpha(Nu, D, a)
t
Take a unit for a extension of degree 2 of the form N[sqrt(x)] and compute the parametrization of the associated conic. Units under composition are isomorphic to GL(2,F)
INPUT:
OUPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: N.<I> = QQ[I]
sage: K.<xul> = N[]
sage: u = random_linear_fraction(K)
sage: unit_to_hc_sqrt(u) == witness([u], name = xul,check = False)
True
Compute an automorphism of P^1 from the images of three points.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import unitfrom3points, is_com_unit
sage: K.<i1, i2, i3, j1, j2, j3> = QQ[]
sage: L.<t> = K.fraction_field()[]
sage: M = [[i1, j1], [i2, j2], [i3, j3]]
sage: u = unitfrom3points(M,t)
sage: is_com_unit(u)
True
sage: u(i1)
j1
sage: u(i2)
j2
sage: u(i3)
j3
sage: K.<t> = NumberField(x^3 - 2, 'a')[]
sage: u = (t+4) / (t + 1)
sage: i1, i2, i3 = randint(0, 10^6), randint(0, 10^6), randint(0, 10^6)
sage: M = [[i1, u(i1)], [i2, u(i2)], [i3, u(i3)] ]
sage: v = unitfrom3points(M,t)
sage: u - v
0
Compute the a witness variety of the parametrization of a rational curve. The witness variety is represented by its standard parametrization.
INPUT:
OUTPUT:
EXAMPLES:
sage: from sage.geometry.hypercircles.hypercircle import *
sage: N.<a> = NumberField(x^3 - 2)
sage: K.<t> = N[]
sage: u = random_linear_fraction(K)
sage: Phi = witness([u])
sage: Phi2 = witness(Phi)
sage: Phi == Phi2
True
sage: v = inverse_unit(u,t)
sage: L = QQ[t].fraction_field()
sage: [simsim(f(v)) in L for f in Phi]
[True, True, True]
sage: Phi[0] + a*Phi[1] + a^2*Phi[2]
t
Relative fields are supported (somehow, if it does not work, please report):
sage: KalphaX.<xl> = NumberField([x^2-2, x^2-3], 's2,s3')[]
sage: Kalpha = KalphaX.base_ring()
sage: K = Kalpha.base_ring()
sage: u = random_linear_fraction(KalphaX, 10, 1)
sage: Phi = witness([u], name = xl)
sage: change = is_hypercircle(Phi, inverse_unit(u,xl)(0))
sage: is_com_unit(change)
True
Non-primitive hypercircle and relative extensions. See Fields of Parametrization and Optimal Affine Reparametrization of Rational Curves:
sage: Kalpha.<a> = NumberField(x^4-4*x^3+12*x^2-16*x+8,'a')
sage: K = Kalpha.base_ring()
sage: KalphaX.<t> = Kalpha['t']
sage: Den = -22+26*a-9*a^2+2*a^3+4*t
sage: Phi = [ (-6 + 18*a -9*a^2 + 6*a^3 + (44 -52*a +18*a^2 -4*a^3)*t -4*t^2) / Den, (-12 -2*a +9*a^2 -a^3 +(4+4*a+4*a^2)*t +(12-16*a+6*a^2-2*a^3)*t^2) / Den]
sage: W = witness(Phi,name = 't')
sage: W[0]+a*W[1]+a^2*W[2]+a^3*W[3]
t
sage: [p.denominator().degree() for p in W]
[1, 1, 1, 1]
It is a conic, so we compute an intermediate field:
sage: point_infinity = [simsim(p(1/t)*t)(0) for p in W]
sage: gamma = point_infinity[0]
sage: gamma.minpoly().degree()
2
gamma is a primitive element of a field such that QQ(gamma) in QQ(a) defines a line as hypercircle, relativize is broken:
sage: gamma *= gamma.lift().denominator() #Bug in relativize
sage: Kgamma_alpha.<new_alpha, new_gamma> = Kalpha.relativize(gamma)
Recompute the original parametrization in Kgamma_alpha:
sage: to_Kgamma_alpha = Kgamma_alpha.structure()[1]
sage: Kgamma_alpha.register_coercion(to_Kgamma_alpha)
sage: NewPhi = [conjugate_pol(p, to_Kgamma_alpha, Kgamma_alpha[t]) for p in Phi]
sage: W2 = witness(NewPhi)
sage: W2
[t + (1/4*new_gamma + 1/2)*new_alpha, -1/4*new_gamma - 1/2]
The new hypercircle is a line. It can be parametrized by [t, -1/4*new_gamma - 1/2] We have found an algebraically optimal affine reparametrization.:
sage: change = t + W2[1] * new_alpha
sage: [simsim(p(change)) for p in NewPhi]
[(-t^2 + (1/6*new_gamma + 5/3)*t - 1/6*new_gamma - 1/6)/(t -
1/12*new_gamma - 5/6), ((-1/6*new_gamma + 1/3)*t^2 + (1/3*new_gamma -
5/3)*t - 1/6*new_gamma + 4/3)/(t - 1/12*new_gamma - 5/6)]
The reparametrization is over QQ[new_gamma].
Now an example not defined over K but over an strict intermediate field:
sage: x = var('x')
sage: N = NumberField(x^4-2,'a')
sage: a = N.gen()
sage: K = N['t']
sage: t = K.gen()
sage: Phi = [(t^2-a^2)/(t^2+(2+a^2)*t-a^2-3), ((4*a^4+3)*t-2)/(t^2+(2+a^2)*t-a^2-3)]
sage: u = random_linear_fraction(K)
sage: Phi = [f(u) for f in Phi]
sage: witness(Phi)
['N', 1, a^2]
sage: witness([t,a*t])
['N', 1, a, a^2, a^3]
sage: N = NumberField(x^9-2, 'a')
sage: a = N.gen()
sage: K = N['t']
sage: t = K.gen()
sage: Phi = [(t^2-a^3)/(t^2+(2+a^3+3*a^6)*t-a^6-3), ((4*a^3+3*a^6+1)*t-2-3*a^6)/(t^2+(2+a^3+3*a^6)*t-a^6-3)]
sage: u = random_linear_fraction(K)
sage: Phi = [f(u) for f in Phi]
sage: witness(Phi)
['N', 1, a^3, a^6]
sage: N = NumberField((x^7-1)/(x-1),'a')
sage: a = N.gen()
sage: K = N['t']
sage: t = K.gen()
sage: Phi = [t^3+1, (3*a+3*a**6+a^3+a^4)*t]
sage: u = random_linear_fraction(K)
sage: Phi = [f(u) for f in Phi]
sage: witness(Phi)
['N', 1, a^5 + a^2, a^4 + a^3]
TODO:
TESTS:
The following used to fail:
sage: N.<I> = QQ[I]
sage: K.<x> = N[]
sage: u = (x-I)/(x+I)
sage: witness([u])
[0, -I*t]
sage: var('x')
x
sage: N.<a> = NumberField(x^4-5, 'a')
sage: K.<x> = N[]
sage: u = (a**2*x-(a**2+1))/((4-a**2)*x+(a**2+3))
sage: witness([u])
[(1/2*t^2 - 1/10*a^2*t + 1/4)/(t - 1/10*a^2 - 1/4), 0, (1/10*a^2*t^2 -
1/20*a^2*t - 1/20*a^2)/(t - 1/10*a^2 - 1/4), 0]
For higher towers there were problems, they should be fixed now with Nbase_to_M morphism:
sage: var('x')
x
sage: KalphaX.<xl> = NumberField([x^2-2, x^2-3, x^2-5], 's2, s3, s5')[]
sage: Kalpha = KalphaX.base_ring()
sage: K = Kalpha.base_ring()
sage: u = random_linear_fraction(KalphaX, 10, 1)
sage: Phi = witness([u], name = xl)
sage: change = is_hypercircle(Phi, inverse_unit(u,xl)(0))
sage: is_com_unit(change)
True