Changeset 8245:81b1756e3f42


Ignore:
Timestamp:
01/26/08 19:38:25 (5 years ago)
Author:
John Cremona <john.cremona@…>
Branch:
default
Message:

trac #740 -- isomorphisms, automorphisms and twists in all characteristics

Location:
sage/schemes/elliptic_curves
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • sage/schemes/elliptic_curves/constructor.py

    r7000 r8245  
    11""" 
    22Elliptic curve constructor 
     3 
     4AUTHOR: 
     5   * William Stein (2005) -- Initial version 
     6   * John Cremona (Jan 2008) -- EllipticCurve(j) fixed for all cases 
    37""" 
    48 
     
    4852 
    4953        -- EllipticCurve(j): Return an elliptic curve with j-invariant 
    50            $j$. (Some mild hypothesis on char of base ring.) 
     54           $j$. 
    5155 
    5256    EXAMPLES: 
     
    116120 
    117121    if rings.is_RingElement(x): 
    118         # TODO: worry about char != 0!!!! 
    119         j = x 
     122        # Fixed for all characteristics and cases by John Cremona 
     123        j=x 
     124        F=j.parent().fraction_field() 
     125        char=F.characteristic() 
     126        if char==2: 
     127            if j==0: 
     128                return EllipticCurve(F, [ 0, 0, 1, 0, 0 ]) 
     129            else: 
     130                return EllipticCurve(F, [ 1, 0, 0, 0, 1/j ]) 
     131        if char==3: 
     132            if j==0: 
     133                return EllipticCurve(F, [ 0, 0, 0, 1, 0 ]) 
     134            else: 
     135                return EllipticCurve(F, [ 0, j, 0, 0, -j**2 ]) 
    120136        if j == 0: 
    121             return EllipticCurve(x.parent(), [ 0, 0, 1, 0, 0 ]) 
    122         elif j == 1728: 
    123             return EllipticCurve(x.parent(), [ 0, 0, 0, 1, 0 ]) 
    124         return EllipticCurve((j/1).parent(), [1,0,0,-36/(j - 1728), -1/(j - 1728)]) 
     137            return EllipticCurve(F, [ 0, 0, 0, 0, 1 ]) 
     138        if j == 1728: 
     139            return EllipticCurve(F, [ 0, 0, 0, 1, 0 ]) 
     140        k=j-1728 
     141        return EllipticCurve(F, [0,0,0,-3*j*k, -2*j*k**2]) 
    125142 
    126143    if not isinstance(x,list): 
  • sage/schemes/elliptic_curves/ell_generic.py

    r7935 r8245  
    1212    sage: latex(E) 
    1313    y^2  = x^3 + x + 3 
     14 
     15AUTHORS: 
     16   * William Stein (2005) -- Initial version 
     17   * Robert Bradshaw et al.... 
     18   * John Cremona (Jan 2008) -- isomorphisms, automorphisms and twists  
     19                                in all characteristics  
    1420""" 
    1521 
     
    2935#***************************************************************************** 
    3036 
    31  
    3237import math 
    3338 
     
    5661import constructor 
    5762import formal_group 
    58 import weierstrass_morphism 
     63import weierstrass_morphism as wm 
     64from constructor import EllipticCurve 
     65 
    5966 
    6067factor = arith.factor 
     
    95102    """ 
    96103    def __init__(self, ainvs, extra=None): 
     104        """ 
     105        Constructor from [a1,a2,a3,a4,a6] or [a4,a6] (see 
     106        constructor.py for more variants) 
     107 
     108        EXAMPLES: 
     109        sage: E = EllipticCurve([1,2,3,4,5]); E 
     110        Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field 
     111        sage: E = EllipticCurve(GF(7),[1,2,3,4,5]); E 
     112        Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Finite Field of size 7 
     113 
     114        Constructor from [a4,a6] sets a1=a2=a3=0: 
     115 
     116        sage: EllipticCurve([4,5]).ainvs() 
     117        [0, 0, 0, 4, 5] 
     118 
     119        Base need not be a field: 
     120 
     121        sage: EllipticCurve(IntegerModRing(91),[1,2,3,4,5]) 
     122        Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Ring of integers modulo 91 
     123 
     124        """ 
    97125        if extra != None:   # possibility of two arguments 
    98126            K, ainvs = ainvs, extra 
     
    682710     
    683711    def c4(self): 
     712        """ 
     713        EXAMPLES: 
     714            sage: E = EllipticCurve([0, -1, 1, -10, -20]) 
     715            sage: E.c4() 
     716            496 
     717        """ 
    684718        try: 
    685719            return self.__c_invariants[0] 
     
    689723 
    690724    def c6(self): 
     725        """ 
     726        EXAMPLES: 
     727            sage: E = EllipticCurve([0, -1, 1, -10, -20]) 
     728            sage: E.c6() 
     729            20008 
     730        """ 
    691731        try: 
    692732            return self.__c_invariants[1] 
     
    697737 
    698738    def base_extend(self, R): 
     739        """ 
     740        EXAMPLES: 
     741        sage: E=EllipticCurve(GF(5),[1,1]); E 
     742        Elliptic Curve defined by y^2  = x^3 + x +1 over Finite Field of size 5 
     743        sage: E1=E.base_extend(GF(125,'a')); E1 
     744        Elliptic Curve defined by y^2  = x^3 + x +1 over Finite Field in a of size 5^3 
     745        """ 
    699746        return constructor.EllipticCurve(R, [R(a) for a in self.a_invariants()]) 
    700747 
     
    772819        return self.gens()[i] 
    773820         
     821 
     822    # Twists: rewritten by John Cremona as follows: 
     823    # 
     824    # Quadratic twist allowed except when char=2, j=0 
     825    # Quartic twist allowed only if j=1728!=0 (so char!=2,3) 
     826    # Sextic  twist allowed only if j=0!=1728 (so char!=2,3) 
     827    # 
     828    # More complicated twists exist in theory for char=2,3 and 
     829    # j=0=1728, but I have never worked them out or seen them used! 
     830    #                                              
     831 
    774832    def quadratic_twist(self, D): 
    775833        """ 
    776834        Return the quadratic twist of this curve by D. 
     835 
     836        In characteristic!=2, D must be nonzero, and the twist is 
     837        isomorphic to self after adjoining sqrt(D) to the base 
     838 
     839        In characteristic==2, D is arbitrary, and the twist is 
     840        isomorphic to self after adjoining a root of x^2+x+D to the 
     841        base 
     842 
     843        In characteristics 2 when j==0 this is not implemented (the 
     844        twists are more complicated than quadratic!) 
    777845 
    778846        EXAMPLES: 
    779847            sage: E = EllipticCurve([GF(1103)(1), 0, 0, 107, 340]); E 
    780848            Elliptic Curve defined by y^2 + x*y  = x^3 + 107*x + 340 over Finite Field of size 1103 
    781             sage: E.quadratic_twist(-1) 
    782             Elliptic Curve defined by y^2  = x^3 + 770*x + 437 over Finite Field of size 1103 
    783         """ 
    784         a = self.ainvs() 
    785         ap = [0, 0, 0] 
    786         ap.append(-27*D**2*a[0]**4 - 216*D**2*a[0]**2*a[1] + 648*D**2*a[0]*a[2] - 432*D**2*a[1]**2 + 1296*D**2*a[3]) 
    787         ap.append(54*D**3*a[0]**6 + 648*D**3*a[0]**4*a[1] 
    788                   - 1944*D**3*a[0]**3*a[2] + 2592*D**3*a[0]**2*a[1]**2 
    789                   - 3888*D**3*a[0]**2*a[3] - 7776*D**3*a[0]*a[1]*a[2] 
    790                   + 3456*D**3*a[1]**3 - 15552*D**3*a[1]*a[3] 
    791                   + 11664*D**3*a[2]**2 + 46656*D**3*a[4]) 
    792         import constructor 
    793         return constructor.EllipticCurve(self.base_ring(), ap) 
     849            sage: F=E.quadratic_twist(-1); F 
     850            Elliptic Curve defined by y^2  = x^3 + 1102*x^2 + 609*x + 300 over Finite Field of size 1103 
     851            sage: E.is_isomorphic(F) 
     852            False 
     853            sage: E.is_isomorphic(F,GF(1103^2,'a')) 
     854            True 
     855 
     856            A characteristic 2 example: 
     857 
     858            sage: E=EllipticCurve(GF(2),[1,0,1,1,1]) 
     859            sage: E1=E.quadratic_twist(1) 
     860            sage: E.is_isomorphic(E1) 
     861            False 
     862            sage: E.is_isomorphic(E1,GF(4,'a')) 
     863            True 
     864 
     865        """ 
     866        K=self.base_ring() 
     867        char=K.characteristic() 
     868        
     869        if char!=2:  
     870            b2,b4,b6,b8=self.b_invariants() 
     871            # E is isomorphic to  [0,b2,0,8*b4,16*b6] 
     872            return EllipticCurve(K,[0,b2*D,0,8*b4*D**2,16*b6*D**3]) 
     873 
     874        # now char==2 
     875        if self.j_invariant() !=0: # iff a1!=0 
     876            a1,a2,a3,a4,a6=self.ainvs() 
     877            E0=self.change_weierstrass_model(a1,a3/a1,0,(a1**2*a4+a3**2)/a1**3) 
     878            # which has the form = [1,A2,0,0,A6] 
     879            assert E0.a1()==K(1) 
     880            assert E0.a3()==K(0) 
     881            assert E0.a4()==K(0) 
     882            return EllipticCurve(K,[1,E0.a2()+D,0,0,E0.a6()]) 
     883        else: 
     884            raise ValueError, "Quadratic twist not implemented in char 2 when j=0" 
     885 
     886    def quartic_twist(self, D): 
     887        """ 
     888        Return the quartic twist of this curve by D. 
     889 
     890        The characteristic must not be 2 or 3 and the j-invariant must be 1728 
     891 
     892        EXAMPLES: 
     893        sage: E=EllipticCurve(GF(13)(1728)); E 
     894        Elliptic Curve defined by y^2  = x^3 + x over Finite Field of size 13 
     895        sage: E1=E.quartic_twist(2); E1 
     896        Elliptic Curve defined by y^2  = x^3 + 5*x over Finite Field of size 13 
     897        sage: E.is_isomorphic(E1) 
     898        False 
     899        sage: E.is_isomorphic(E1,GF(13^2,'a')) 
     900        False 
     901        sage: E.is_isomorphic(E1,GF(13^4,'a')) 
     902        True 
     903        """ 
     904        K=self.base_ring() 
     905        char=K.characteristic() 
     906        
     907        if char==2 or char==3:  
     908            raise ValueError, "Quartic twist not defined in chars 2,3" 
     909 
     910        if self.j_invariant() !=K(1728): 
     911            raise ValueError, "Quartic twist not defined when j!=1728" 
     912 
     913        c4,c6=self.c_invariants() 
     914        # E is isomorphic to  [0,0,0,-27*c4,0] 
     915        assert c6==0 
     916        return EllipticCurve(K,[0,0,0,-27*c4*D,0]) 
     917 
     918    def sextic_twist(self, D): 
     919        """ 
     920        Return the sextic twist of this curve by D. 
     921 
     922        The characteristic must not be 2 or 3 and the j-invariant must be 0 
     923 
     924        EXAMPLES: 
     925        sage: E=EllipticCurve(GF(13)(0)); E 
     926        Elliptic Curve defined by y^2  = x^3 +1 over Finite Field of size 13 
     927        sage: E1=E.sextic_twist(2); E1 
     928        Elliptic Curve defined by y^2  = x^3 + 11 over Finite Field of size 13 
     929        sage: E.is_isomorphic(E1) 
     930        False 
     931        sage: E.is_isomorphic(E1,GF(13^2,'a')) 
     932        False 
     933        sage: E.is_isomorphic(E1,GF(13^4,'a')) 
     934        False 
     935        sage: E.is_isomorphic(E1,GF(13^6,'a')) 
     936        True 
     937        """ 
     938        K=self.base_ring() 
     939        char=K.characteristic() 
     940        
     941        if char==2 or char==3:  
     942            raise ValueError, "Sextic twist not defined in chars 2,3" 
     943 
     944        if self.j_invariant() !=K(0): 
     945            raise ValueError, "Sextic twist not defined when j!=1728" 
     946 
     947        c4,c6=self.c_invariants() 
     948        # E is isomorphic to  [0,0,0,0,-54*c6] 
     949        assert c4==0 
     950        return EllipticCurve(K,[0,0,0,0,-54*c6*D]) 
    794951 
    795952    def rst_transform(self, r, s, t): 
    796953        """ 
    797         Transforms the elliptic curve using the unimodular (u=1) transform with standard parameters [r,s,t]. 
     954        Transforms the elliptic curve using the unimodular (u=1) 
     955        transform with standard parameters [r,s,t].  This is just a 
     956        special case of change_weierstrass_model(). 
    798957         
    799958        Returns the transformed curve. 
    800         """ 
    801         ##Ported from John Cremona's code implementing Tate's algorithm. 
    802         (a1, a2, a3, a4, a6) = self.a_invariants() 
    803         a6 += r*(a4 + r*(a2 + r)) - t*(a3 + r*a1 + t) 
    804         a4 += -s*a3 + 2*r*a2 - (t + r*s)*a1 + 3*r*r - 2*s*t 
    805         a3 += r*a1 + 2*t 
    806         a2 += -s*a1 + 3*r - s*s 
    807         a1 += 2*s 
    808         return constructor.EllipticCurve([a1, a2, a3, a4, a6]) 
     959        EXAMPLES: 
     960        sage: R.<r,s,t>=QQ[] 
     961        sage: E=EllipticCurve([1,2,3,4,5]) 
     962        sage: E.rst_transform(r,s,t) 
     963        Elliptic Curve defined by y^2 + (2*s+1)*x*y + (r+2*t+3)*y = x^3 + (-s^2+3*r-s+2)*x^2 + (3*r^2-r*s-2*s*t+4*r-3*s-t+4)*x + (r^3+2*r^2-r*t-t^2+4*r-3*t+5) over Multivariate Polynomial Ring in r, s, t over Rational Field 
     964 
     965        """ 
     966        return self.change_weierstrass_model(1,r,s,t) 
    809967     
    810968    def scale_curve(self, u): 
    811969        """ 
    812970        Transforms the elliptic curve using scale factor $u$, 
    813         i.e. multiplies $c_i$ by $u^i$. 
     971        i.e. multiplies $c_i$ by $u^i$.  This is  another 
     972        special case of change_weierstrass_model(). 
    814973         
    815974        Returns the transformed curve. 
    816         """ 
    817         ##Ported from John Cremona's code implementing Tate's algorithm. 
    818         (a1,a2,a3,a4,a6) = self.a_invariants() 
    819         a1 *= u 
    820         a2 *= u^2 
    821         a3 *= u^3 
    822         a4 *= u^4 
    823         a6 *= u^6 
    824         return constructor.EllipticCurve([a1, a2, a3, a4, a6]) 
     975 
     976        EXAMPLES: 
     977        sage: K=Frac(PolynomialRing(QQ,'u')) 
     978        sage: u=K.gen() 
     979        sage: E=EllipticCurve([1,2,3,4,5])        
     980        sage: E.scale_curve(u) 
     981        Elliptic Curve defined by y^2 + u*x*y + 3*u^3*y = x^3 + 2*u^2*x^2 + 4*u^4*x + 5*u^6 over Fraction Field of Univariate Polynomial Ring in u over Rational Field 
     982 
     983        """ 
     984        if isinstance(u, (int,long)): 
     985            u=self.base_ring()(u)       # because otherwise 1/u would round! 
     986        return self.change_weierstrass_model(1/u,0,0,0) 
    825987 
    826988    def discriminant(self): 
     
    12981460              From: Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field 
    12991461              To:   Abelian group of points on Elliptic Curve defined by y^2  = x^3 - x + 1/4 over Rational Field 
     1462              Via:  (u,r,s,t) = (1, 0, 0, -1/2) 
    13001463            sage: P = E(0,-1,1) 
    13011464            sage: w(P) 
     
    13141477                From: Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field 
    13151478                To:   Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 + (-1)*x over Number Field in a with defining polynomial x^3 - 7 
    1316         """ 
    1317         return weierstrass_morphism.WeierstrassIsomorphism(self, other) 
    1318      
    1319     def is_isomorphic(self, other): 
    1320         """ 
    1321         Returns whether or not self is isomorphic to other, i.e. they define the same curve over  
    1322         the same basering.  
     1479                Via:  (u,r,s,t) = (1, 0, 0, 0) 
     1480        """ 
     1481        return wm.WeierstrassIsomorphism(self, None, other) 
     1482 
     1483    def automorphisms(self, field=None): 
     1484        """ 
     1485        Return the set of isomorphisms from self to itself (as a list). 
     1486         
     1487        EXAMPLES:  
     1488        sage: E = EllipticCurve(QQ(0)) # a curve with j=0 over QQ 
     1489        sage: E.automorphisms();        
     1490        [Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Rational Field 
     1491        Via:  (u,r,s,t) = (-1, 0, 0, 0), Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Rational Field 
     1492        Via:  (u,r,s,t) = (1, 0, 0, 0)] 
     1493 
     1494        We can also find automorphisms defined over extension fields:  
     1495        sage: K.<a> = NumberField(x^2+3) # adjoin roots of unity 
     1496        sage: E.automorphisms(K) 
     1497        [Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Number Field in a with defining polynomial x^2 + 3 
     1498        Via:  (u,r,s,t) = (1, 0, 0, 0), 
     1499        ... 
     1500        Generic endomorphism of Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Number Field in a with defining polynomial x^2 + 3 
     1501        Via:  (u,r,s,t) = (-1/2*a - 1/2, 0, 0, 0)] 
     1502 
     1503        sage: [ len(EllipticCurve(GF(q,'a')(0)).automorphisms()) for q in [2,4,3,9,5,25,7,49]] 
     1504        [2, 24, 2, 12, 2, 6, 6, 6] 
     1505        """ 
     1506        if field==None: 
     1507            return [wm.WeierstrassIsomorphism(self, urst, self) 
     1508                    for urst in wm.isomorphisms(self,self)] 
     1509        E=self.change_ring(field) 
     1510        return [wm.WeierstrassIsomorphism(E, urst, E) 
     1511                for urst in wm.isomorphisms(E,E)] 
     1512         
     1513    def isomorphisms(self, other, field=None): 
     1514        """ 
     1515        Return the set of isomorphisms from self to other (as a list). 
     1516         
     1517        EXAMPLES:  
     1518        sage: E = EllipticCurve(QQ(0)) # a curve with j=0 over QQ 
     1519        sage: F = EllipticCurve('36a1') # should be the same one 
     1520        sage: E.isomorphisms(F);        
     1521        [Generic morphism: 
     1522        From: Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Rational Field 
     1523        To:   Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Rational Field 
     1524        Via:  (u,r,s,t) = (-1, 0, 0, 0), Generic morphism: 
     1525        From: Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Rational Field 
     1526        To:   Abelian group of points on Elliptic Curve defined by y^2  = x^3 +1 over Rational Field 
     1527        Via:  (u,r,s,t) = (1, 0, 0, 0)] 
     1528 
     1529        We can also find istomorphisms defined over extension fields:  
     1530        sage: E=EllipticCurve(GF(7),[0,0,0,1,1]) 
     1531        sage: F=EllipticCurve(GF(7),[0,0,0,1,-1]) 
     1532        sage: E.isomorphisms(F) 
     1533        [] 
     1534        sage: E.isomorphisms(F,GF(49,'a')) 
     1535        [Generic morphism: 
     1536        From: Abelian group of points on Elliptic Curve defined by y^2  = x^3 + x +1 over Finite Field in a of size 7^2 
     1537        To:   Abelian group of points on Elliptic Curve defined by y^2  = x^3 + x + 6 over Finite Field in a of size 7^2 
     1538        Via:  (u,r,s,t) = (a + 3, 0, 0, 0), Generic morphism: 
     1539        From: Abelian group of points on Elliptic Curve defined by y^2  = x^3 + x +1 over Finite Field in a of size 7^2 
     1540        To:   Abelian group of points on Elliptic Curve defined by y^2  = x^3 + x + 6 over Finite Field in a of size 7^2 
     1541        Via:  (u,r,s,t) = (6*a + 4, 0, 0, 0)] 
     1542        """ 
     1543        if field==None: 
     1544            return [wm.WeierstrassIsomorphism(self, urst, other) 
     1545                    for urst in wm.isomorphisms(self,other)] 
     1546        E=self.change_ring(field) 
     1547        F=other.change_ring(field) 
     1548        return [wm.WeierstrassIsomorphism(E, urst, F) 
     1549                for urst in wm.isomorphisms(E,F)] 
     1550         
     1551    def is_isomorphic(self, other, field=None): 
     1552        """ 
     1553        Returns whether or not self is isomorphic to other, i.e. they 
     1554        define the same curve over the same basering. 
     1555 
     1556        If field!=None then both curves must base_extend-able to it 
     1557        and the isomorphism is then checked over that field 
    13231558         
    13241559        EXAMPLES:  
    13251560            sage: E = EllipticCurve('389a') 
    13261561            sage: F = E.change_weierstrass_model([2,3,4,5]); F 
    1327             Elliptic Curve defined by y^2 - 8*x*y + 22*y = x^3 - 21*x^2 + 59*x over Rational Field 
     1562            Elliptic Curve defined by y^2 + 4*x*y + 11/8*y = x^3 - 3/2*x^2 - 13/16*x over Rational Field 
    13281563            sage: E.is_isomorphic(F) 
    13291564            True 
     
    13331568        if not is_EllipticCurve(other): 
    13341569            return False 
    1335         elif self.base_ring() != other.base_ring(): 
    1336             return False 
     1570        if field==None: 
     1571            if self.base_ring() != other.base_ring(): 
     1572                return False 
     1573            elif self.j_invariant() != other.j_invariant(): ## easy check 
     1574                return False 
     1575            else: 
     1576                return wm.isomorphisms(self,other,True) != None 
    13371577        else: 
    1338             try: 
    1339                 phi = self.isomorphism_to(other) 
    1340                 return True 
    1341             except ValueError: 
     1578            E=self.base_extend(field) 
     1579            F=other.base_extend(field) 
     1580            if E.j_invariant() != F.j_invariant(): ## easy check 
    13421581                return False 
     1582            else: 
     1583                return wm.isomorphisms(E,other,F) != None 
    13431584         
    13441585    def change_weierstrass_model(self, *urst): 
     
    13491590        EXAMPLES:  
    13501591            sage: E = EllipticCurve('15a') 
    1351             sage: F = E.change_weierstrass_model([1/2,0,0,0]); F 
    1352             Elliptic Curve defined by y^2 + 1/2*x*y + 1/8*y = x^3 + 1/4*x^2 - 5/8*x - 5/32 over Rational Field 
    1353             sage: F = E.change_weierstrass_model([7,2,1/3,5]); F 
    1354             Elliptic Curve defined by y^2 + 19/3*x*y + 961/3*y = x^3 + 407/9*x^2 - 216512/9*x - 10141876/9 over Rational Field 
    1355             sage: E.is_isomorphic(F) 
     1592            sage: F1 = E.change_weierstrass_model([1/2,0,0,0]); F1 
     1593            Elliptic Curve defined by y^2 + 2*x*y + 8*y = x^3 + 4*x^2 - 160*x - 640 over Rational Field 
     1594            sage: F2 = E.change_weierstrass_model([7,2,1/3,5]); F2 
     1595            Elliptic Curve defined by y^2 + 5/21*x*y + 13/343*y = x^3 + 59/441*x^2 - 10/7203*x - 58/117649 over Rational Field 
     1596            sage: F1.is_isomorphic(F2) 
    13561597            True 
    13571598        """ 
    13581599        if isinstance(urst[0], (tuple, list)): 
    13591600            urst = urst[0] 
    1360         u, r, s, t = urst 
    1361         a1, a2, a3, a4, a6 = self.ainvs() 
    1362         b1 = a1*u - 2*s 
    1363         b3 = a3*u**3 - a1*r*u - 2*t + 2*r*s 
    1364         b2 = a2*u**2 + a1*s*u - s**2 - 3*r 
    1365         b4 = a4*u**4 + a3*s*u**3 - 2*a2*r*u**2 + a1*t*u - 2*a1*r*s*u - 2*s*t + 2*r*s**2 + 3*r**2 
    1366         b6 = a6*u**6 - a4*r*u**4 + a3*t*u**3 - a3*r*s*u**3 + a2*r**2*u**2 - a1*r*t*u + a1*r**2*s*u - t**2 + 2*r*s*t - r**2*s**2 - r**3 
    1367         return constructor.EllipticCurve(self.base_ring(), [b1,b2,b3,b4,b6]) 
     1601        return constructor.EllipticCurve((wm.baseWI(*urst))(self.ainvs())) 
    13681602         
    13691603    def weierstrass_model(self): 
     
    13891623        """ 
    13901624        import constructor 
     1625        K = self.base_ring() 
     1626        if K.characteristic() == 2 or K.characteristic() == 3: 
     1627            raise ValueError, "weierstrass_model(): no short model for %s (characteristic is %s)"%(self,K.characteristic()) 
     1628 
    13911629        c4, c6 = self.c_invariants() 
    13921630        return constructor.EllipticCurve([-c4/(2**4*3), -c6/(2**5*3**3)]) 
    13931631         
     1632### The following functions should not be in ell_generic.py but in ell_rational_field.py!  John Cremona 
    13941633    def integral_weierstrass_model(self): 
    13951634        raise NotImplementedError 
     
    14011640 
    14021641        EXAMPLES: 
    1403             sage: E = EllipticCurve([1,2,3,4,5]); E 
    1404             Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field 
     1642            sage: E = EllipticCurve([1/2,2/3,3/5,4/7,5/11]); E 
     1643            Elliptic Curve defined by y^2 + 1/2*x*y + 3/5*y = x^3 + 2/3*x^2 + 4/7*x + 5/11 over Rational Field 
    14051644            sage: E.integral_model() 
    1406             (Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field, Generic morphism: 
    1407               From: Abelian group of points on Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field 
    1408               To:   Abelian group of points on Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field)         
     1645            (Elliptic Curve defined by y^2 + 1155*x*y + 7395834600*y = x^3 + 3557400*x^2 + 16270836120000*x + 69063597765855000000 over Rational Field, 
     1646            Generic morphism: 
     1647             From: Abelian group of points on Elliptic Curve defined by y^2 + 1/2*x*y + 3/5*y = x^3 + 2/3*x^2 + 4/7*x + 5/11 over Rational Field 
     1648             To:   Abelian group of points on Elliptic Curve defined by y^2 + 1155*x*y + 7395834600*y = x^3 + 3557400*x^2 + 16270836120000*x + 69063597765855000000 over Rational Field 
     1649             Via:  (u,r,s,t) = (1/2310, 0, 0, 0)) 
     1650            sage: E.integral_model()[0].minimal_model() 
     1651            Elliptic Curve defined by y^2 + x*y + y = x^3 - x^2 + 15495546786070*x + 60459278934205438697 over Rational Field 
    14091652        """ 
    14101653        denom = lcm([a.denominator() for a in self.ainvs()]) 
  • sage/schemes/elliptic_curves/weierstrass_morphism.py

    r7687 r8245  
     1""" 
     2Isomorphisms between Weierstrass models of elliptic curves 
     3 
     4AUTHORS: 
     5   * Robert Bradshaw (2007): initial version 
     6   * John Cremona (Jan 2008): isomorphisms, automorphisms and twists  
     7                              in all characteristics  
     8""" 
    19 
    210#***************************************************************************** 
     
    2028from sage.categories.homset import Hom 
    2129 
    22 class WeierstrassIsomorphism(Morphism): 
     30class baseWI: 
     31    r""" This class implements the basic arithmetic of isomorphisms 
     32    between Weierstrass models of elliptic curves.  These are 
     33    specified by lists of the form [u,r,s,t] (with u!=0) which 
     34    specifies a transformation $(x,y) \mapsto (x',y')$ where 
     35         
     36            $$ (x,y) = (u^2x'+r , u^3y' + su^2x' + t) $$ 
     37    """ 
     38 
     39    def __init__(self, u=1, r=0, s=0, t=0): 
     40        r"""  
     41        Constructor: check for valid parameters (defaults to identity) 
     42 
     43        EXAMPLES: 
     44        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     45        sage: baseWI() 
     46        (1, 0, 0, 0) 
     47        sage: baseWI(2,3,4,5) 
     48        (2, 3, 4, 5) 
     49        sage: w=var('u r s t'); baseWI(u,r,s,t) 
     50        (u, r, s, t) 
     51        """     
     52        if u==0: 
     53            raise ValueError, "u!=0 required for baseWI" 
     54        self.u=u; self.r=r; self.s=s; self.t=t 
     55 
     56    def __cmp__(self, other): 
     57        """ 
     58        The ordering is just lexicographic on the tuple (u,r,s,t)  
     59 
     60        NB In a list of automprhisms there's no guarantee that the 
     61        identity will be first! 
     62 
     63        EXAMPLE: 
     64 
     65        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     66        sage: baseWI(1,2,3,4)==baseWI(1,2,3,4) 
     67        True 
     68        sage: baseWI(1,2,3,4)<baseWI(1,2,3,5) 
     69        True 
     70        sage: baseWI(1,2,3,4)>baseWI(1,2,3,4) 
     71        False 
     72 
     73        It will never return equaity if other is of another type: 
     74        sage: baseWI() == 1 
     75        False 
     76 
     77        """ 
     78        if not isinstance(other, baseWI): 
     79            return cmp(type(self), type(other)) 
     80        return cmp(self.tuple(), other.tuple()) 
     81 
     82    def tuple(self): 
     83        r""" 
     84        Returns the parameters u,r,s,t as a tuple 
     85 
     86        EXAMPLES: 
     87        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     88        sage: u,r,s,t=baseWI(2,3,4,5).tuple() 
     89        sage: w=baseWI(2,3,4,5) 
     90        sage: u,r,s,t=w.tuple() 
     91        sage: u 
     92        2 
     93 
     94        """ 
     95        return (self.u,self.r,self.s,self.t) 
     96 
     97    def __mul__(self, other): 
     98        r"""  
     99        Composition of isomorphisms 
     100 
     101        EXAMPLES: 
     102        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     103        sage: baseWI(1,2,3,4)*baseWI(5,6,7,8) 
     104        (5, 56, 22, 858) 
     105        sage: baseWI()*baseWI(1,2,3,4)*baseWI() 
     106        (1, 2, 3, 4) 
     107        """     
     108        u1,r1,s1,t1=other.tuple() 
     109        u2,r2,s2,t2=self.tuple() 
     110        return baseWI(u1*u2,(u1**2)*r2+r1,u1*s2+s1,(u1**3)*t2+s1*(u1**2)*r2+t1) 
     111 
     112    def __invert__(self): 
     113        r"""  
     114        Inverse of an isomporphism 
     115 
     116        EXAMPLES: 
     117        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     118        sage: w=baseWI(2,3,4,5) 
     119        sage: ~w 
     120        (1/2, -3/4, -2, 7/8) 
     121        sage: w*~w 
     122        (1, 0, 0, 0) 
     123        sage: ~w*w 
     124        (1, 0, 0, 0) 
     125        sage: var('u r s t'); w=baseWI(u,r,s,t) 
     126        (u, r, s, t) 
     127        sage: ~w 
     128        (1/u, -r/u^2, -s/u, (r*s - t)/u^3) 
     129        sage: ~w*w 
     130        (1, 0, 0, 0) 
     131        """     
     132        u,r,s,t=self.tuple() 
     133        return baseWI(1/u,-r/(u**2),-s/u,(r*s-t)/(u**3)) 
     134 
     135    def __repr__(self): 
     136        r"""  
     137        String representation  
     138 
     139        EXAMPLES: 
     140        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     141        sage: baseWI(2,3,4,5) 
     142        (2, 3, 4, 5) 
     143        """     
     144        return self.tuple().__repr__() 
     145 
     146    def is_identity(self): 
     147        r"""  
     148        Tests for identity 
     149 
     150        EXAMPLES: 
     151        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     152        sage: w=baseWI(); w.is_identity() 
     153        True 
     154        sage: w=baseWI(2,3,4,5); w.is_identity() 
     155        False 
     156        """     
     157        return self.tuple()==(1,0,0,0) 
     158 
     159    def __call__(self, EorP): 
     160        r""" Base application of isomorphisms to curves and points: a baseWI w 
     161        may be applied to a list [a1,a2,a3,a4,a6] representing the 
     162        a-invariants of an elliptic curve E, returning the 
     163        a-invariants of w(E); or to P=[x,y] or P=[x,y,z] represeting a 
     164        point in A^2 or P^2, returning the transformed point. 
     165 
     166        EXAMPLES: 
     167        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     168        sage: E=EllipticCurve([0,0,1,-7,6]) 
     169        sage: w=baseWI(2,3,4,5); 
     170        sage: w(E.ainvs()) 
     171        [4, -7/4, 11/8, -3/2, -9/32] 
     172        sage: P=E(-2,3) 
     173        sage: w(P.xy()) 
     174        [-5/4, 9/4] 
     175        sage: EllipticCurve(w(E.ainvs()))(w(P.xy())) 
     176        (-5/4 : 9/4 : 1) 
     177        """     
     178        u,r,s,t=self.tuple() 
     179        if len(EorP)==5: 
     180           a1,a2,a3,a4,a6=EorP 
     181           a6 += r*(a4 + r*(a2 + r)) - t*(a3 + r*a1 + t); 
     182           a4 += -s*a3 + 2*r*a2 - (t + r*s)*a1 + 3*r*r - 2*s*t; 
     183           a3 += r*a1 +t+t; 
     184           a2 += -s*a1 + 3*r - s*s; 
     185           a1 += 2*s; 
     186           return [a1/u,a2/u**2,a3/u**3,a4/u**4,a6/u**6] 
     187        if len(EorP)==2: 
     188           x,y=EorP 
     189           x-=r 
     190           y-=(s*x+t) 
     191           return [x/u**2,y/u**3] 
     192        if len(EorP)==3: 
     193           x,y,z=EorP 
     194           x-=r*z 
     195           y-=(s*x+t*z) 
     196           return [x/u**2,y/u**3,z] 
     197        raise ValueError, "baseWI(a) only for a=(x,y), (x:y:z) or (a1,a2,a3,a4,a6)" 
     198            
     199def isomorphisms(E,F,JustOne=False): 
     200    r""" function to find one or all isomorphisms between two elliptic 
     201    curves, as simple (u,r,s,t) tuples. 
     202 
     203    if JustOne==False it returns a list of these, which is empty 
     204    iff the curves are not isomorphic.  The length of the list is 
     205    very often 0 or 2. 
    23206     
    24     def __init__(self, E, F): 
    25         r""" 
    26         Given two Weierstrass models $E$ and $F$ of the same elliptic curve,  
    27         construct an isomorphism from $E$ to $F$.  
     207    if JustOne==True it returns one isomorphism, or None if the 
     208    curves are not isomorphic. 
     209 
     210    This function is not intended for users, who should use the 
     211    interface provided by ell_generic. 
     212 
     213    EXAMPLES: 
     214    sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     215    sage: isomorphisms(EllipticCurve(0),EllipticCurve('36a1')) 
     216    [(-1, 0, 0, 0), (1, 0, 0, 0)] 
     217    sage: isomorphisms(EllipticCurve(0),EllipticCurve('36a1'),JustOne=True) 
     218    (1, 0, 0, 0) 
     219    sage: isomorphisms(EllipticCurve(0),EllipticCurve('36a2')) 
     220    [] 
     221    sage: isomorphisms(EllipticCurve(0),EllipticCurve('36a2'),JustOne=True) 
     222    """     
     223    from ell_generic import is_EllipticCurve 
     224    if not is_EllipticCurve(E) or not is_EllipticCurve(F): 
     225        raise ValueError, "arguments are not elliptic curves" 
     226    K = E.base_ring() 
     227#   if not K == F.base_ring(): return [] 
     228    j=E.j_invariant() 
     229    if  j != F.j_invariant(): 
     230        if JustOne: return None 
     231        return [] 
     232     
     233    from sage.rings.all import PolynomialRing 
     234    x=PolynomialRing(K,'x').gen() 
     235 
     236    a1E, a2E, a3E, a4E, a6E = E.ainvs()             
     237    a1F, a2F, a3F, a4F, a6F = F.ainvs() 
     238 
     239    char=K.characteristic() 
    28240         
    29         Alternatively, the second argument may be a vector of the form [u,r,s,t] 
    30         which specifies a transformation 
     241    if char==2: 
     242        if j==0: 
     243            ulist=[u[0] for u in (x**3-(a3E/a3F)).roots()] 
     244            ans=[] 
     245            for u in ulist: 
     246                slist=[s[0] for s in 
     247                       (x**4+a3E*x+(a2F**2+a4F)*u**4+a2E**2+a4E).roots()] 
     248                for s in slist: 
     249                    r=s**2+a2E+a2F*u**2 
     250                    tlist=[t[0] for t in 
     251                    (x**2 + a3E*x + r**3 + a2E*r**2 + a4E*r + a6E + a6F*u**6).roots()] 
     252                    for t in tlist: 
     253                        if JustOne: return (u,r,s,t) 
     254                        ans.append((u,r,s,t)) 
     255            if JustOne: return None 
     256            ans.sort() 
     257            return ans 
     258        else: 
     259            ans=[] 
     260            u=a1E/a1F 
     261            r=(a3E+a3F*u**3)/a1E 
     262            slist=[s[0] for s in (x**2+a1E*x+(r+a2E+a2F*u**2)).roots()] 
     263            for s in slist: 
     264                t = (a4E+a4F*u**4 + s*a3E + r*s*a1E + r**2) 
     265                if JustOne: return (u,r,s,t) 
     266                ans.append((u,r,s,t)) 
     267            if JustOne: return None 
     268            ans.sort() 
     269            return ans 
     270 
     271    b2E, b4E, b6E, b8E      = E.b_invariants() 
     272    b2F, b4F, b6F, b8F      = F.b_invariants() 
     273 
     274    if char==3: 
     275        if j==0: 
     276            ulist=[u[0] for u in (x**4-(b4E/b4F)).roots()] 
     277            ans=[] 
     278            for u in ulist: 
     279                s=a1E-a1F*u 
     280                t=a3E-a3F*u**3 
     281                rlist=[r[0] for r in (x**3-b4E*x+(b6E-b6F*u**6)).roots()] 
     282                for r in rlist: 
     283                    if JustOne: return (u,r,s,t+r*a1E) 
     284                    ans.append((u,r,s,t+r*a1E)) 
     285            if JustOne: return None 
     286            ans.sort() 
     287            return ans 
     288        else: 
     289            ulist=[u[0] for u in (x**2-(b2E/b2F)).roots()] 
     290            ans=[] 
     291            for u in ulist: 
     292                r = (b4F*u**4 -b4E)/b2E 
     293                s = (a1E-a1F*u) 
     294                t = (a3E-a3F*u**3 + a1E*r) 
     295                if JustOne: return (u,r,s,t) 
     296                ans.append((u,r,s,t)) 
     297            if JustOne: return None 
     298            ans.sort() 
     299            return ans 
     300             
     301# now char!=2,3: 
     302    c4E,c6E = E.c_invariants() 
     303    c4F,c6F = F.c_invariants() 
     304 
     305    if j==0: 
     306        m,um = 6,c6E/c6F 
     307    elif j==1728: 
     308        m,um=4,c4E/c4F 
     309    else: 
     310        m,um=2,(c6E*c4F)/(c6F*c4E) 
     311    ulist=[u[0] for u in (x**m-um).roots()] 
     312    ans=[] 
     313    for u in ulist: 
     314        s = (a1F*u - a1E)/2 
     315        r = (a2F*u**2 + a1E*s + s**2 - a2E)/3 
     316        t = (a3F*u**3 - a1E*r - a3E)/2 
     317        if JustOne: return (u,r,s,t) 
     318        ans.append((u,r,s,t)) 
     319    if JustOne: return None 
     320    ans.sort() 
     321    return ans 
    31322         
    32             $$ (x,y) \mapsto (x',y') = (u^2*x+r , u^3*y + s*u^2*x' + t) $$ 
     323class WeierstrassIsomorphism(baseWI,Morphism): 
     324     
     325    def __init__(self, E=None, urst=None, F=None): 
     326        r"""  
     327        Given two Weierstrass models $E$ and $F$ of elliptic curves, and 
     328        a transformation urst from E to F, construct an isomorphism 
     329        from $E$ to $F$.  An exception is raised if urst(E)!=F. 
     330        At most one of E, F, urst can be None.  If F==None then it 
     331        constructed as urst(E).  If E==None then it is constructed as 
     332        urst^-1(F).  (This is needed surprisingly often).  If 
     333        urst==None then an isomorphism from E to F is constructed, and 
     334        an exception is raised if they are not isomorphic.  Otherwise 
     335        urst can be a tuple of length 4 or a baseWI. 
    33336         
    34         It is usually easier to use the \code{isomorphism_to} method of elliptic curves.  
    35         """ 
     337        Users will not usually need to use this class directly, but instead use 
     338        methods such as \code{isomorphism} of elliptic curves.    
     339 
     340        EXAMPLES: 
     341        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     342        sage: WeierstrassIsomorphism(EllipticCurve([0,1,2,3,4]),(-1,2,3,4)) 
     343        Generic morphism: 
     344        From: Abelian group of points on Elliptic Curve defined by y^2 + 2*y = x^3 + x^2 + 3*x + 4 over Rational Field 
     345        To:   Abelian group of points on Elliptic Curve defined by y^2 - 6*x*y - 10*y = x^3 - 2*x^2 - 11*x - 2 over Rational Field 
     346        Via:  (u,r,s,t) = (-1, 2, 3, 4) 
     347        sage: E=EllipticCurve([0,1,2,3,4]) 
     348        sage: F=EllipticCurve(E.cremona_label()) 
     349        sage: WeierstrassIsomorphism(E,None,F) 
     350        Generic morphism: 
     351        From: Abelian group of points on Elliptic Curve defined by y^2 + 2*y = x^3 + x^2 + 3*x + 4 over Rational Field 
     352        To:   Abelian group of points on Elliptic Curve defined by y^2  = x^3 + x^2 + 3*x + 5 over Rational Field 
     353        Via:  (u,r,s,t) = (1, 0, 0, -1) 
     354        sage: w=WeierstrassIsomorphism(None,(1,0,0,-1),F) 
     355        sage: w._domain_curve==E 
     356        True 
     357        """     
    36358        from ell_generic import is_EllipticCurve 
    37         if is_EllipticCurve(E): 
    38             K = F.base_ring() 
    39             a1, a2, a3, a4, a6 = E.ainvs() 
    40             b1, b2, b3, b4, b5 = F.ainvs() 
    41             try: 
    42                 D = (F.discriminant() / E.discriminant()) 
    43                 u = D.nth_root(12) 
    44             except ValueError: 
     359         
     360        if E!=None: 
     361            if not is_EllipticCurve(E): 
     362                raise ValueError, "First argument must be an elliptic curve or None" 
     363        if F!=None: 
     364            if not is_EllipticCurve(F): 
     365                raise ValueError, "Third argument must be an elliptic curve or None" 
     366        if urst!=None: 
     367            if len(urst)!=4: 
     368                raise ValueError, "Second argument must be [u,r,s,t] or None" 
     369        if len([par for par in [E,urst,F] if par!=None])<2: 
     370            raise ValueError, "At most 1 argument can be None" 
     371             
     372        if F==None:  # easy case 
     373            baseWI.__init__(self,*urst) 
     374            F=EllipticCurve(baseWI.__call__(self,E.a_invariants())) 
     375            Morphism.__init__(self, Hom(E(0).parent(), F(0).parent())) 
     376            self._domain_curve = E 
     377            self._codomain_curve = F 
     378            return 
     379 
     380        if E==None:  # easy case in reverse 
     381            baseWI.__init__(self,*urst) 
     382            inv_urst=baseWI.__invert__(self) 
     383            E=EllipticCurve(baseWI.__call__(inv_urst,F.a_invariants())) 
     384            Morphism.__init__(self, Hom(E(0).parent(), F(0).parent())) 
     385            self._domain_curve = E 
     386            self._codomain_curve = F 
     387            return 
     388 
     389        if urst==None: # try to construct the morphism 
     390            urst=isomorphisms(E,F,True) 
     391            if urst==None: 
    45392                raise ValueError, "Elliptic curves not isomorphic." 
    46             if u.parent() is not K or K.is_exact() and u**12 != D: 
    47                 raise ValueError, "Elliptic curves not isomorphic." 
    48             s = (a1*u - b1)/2 
    49             r = (a2*u**2 + a1*s*u - s**2 - b2)/3 
    50             t = (a3*u**3 - a1*r*u + 2*r*s - b3)/2 
    51             urst = u, r, s, t 
    52             if F != E.change_ring(K).change_weierstrass_model(urst): 
    53                 raise ValueError, "Elliptic curves not isomorphic." 
     393            baseWI.__init__(self, *urst) 
     394            Morphism.__init__(self, Hom(E(0).parent(), F(0).parent())) 
     395            self._domain_curve = E 
     396            self._codomain_curve = F 
     397            return 
     398                 
     399 
     400        # none of the parameters is None: 
     401        baseWI.__init__(self,*urst) 
     402        if F!=EllipticCurve(baseWI.__call__(self,E.a_invariants())): 
     403            raise ValueError, "second argument is not an isomorphism from first argument to third argument" 
    54404        else: 
    55             urst = F 
    56             F = E.change_weierstrass_model(urst) 
    57         Morphism.__init__(self, Hom(E(0).parent(), F(0).parent())) 
    58         self._urst = urst 
    59         self._domain_curve = E 
    60         self._codomain_curve = F 
     405            Morphism.__init__(self, Hom(E(0).parent(), F(0).parent())) 
     406            self._domain_curve = E 
     407            self._codomain_curve = F 
     408        return 
    61409     
    62     def _call_(self, P): 
     410    def __cmp__(self, other): 
     411        """ 
     412        EXAMPLE: 
     413        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     414        sage: E=EllipticCurve('389a1') 
     415        sage: F=E.change_weierstrass_model(1,2,3,4) 
     416        sage: w1=E.isomorphism_to(F) 
     417        sage: w1==w1 
     418        True 
     419        sage: w2 = F.automorphisms()[0] *w1 
     420        sage: w1==w2 
     421        False 
     422 
     423        sage: E=EllipticCurve(GF(7)(0)) 
     424        sage: F=E.change_weierstrass_model(2,3,4,5) 
     425        sage: a=E.isomorphisms(F) 
     426        sage: b=[w*a[0] for w in F.automorphisms()] 
     427        sage: b.sort() 
     428        sage: a==b 
     429        True 
     430        sage: c=[a[0]*w for w in E.automorphisms()] 
     431        sage: c.sort() 
     432        sage: a==c 
     433        True 
     434        """ 
     435        if not isinstance(other, WeierstrassIsomorphism): 
     436            return cmp(type(self), type(other)) 
     437        t = cmp(self._domain_curve, other._domain_curve) 
     438        if t: return t 
     439        t = cmp(self._codomain_curve, other._codomain_curve) 
     440        if t: return t 
     441        return baseWI.__cmp__(self,other) 
     442 
     443    def __call__(self, P): 
     444        r""" Allow a WeierstrassIsomorphism to be called as a function on 
     445        points of the domain curve to points of the codomain 
     446 
     447        EXAMPLES:  
     448        sage: from sage.schemes.elliptic_curves.weierstrass_morphism import * 
     449        sage: E=EllipticCurve('37a1') 
     450        sage: w=WeierstrassIsomorphism(E,(2,3,4,5)) 
     451        sage: P=E.gens()[0] 
     452        sage: w(P) 
     453        (-3/4 : 3/4 : 1) 
     454        sage: w(P).codomain()==E.change_weierstrass_model((2,3,4,5)) 
     455        True 
     456 
     457        In the last example, note that the curve a point lies on is 
     458        the "codomain" of the point since points are schemes... 
     459        """ 
    63460        if P[2] == 0: 
    64461            return self._codomain_curve(0) 
    65462        else: 
    66             x, y, z = P 
    67             u, r, s, t = self._urst 
    68             new_x = u**2*x+r 
    69             new_y = u**3*y + s*u**2*x + t 
    70             return self._codomain_curve((new_x, new_y, 1)) 
     463            return self._codomain_curve(baseWI.__call__(self,P.xy())) 
    71464 
    72465    def __invert__(self): 
    73         """ 
     466        """        
    74467        EXAMPLES:  
    75468            sage: E = EllipticCurve('5077') 
    76469            sage: F = E.change_weierstrass_model([2,3,4,5]); F 
    77             Elliptic Curve defined by y^2 - 8*x*y + 22*y = x^3 - 25*x^2 + 3*x + 588 over Rational Field 
     470            Elliptic Curve defined by y^2 + 4*x*y + 11/8*y = x^3 - 7/4*x^2 - 3/2*x - 9/32 over Rational Field 
    78471            sage: w = E.isomorphism_to(F) 
    79472            sage: P = E(-2,3,1) 
    80473            sage: w(P) 
    81             (-5 : -3 : 1) 
     474            (-5/4 : 9/4 : 1) 
    82475            sage: ~w 
    83476            Generic morphism: 
    84             From: Abelian group of points on Elliptic Curve defined by y^2 - 8*x*y + 22*y = x^3 - 25*x^2 + 3*x + 588 over Rational Field 
    85             To:   Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - 7*x + 6 over Rational Field 
     477                    From: Abelian group of points on Elliptic Curve defined by y^2 + 4*x*y + 11/8*y = x^3 - 7/4*x^2 - 3/2*x - 9/32 over Rational Field 
     478              To:   Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - 7*x + 6 over Rational Field 
     479              Via:  (u,r,s,t) = (1/2, -3/4, -2, 7/8) 
    86480            sage: Q = w(P); Q 
    87             (-5 : -3 : 1) 
     481            (-5/4 : 9/4 : 1) 
    88482            sage: (~w)(Q) 
    89483            (-2 : 3 : 1) 
    90484        """ 
    91         return WeierstrassIsomorphism(self._codomain_curve, self._domain_curve) 
    92                      
     485        winv=baseWI.__invert__(self).tuple() 
     486        return WeierstrassIsomorphism(self._codomain_curve, winv, self._domain_curve) 
     487         
     488    def __mul__(self,other): 
     489        """ 
     490        WeierstrassMorphisms can be composed using * if the codomain & 
     491        domain match: (w1*w2)()=w1(w2()) so require domain(w1)==codomain(w2) 
     492 
     493        EXAMPLES:  
     494            sage: E1 = EllipticCurve('5077') 
     495            sage: E2 = E1.change_weierstrass_model([2,3,4,5]) 
     496            sage: w1 = E1.isomorphism_to(E2) 
     497            sage: E3 = E2.change_weierstrass_model([6,7,8,9]) 
     498            sage: w2 = E2.isomorphism_to(E3) 
     499            sage: P = E1(-2,3,1) 
     500            sage: (w2*w1)(P)==w2(w1(P)) 
     501            True 
     502        """ 
     503        if self._domain_curve==other._codomain_curve: 
     504            w=baseWI.__mul__(self,other) 
     505            return WeierstrassIsomorphism(other._domain_curve, w.tuple(), self._codomain_curve) 
     506        else: 
     507            raise ValueError, "Domain of first argument must equal codomain of second" 
     508 
     509    def __repr__(self): 
     510        """ 
     511        Output of WeierstrassMorphisms outputs the underlying Morphism 
     512        with an extra line showing the (u,r,s,t) parameters 
     513        EXAMPLES:  
     514            sage: E1 = EllipticCurve('5077') 
     515            sage: E2 = E1.change_weierstrass_model([2,3,4,5]) 
     516            sage: E1.isomorphism_to(E2) 
     517            Generic morphism: 
     518            From: Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - 7*x + 6 over Rational Field 
     519            To:   Abelian group of points on Elliptic Curve defined by y^2 + 4*x*y + 11/8*y = x^3 - 7/4*x^2 - 3/2*x - 9/32 over Rational Field 
     520            Via:  (u,r,s,t) = (2, 3, 4, 5) 
     521        """ 
     522        return Morphism.__repr__(self)+"\n  Via:  (u,r,s,t) = "+baseWI.__repr__(self) 
     523 
     524 
Note: See TracChangeset for help on using the changeset viewer.