source: sage/rings/multi_polynomial_libsingular.pyx @ 4109:13b9e9dfbe0b

Revision 4109:13b9e9dfbe0b, 72.7 KB checked in by William Stein <wstein@…>, 6 years ago (diff)

Fix bug in new unpickle method for multivariate polynomial ring parents.

Line 
1"""
2Multivariate polynomials over QQ and GF(p) implemented using SINGULAR as backend.
3
4AUTHORS:
5    Martin Albrecht <malb@informatik.uni-bremen.de>
6
7
8
9TODO:
10   * implement every method from multi_polynomial_ring and multi_polynomial_element
11   * check SINGULAR code base for 'interesting' methods to add
12   * add required methods for F4, Buchberger etc.
13   * implement GF(p^n)
14   * implement block orderings
15   * implement Real, Complex
16   * test under CYGWIN (link error) 
17
18"""
19# We do this as we get a link error for init_csage(). However, on
20# obscure plattforms (Windows) we might need to link to csage anyway.
21
22cdef extern from "stdsage.h":
23    ctypedef void PyObject
24    object PY_NEW(object t)
25    int PY_TYPE_CHECK(object o, object t)
26    PyObject** FAST_SEQ_UNSAFE(object o)
27    void init_csage()
28
29    void  sage_free(void *p)
30    void* sage_realloc(void *p, size_t n)
31    void* sage_malloc(size_t)
32
33import os
34import sage.rings.memory
35
36from sage.libs.singular.singular import Conversion
37from sage.libs.singular.singular cimport Conversion
38
39cdef Conversion co
40co = Conversion()
41
42from sage.rings.multi_polynomial_ring import singular_name_mapping, TermOrder
43from sage.rings.multi_polynomial_ideal import MPolynomialIdeal
44from sage.rings.polydict import ETuple
45
46from sage.rings.rational_field import RationalField
47from sage.rings.finite_field import FiniteField_prime_modn
48
49from  sage.rings.rational cimport Rational
50
51from sage.interfaces.singular import singular as singular_default, is_SingularElement, SingularElement
52from sage.interfaces.macaulay2 import macaulay2 as macaulay2_default, is_Macaulay2Element
53from sage.structure.factorization import Factorization
54
55from complex_field import is_ComplexField
56from real_mpfr import is_RealField
57
58
59from sage.rings.integer_ring import IntegerRing
60from sage.structure.element cimport EuclideanDomainElement, \
61     RingElement, \
62     ModuleElement, \
63     Element, \
64     CommutativeRingElement
65
66from sage.rings.integer cimport Integer
67
68from sage.structure.parent cimport Parent
69from sage.structure.parent_base cimport ParentWithBase
70from sage.structure.parent_gens cimport ParentWithGens
71
72from sage.misc.sage_eval import sage_eval
73
74# shared library loading
75cdef extern from "dlfcn.h":
76    void *dlopen(char *, long)
77    char *dlerror()
78
79cdef extern from "string.h":
80    char *strdup(char *s)
81
82
83cdef init_singular():
84    """
85    This initializes the Singular library. Right now, this is a hack.
86
87    SINGULAR has a concept of compiled extension modules similar to
88    SAGE. For this, the compiled modules need to see the symbols from
89    the main programm. However, SINGULAR is a shared library in this
90    context these symbols are not known globally. The work around so
91    far is to load the library again and to specifiy RTLD_GLOBAL.
92    """
93
94    cdef void *handle
95
96
97    lib = os.environ['SAGE_LOCAL']+"/lib/libsingular.so"
98    if os.path.exists(lib):
99       handle = dlopen(lib, 256+1)
100    else:
101        lib = os.environ['SAGE_LOCAL']+"/lib/libsingular.dylib"
102        if os.path.exists(lib):
103            handle = dlopen(lib, 256+1)
104        else:
105            raise ImportError, "cannot load libSINGULAR library"
106
107    if handle == NULL:
108        print dlerror()
109
110    # Load Singular
111    siInit(lib)
112
113    # Steal Memory Manager back or weird things may happen
114    sage.rings.memory.pmem_malloc()
115
116# call it
117init_singular()
118
119cdef class MPolynomialRing_libsingular(MPolynomialRing_generic):
120    """
121    A multivariate polynomial ring over QQ or GF(p) implemented using SINGULAR.
122
123    """
124    def __init__(self, base_ring, n, names, order='degrevlex'):
125        """
126
127        Construct a multivariate polynomial ring subject to the following conditions:
128
129        INPUT:
130            base_ring -- base ring (must be either GF(p) (p prime) or QQ)
131            n -- number of variables (must be at least 1)
132            names -- names of ring variables, may be string of list/tuple
133            order -- term order (default: degrevlex)
134
135        EXAMPLES:
136            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
137            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
138            sage: P
139            Polynomial Ring in x, y, z over Rational Field
140
141            sage: f = 27/113 * x^2 + y*z + 1/2; f
142            27/113*x^2 + y*z + 1/2
143
144            sage: P.term_order()
145            Degree reverse lexicographic term order
146
147            sage: P = MPolynomialRing_libsingular(GF(127),3,names='abc', order='lex')
148            sage: P
149            Polynomial Ring in a, b, c over Finite Field of size 127
150
151            sage: a,b,c = P.gens()
152            sage: f = 57 * a^2*b + 43 * c + 1; f
153            57*a^2*b + 43*c + 1
154
155            sage: P.term_order()
156            Lexicographic term order
157       
158        """
159        cdef char **_names
160        cdef char *_name
161        cdef int i
162        cdef int characteristic
163
164        n = int(n)
165        if n<1:
166            raise ArithmeticError, "number of variables must be at least 1"
167
168        self.__ngens = n
169
170        MPolynomialRing_generic.__init__(self, base_ring, n, names, TermOrder(order))
171
172        self._has_singular = True
173
174        _names = <char**>sage_malloc(sizeof(char*)*len(self._names))
175
176        for i from 0 <= i < n:
177            _name = self._names[i]
178            _names[i] = strdup(_name)
179
180        if PY_TYPE_CHECK(base_ring, FiniteField_prime_modn):
181            characteristic = base_ring.characteristic()
182           
183        elif PY_TYPE_CHECK(base_ring, RationalField):
184            characteristic = 0
185
186        else:
187            raise NotImplementedError, "Only GF(p) and QQ are supported right now, sorry"
188
189        try:
190            order = singular_name_mapping[order]
191        except KeyError:
192            pass
193
194        self._ring = rDefault(characteristic, n, _names)
195        if(self._ring != currRing): rChangeCurrRing(self._ring)
196
197        rUnComplete(self._ring)
198
199        omFree(self._ring.wvhdl)
200        omFree(self._ring.order)
201        omFree(self._ring.block0)
202        omFree(self._ring.block1)
203
204        self._ring.wvhdl  = <int **>omAlloc0(3 * sizeof(int*))
205        self._ring.order  = <int *>omAlloc0(3* sizeof(int *))
206        self._ring.block0 = <int *>omAlloc0(3 * sizeof(int *))
207        self._ring.block1 = <int *>omAlloc0(3 * sizeof(int *))
208
209        if order == "dp":
210            self._ring.order[0] = ringorder_dp
211        elif order == "Dp":
212            self._ring.order[0] = ringorder_Dp
213        elif order == "lp":
214            self._ring.order[0] = ringorder_lp
215        elif order == "rp":
216            self._ring.order[0] = ringorder_rp
217        else:
218            self._ring.order[0] = ringorder_lp
219
220        self._ring.order[1] = ringorder_C
221
222        self._ring.block0[0] = 1
223        self._ring.block1[0] = n
224
225        rComplete(self._ring, 1)
226        self._ring.ShortOut = 0
227
228        self._zero = <MPolynomial_libsingular>new_MP(self,NULL)
229
230        sage_free(_names)
231
232    def __dealloc__(self):
233        """
234        """
235        rDelete(self._ring)
236
237    cdef _coerce_c_impl(self, element):
238        """
239
240        Coerces elements to self.
241
242        EXAMPLES:
243            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
244            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
245
246            We can coerce elements of self to self
247           
248            sage: P._coerce_(x*y + 1/2)
249            x*y + 1/2
250
251            We can coerce elements for a ring with the same algebraic properties
252
253            sage: R.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
254            sage: P == R
255            True
256
257            sage: P is R
258            False
259
260            sage: P._coerce_(x*y + 1)
261            x*y + 1
262
263            We can coerce base ring elements
264
265            sage: P._coerce_(3/2)
266            3/2
267
268            sage: P._coerce_(ZZ(1))
269            1
270
271            sage: P._coerce_(int(1))
272            1
273       
274        """
275        cdef poly *_p
276        cdef ring *_ring
277        cdef number *_n
278        cdef poly *rm
279        cdef int ok
280        cdef int i
281        cdef ideal *from_id, *to_id, *res_id
282           
283        _ring = self._ring
284
285        if(_ring != currRing): rChangeCurrRing(_ring)
286
287        if PY_TYPE_CHECK(element, MPolynomial_libsingular):
288            if element.parent() is <object>self:
289                return element
290            elif element.parent() == self:
291                # is this safe?
292                _p = p_Copy((<MPolynomial_libsingular>element)._poly, _ring)
293            else:
294                raise TypeError, "parents do not match"
295
296        elif PY_TYPE_CHECK(element, CommutativeRingElement):
297            # Accepting ZZ
298            if element.parent() is IntegerRing():
299                _p = p_ISet(int(element), _ring)
300
301            elif  <Parent>element.parent() is self._base:
302                # Accepting GF(p)
303                if PY_TYPE_CHECK(self._base, FiniteField_prime_modn):
304                    _p = p_ISet(int(element), _ring)
305
306                # Accepting QQ
307                elif PY_TYPE_CHECK(self._base, RationalField):
308                    _n = co.sa2si_QQ(element,_ring)
309                    _p = p_NSet(_n, _ring)
310                else:
311                    raise NotImplementedError
312            else:
313                raise TypeError, "base rings must be identical"
314
315        # Accepting int
316        elif PY_TYPE_CHECK(element, int):
317            _p = p_ISet(int(element), _ring)
318        else:
319            raise TypeError, "Cannot coerce element"
320
321        return new_MP(self,_p)
322
323    def __call__(self, element):
324        """
325        Construct a new element in self.
326
327        INPUT:
328            element -- several types are supported, see below
329       
330        EXAMPLE:
331            Call supports all conversions _coerce_ supports, plus:
332
333            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
334            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
335            sage: P('x+y + 1/4')
336            x + y + 1/4
337
338            sage: P._singular_()
339            //   characteristic : 0
340            //   number of vars : 3
341            //        block   1 : ordering dp
342            //                  : names    x y z
343            //        block   2 : ordering C
344
345            sage: P(singular('x + 3/4'))
346            x + 3/4
347        """
348        cdef poly *_m, *_p, *_tmp
349        cdef ring *_r
350
351        if PY_TYPE_CHECK(element, SingularElement):
352            element = str(element)
353       
354        if PY_TYPE_CHECK(element,str):
355            # let python do the the parsing
356            return sage_eval(element,self.gens_dict())
357
358               # this almost does what I want, besides variables with 0-9 in their names       
359##             _r = self._ring
360##             if(_r != currRing): rChangeCurrRing(_r)
361##             # improve this
362##             element = element.replace('^','').replace('**','').replace(' ','').strip()
363##             monomials = element.split('+')
364##             _p = p_ISet(0, _r)
365##             # wrong for e.g. x0,..,x7
366##             for m in monomials:
367##                 _m = p_ISet(1,_r)
368##                 for var in m.split('*'):
369##                     p_Read(var, _tmp, _r)
370##                     _m = p_Mult_q(_m,_tmp, _r)
371
372##                 _p = p_Add_q(_p,_m,_r)
373
374##             return new_MP(self,_p)           
375
376        return self._coerce_c_impl(element)
377
378    def _repr_(self):
379        """
380        EXAMPLE:
381            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
382            sage: P.<x,y> = MPolynomialRing_libsingular(QQ, 2)
383            sage: P
384            Polynomial Ring in x, y over Rational Field
385
386        """
387        varstr = ", ".join([ rRingVar(i,self._ring)  for i in range(self.__ngens) ])
388        return "Polynomial Ring in %s over %s"%(varstr,self._base)
389
390    def ngens(self):
391        """
392        Returns the number of variables in self.
393
394        EXAMPLES:
395            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
396            sage: P.<x,y> = MPolynomialRing_libsingular(QQ, 2)
397            sage: P.ngens()
398            2
399
400            sage: P = MPolynomialRing_libsingular(GF(127),1000,'x')
401            sage: P.ngens()
402            1000
403       
404        """
405        return int(self.__ngens)
406
407    def gens(self):
408        """
409        Return the tuple of variables in self.
410
411        EXAMPLES:
412            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
413            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
414            sage: P.gens()
415            (x, y, z)
416
417            sage: P = MPolynomialRing_libsingular(QQ,10,'x')
418            sage: P.gens()
419            (x0, x1, x2, x3, x4, x5, x6, x7, x8, x9)
420
421            sage: P.<SAGE,SINGULAR> = MPolynomialRing_libsingular(QQ,2) # weird names
422            sage: P.gens()
423            (SAGE, SINGULAR)
424
425        """
426        return tuple([self.gen(i) for i in range(self.__ngens)  ])
427
428    def gen(self, int n=0):
429        """
430        Returns the n-th generator of self.
431
432        EXAMPLES:
433            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
434            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
435            sage: P.gen(),P.gen(1)
436            (x, y)         
437
438            sage: P = MPolynomialRing_libsingular(GF(127),1000,'x')
439            sage: P.gen(500)
440            x500
441
442            sage: P.<SAGE,SINGULAR> = MPolynomialRing_libsingular(QQ,2) # weird names
443            sage: P.gen(1)
444            SINGULAR
445
446        """
447        cdef poly *_p
448        cdef ring *_ring
449
450        if n < 0 or n >= self.__ngens:
451            raise ValueError, "Generator not defined."
452
453        _ring = self._ring
454        _p = p_ISet(1,_ring)
455
456        # oddly enough, Singular starts counting a 1!!!
457        p_SetExp(_p, n+1, 1, _ring)
458        p_Setm(_p, _ring);
459
460        return new_MP(self,_p)
461
462    def ideal(self, gens, coerce=True):
463        """
464        Create an ideal in this polynomial ring.
465
466        INPUT:
467            gens -- generators of the ideal
468            coerce -- shall the generators be coerced first (default:True)
469
470        EXAMPLE:
471            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
472            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
473            sage: sage.rings.ideal.Katsura(P)
474            Ideal (x + 2*y + 2*z - 1, x^2 + 2*y^2 + 2*z^2 - x, 2*x*y + 2*y*z - y) of Polynomial Ring in x, y, z over Rational Field
475
476            sage: P.ideal([x + 2*y + 2*z-1, 2*x*y + 2*y*z-y, x^2 + 2*y^2 + 2*z^2-x])
477            Ideal (x + 2*y + 2*z - 1, 2*x*y + 2*y*z - y, x^2 + 2*y^2 + 2*z^2 - x) of Polynomial Ring in x, y, z over Rational Field
478
479        """
480        if is_SingularElement(gens):
481            gens = list(gens)
482            coerce = True
483        if is_Macaulay2Element(gens):
484            gens = list(gens)
485            coerce = True
486        elif not isinstance(gens, (list, tuple)):
487            gens = [gens]
488        if coerce:
489            gens = [self(x) for x in gens]  # this will even coerce from singular ideals correctly!
490        return MPolynomialIdeal(self, gens, coerce=False)
491       
492    def _macaulay2_(self, macaulay2=macaulay2_default):
493        """
494        Create a M2 representation of self if Macaulay2 is installed.
495
496        INPUT:
497            macaulay2 -- M2 interpreter (default: macaulay2_default)
498        """
499        try:
500            R = self.__macaulay2
501            if R is None or not (R.parent() is macaulay2):
502                raise ValueError
503            R._check_valid()
504            return R
505        except (AttributeError, ValueError):
506            if self.base_ring().is_prime_field():
507                if self.characteristic() == 0:
508                    base_str = "QQ"
509                else:
510                    base_str = "ZZ/" + str(self.characteristic())
511            elif is_IntegerRing(self.base_ring()):
512                base_str = "ZZ"
513            else:
514                raise TypeError, "no conversion of to a Macaulay2 ring defined"
515            self.__macaulay2 = macaulay2.ring(base_str, str(self.gens()), \
516                                              self.term_order().macaulay2_str())
517        return self.__macaulay2
518
519    def _singular_(self, singular=singular_default):
520        """
521        Create a SINGULAR (as in the CAS) representation of self. The
522        result is cached.
523
524        INPUT:
525            singular -- SINGULAR interpreter (default: singular_default)
526
527        EXAMPLES:
528            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
529            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
530            sage: P._singular_()
531            //   characteristic : 0
532            //   number of vars : 3
533            //        block   1 : ordering dp
534            //                  : names    x y z
535            //        block   2 : ordering C
536
537            sage: P._singular_() is P._singular_()
538            True
539
540            sage: P._singular_().name() == P._singular_().name()
541            True
542
543        TESTS:
544            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
545            sage: P.<x> = MPolynomialRing_libsingular(QQ,1)
546            sage: P._singular_()
547            //   characteristic : 0
548            //   number of vars : 1
549            //        block   1 : ordering lp
550            //                  : names    x
551            //        block   2 : ordering C
552       
553        """
554        try:
555            R = self.__singular
556            if R is None or not (R.parent() is singular):
557                raise ValueError
558            R._check_valid()
559            if self.base_ring().is_prime_field():
560                return R
561            if self.base_ring().is_finite():
562                R.set_ring() #sorry for that, but needed for minpoly
563                if  singular.eval('minpoly') != self.__minpoly:
564                    singular.eval("minpoly=%s"%(self.__minpoly))
565            return R
566        except (AttributeError, ValueError):
567            return self._singular_init_(singular)
568
569    def _singular_init_(self, singular=singular_default):
570        """
571        Create a SINGULAR (as in the CAS) representation of self. The
572        result is NOT cached.
573
574        INPUT:
575            singular -- SINGULAR interpreter (default: singular_default)
576
577        EXAMPLES:
578            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
579            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
580            sage: P._singular_init_()
581            //   characteristic : 0
582            //   number of vars : 3
583            //        block   1 : ordering dp
584            //                  : names    x y z
585            //        block   2 : ordering C
586            sage: P._singular_init_() is P._singular_init_()
587            False
588
589            sage: P._singular_init_().name() == P._singular_init_().name()
590            False
591
592        TESTS:
593            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
594            sage: P.<x> = MPolynomialRing_libsingular(QQ,1)
595            sage: P._singular_init_()
596            //   characteristic : 0
597            //   number of vars : 1
598            //        block   1 : ordering lp
599            //                  : names    x
600            //        block   2 : ordering C
601       
602        """
603        if self.ngens()==1:
604            _vars = str(self.gen())
605            if "*" in _vars: # 1.000...000*x
606                _vars = _vars.split("*")[1]
607            order = 'lp'
608        else:
609            _vars = str(self.gens())
610            order = self.term_order().singular_str()
611           
612        if is_RealField(self.base_ring()):
613            # singular converts to bits from base_10 in mpr_complex.cc by:
614            #  size_t bits = 1 + (size_t) ((float)digits * 3.5);
615            precision = self.base_ring().precision()
616            digits = sage.rings.arith.ceil((2*precision - 2)/7.0)
617            self.__singular = singular.ring("(real,%d,0)"%digits, _vars, order=order)
618
619        elif is_ComplexField(self.base_ring()):
620            # singular converts to bits from base_10 in mpr_complex.cc by:
621            #  size_t bits = 1 + (size_t) ((float)digits * 3.5);
622            precision = self.base_ring().precision()
623            digits = sage.rings.arith.ceil((2*precision - 2)/7.0)
624            self.__singular = singular.ring("(complex,%d,0,I)"%digits, _vars,  order=order)
625
626        elif self.base_ring().is_prime_field():
627            self.__singular = singular.ring(self.characteristic(), _vars, order=order)
628       
629        elif self.base_ring().is_finite(): #must be extension field
630            gen = str(self.base_ring().gen())
631            r = singular.ring( "(%s,%s)"%(self.characteristic(),gen), _vars, order=order)
632            self.__minpoly = "("+(str(self.base_ring().modulus()).replace("x",gen)).replace(" ","")+")"
633            singular.eval("minpoly=%s"%(self.__minpoly) )
634
635            self.__singular = r
636        else:   
637            raise TypeError, "no conversion to a Singular ring defined"
638
639        return self.__singular
640
641    def __hash__(self):
642        """
643        Return a hash for self, that is, a hash of the string representation of self
644
645        EXAMPLE:
646            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
647            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
648            sage: hash(P)
649            -6257278808099690586 # 64-bit
650            -1767675994 # 32-bit
651        """
652        return hash(self.__repr__())
653
654    def __richcmp__(left, right, int op):
655        return (<Parent>left)._richcmp(right, op)
656   
657    cdef int _cmp_c_impl(left, Parent right) except -2:
658        """
659        Multivariate polynomial rings are said to be equal if:
660         * their base rings match
661         * their generator names match
662         * their term orderings match
663
664        EXAMPLES:
665            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
666            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
667            sage: R.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
668            sage: P == R
669            True
670
671            sage: R.<x,y,z> = MPolynomialRing_libsingular(GF(127),3)
672            sage: P == R
673            False
674
675            sage: R.<x,y> = MPolynomialRing_libsingular(QQ,2)
676            sage: P == R
677            False
678
679            sage: R.<x,y,z> = MPolynomialRing_libsingular(QQ,3,order='revlex')
680            sage: P == R
681            False
682
683         
684        """
685        if PY_TYPE_CHECK(right, MPolynomialRing_libsingular):
686            return cmp( (left.base_ring(), map(str, left.gens()), left.term_order()),
687                        (right.base_ring(), map(str, right.gens()), right.term_order())
688                        )
689        else:
690            return cmp(type(left),type(right))
691
692    def __reduce__(self):
693        """
694        Serializes self.
695       
696        EXAMPLES:
697            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
698            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3, order='degrevlex')
699            sage: P == loads(dumps(P))
700            True
701
702            sage: P = MPolynomialRing_libsingular(GF(127),3,names='abc')
703            sage: P == loads(dumps(P))
704            True
705           
706        """
707        return sage.rings.multi_polynomial_libsingular.unpickle_MPolynomialRing_libsingular, ( self.base_ring(),
708                                                                                               map(str, self.gens()),
709                                                                                               self.term_order() )
710
711    ### The following methods are handy for implementing e.g. F4. They
712    ### do only superficial type/sanity checks and should be called
713    ### carefully.
714
715    def monomial_m_div_n(self, MPolynomial_libsingular f, MPolynomial_libsingular g, coeff=False):
716        """
717        Return f/g, where both f and g are treated as
718        monomials. Coefficients are ignored by default.
719
720        INPUT:
721            f -- monomial
722            g -- monomial
723            coeff -- divide coefficents as well (default: False)
724
725        EXAMPLE:
726            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
727            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
728            sage: P.monomial_m_div_n(3/2*x*y,x)
729            y
730
731            sage: P.monomial_m_div_n(3/2*x*y,x,coeff=True)
732            3/2*y
733
734        TESTS:
735            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
736            sage: R.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
737            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
738            sage: P.monomial_m_div_n(x*y,x)
739            y
740
741            sage: P.monomial_m_div_n(x*y,R.gen())
742            y
743
744            sage: P.monomial_m_div_n(P(0),P(1))
745            0
746
747            sage: P.monomial_m_div_n(P(1),P(0))
748            Traceback (most recent call last):
749            ...
750            ZeroDivisionError
751
752            sage: P.monomial_m_div_n(P(3/2),P(2/3), coeff=True)
753            9/4
754
755            sage: P.monomial_m_div_n(x,y) # Note the wrong result
756            x*y^1048575*z^1048575 # 64-bit
757            x*y^65535*z^65535 # 32-bit 
758
759            sage: P.monomial_m_div_n(x,P(1))
760            x
761
762        NOTE: Assumes that the head term of f is a multiple of the
763        head term of g and return the multiplicant m. If this rule is
764        violated, funny things may happen.
765
766       
767
768        """
769        cdef poly *res
770        cdef ring *r = self._ring
771       
772        if not <ParentWithBase>self is f._parent:
773            f = self._coerce_c(f)
774        if not <ParentWithBase>self is g._parent:
775            g = self._coerce_c(g)
776
777        if(r != currRing): rChangeCurrRing(r)
778
779        if not f._poly:
780            return self._zero
781        if not g._poly:
782            raise ZeroDivisionError
783       
784        res = pDivide(f._poly,g._poly)
785        if coeff:
786            p_SetCoeff(res, n_Div( p_GetCoeff(f._poly, r) , p_GetCoeff(g._poly, r), r), r)
787        else:
788            p_SetCoeff(res, n_Init(1, r), r)
789        return new_MP(self, res)
790   
791    def monomial_lcm(self, MPolynomial_libsingular f, MPolynomial_libsingular g):
792        """
793        LCM for monomials. Coefficients are ignored.
794       
795        INPUT:
796            f -- monomial
797            g -- monomial
798
799        EXAMPLE:
800            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
801            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
802            sage: P.monomial_lcm(3/2*x*y,x)
803            x*y
804
805        TESTS:
806            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
807            sage: R.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
808            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
809            sage: P.monomial_lcm(x*y,R.gen())
810            x*y
811
812            sage: P.monomial_lcm(P(3/2),P(2/3))
813            1
814
815            sage: P.monomial_lcm(x,P(1))
816            x
817
818        """
819        cdef poly *m = p_ISet(1,self._ring)
820       
821        if not <ParentWithBase>self is f._parent:
822            f = self._coerce_c(f)
823        if not <ParentWithBase>self is g._parent:
824            g = self._coerce_c(g)
825
826        if f._poly == NULL:
827            if g._poly == NULL:
828                return self._zero
829            else:
830                raise ArithmeticError, "cannot compute lcm of zero and nonzero element"
831        if g._poly == NULL:
832            raise ArithmeticError, "cannot compute lcm of zero and nonzero element"
833
834        if(self._ring != currRing): rChangeCurrRing(self._ring)
835       
836        pLcm(f._poly, g._poly, m)
837        return new_MP(self,m)
838       
839    def monomial_reduce(self, MPolynomial_libsingular f, G):
840        """
841        Try to find a g in G where g.lm() divides f. If found (g,flt)
842        is returned, (0,0) otherwise, where flt is f/g.lm().
843
844        It is assumed that G is iterable and contains ONLY elements in
845        self.
846       
847        INPUT:
848            f -- monomial
849            G -- list/set of mpolynomials
850           
851        EXAMPLES:
852            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
853            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
854            sage: f = x*y^2
855            sage: G = [ 3/2*x^3 + y^2 + 1/2, 1/4*x*y + 2/7, 1/2  ]
856            sage: P.monomial_reduce(f,G)
857            (1/4*x*y + 2/7, y)
858
859        TESTS:
860            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
861            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
862            sage: f = x*y^2
863            sage: G = [ 3/2*x^3 + y^2 + 1/2, 1/4*x*y + 2/7, 1/2  ]
864
865            sage: P.monomial_reduce(P(0),G)
866            (0, 0)
867
868            sage: P.monomial_reduce(f,[P(0)])
869            (0, 0)
870           
871        """
872        cdef poly *m = f._poly
873        cdef ring *r = self._ring
874        cdef poly *flt
875
876        if not m:
877            return f,f
878       
879        for g in G:
880            if PY_TYPE_CHECK(g, MPolynomial_libsingular) \
881                   and (<MPolynomial_libsingular>g) \
882                   and p_LmDivisibleBy((<MPolynomial_libsingular>g)._poly, m, r):
883                flt = pDivide(f._poly, (<MPolynomial_libsingular>g)._poly)
884                #p_SetCoeff(flt, n_Div( p_GetCoeff(f._poly, r) , p_GetCoeff((<MPolynomial_libsingular>g)._poly, r), r), r)
885                p_SetCoeff(flt, n_Init(1, r), r)
886                return g, new_MP(self,flt)
887        return self._zero,self._zero
888
889    def monomial_pairwise_prime(self, MPolynomial_libsingular g, MPolynomial_libsingular h):
890        """
891        Return True if h and g are pairwise prime. Both are treated as monomials.
892
893        INPUT:
894            h -- monomial
895            g -- monomial
896
897        EXAMPLES:
898            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
899            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
900            sage: P.monomial_pairwise_prime(x^2*z^3, y^4)
901            True
902
903            sage: P.monomial_pairwise_prime(1/2*x^3*y^2, 3/4*y^3)
904            False
905
906        TESTS:
907            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
908            sage: Q.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
909            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
910            sage: P.monomial_pairwise_prime(x^2*z^3, Q('y^4'))
911            True
912
913            sage: P.monomial_pairwise_prime(1/2*x^3*y^2, Q(0))
914            True
915
916            sage: P.monomial_pairwise_prime(P(1/2),x)
917            False
918
919       
920        """
921        cdef int i
922        cdef ring *r
923        cdef poly *p, *q
924
925        if h._parent is not g._parent:
926            g = (<MPolynomialRing_libsingular>h._parent)._coerce_c(g)
927
928        r = (<MPolynomialRing_libsingular>h._parent)._ring
929        p = g._poly
930        q = h._poly
931
932        if p == NULL:
933            if q == NULL:
934                return False #GCD(0,0) = 0
935            else:
936                return True #GCD(x,0) = 1
937
938        elif q == NULL:
939            return True # GCD(0,x) = 1
940
941        elif p_IsConstant(p,r) or p_IsConstant(q,r): # assuming a base field
942            return False
943
944        for i from 1 <= i <= r.N:
945            if p_GetExp(p,i,r) and p_GetExp(q,i,r):
946                return False
947        return True
948
949    def monomial_is_divisible_by(self, MPolynomial_libsingular a, MPolynomial_libsingular b):
950        """
951        Return 0 if b does not divide a and the factor
952        otherwise. Coefficients are ignored.
953       
954        INPUT:
955            a -- monomial
956            b -- monomial
957
958        EXAMPLES:
959            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
960            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
961            sage: P.monomial_is_divisible_by(x^3*y^2*z^4, x*y*z)
962            x^2*y*z^3
963            sage: P.monomial_is_divisible_by(x*y*z, x^3*y^2*z^4)
964            0
965
966        TESTS:
967            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
968            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
969            sage: P.monomial_is_divisible_by(P(0),P(1))
970            0
971            sage: P.monomial_is_divisible_by(x,P(1))
972            x
973
974        """
975        cdef poly *_a
976        cdef poly *_b
977        cdef ring *_r
978        cdef poly *ret
979        if a._parent is not b._parent:
980            b = (<MPolynomialRing_libsingular>a._parent)._coerce_c(b)
981
982        _a = a._poly
983        _b = b._poly
984        _r = (<MPolynomialRing_libsingular>a._parent)._ring
985        if(_r != currRing): rChangeCurrRing(_r)
986
987        if _b == NULL:
988            raise ZeroDivisionError
989        if _a == NULL:
990            return self._zero
991       
992        if not p_DivisibleBy(_b, _a, _r):
993           return self._zero
994        else:
995            ret = pDivide(_a,_b)
996            p_SetCoeff(ret, n_Init(1, _r), _r)
997            return new_MP(self,ret)
998   
999def unpickle_MPolynomialRing_libsingular(base_ring, names, term_order):
1000    """
1001    inverse function for MPolynomialRing_libsingular.__reduce__
1002
1003    """
1004    return MPolynomialRing_libsingular(base_ring, len(names), names, term_order)
1005
1006cdef MPolynomial_libsingular new_MP(MPolynomialRing_libsingular parent, poly *juice):
1007    """
1008    Construct a new MPolynomial_libsingular element
1009    """
1010    cdef MPolynomial_libsingular p
1011    p = PY_NEW(MPolynomial_libsingular)
1012    p._parent = <ParentWithBase>parent
1013    p._poly = juice
1014    return p
1015
1016cdef class MPolynomial_libsingular(sage.rings.multi_polynomial.MPolynomial):
1017    """
1018    A multivariate polynomial implemented using libSINGULAR.
1019    """
1020    def __init__(self, MPolynomialRing_libsingular parent):
1021        """
1022        Construct a zero element in parent.
1023        """
1024        self._parent = <ParentWithBase>parent
1025
1026    def __dealloc__(self):
1027        if self._poly:
1028            p_Delete(&self._poly, (<MPolynomialRing_libsingular>self._parent)._ring)
1029
1030    def __call__(self, *x):
1031        """
1032        Evaluate a polynomial at the given point x
1033
1034        INPUT:
1035            x -- a list of elements in self.parent()
1036
1037        EXAMPLES:
1038            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1039            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
1040            sage: f = 3/2*x^2*y + 1/7 * y^2 + 13/27
1041            sage: f(0,0,0)
1042            13/27
1043
1044            sage: f(1,1,1)
1045            803/378
1046            sage: 3/2 + 1/7 + 13/27
1047            803/378
1048
1049            sage: f(45/2,19/3,1)
1050            7281167/1512
1051
1052        TESTS:
1053            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1054            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
1055            sage: P(0)(1,2,3)
1056            0
1057            sage: P(3/2)(1,2,3)
1058            3/2
1059        """
1060        cdef int l = len(x)
1061        cdef MPolynomialRing_libsingular parent = (<MPolynomialRing_libsingular>self._parent)
1062        cdef ring *_ring = parent._ring
1063
1064        cdef poly *_p
1065
1066        if l != parent._ring.N:
1067            raise TypeError, "number of arguments does not match number of variables in parent"
1068       
1069        cdef ideal *to_id = idInit(l,1)
1070
1071        try:
1072            for i from 0 <= i < l:
1073                e = x[i] # TODO: optimize this line
1074                to_id.m[i]= p_Copy( (<MPolynomial_libsingular>(<MPolynomialRing_libsingular>parent._coerce_c(x[i])))._poly, _ring)
1075
1076        except TypeError:
1077            id_Delete(&to_id, _ring)
1078            raise TypeError, "cannot coerce in arguments"
1079
1080        cdef ideal *from_id=idInit(1,1)
1081        from_id.m[0] = p_Copy(self._poly, _ring)
1082
1083        cdef ideal *res_id = fast_map(from_id, _ring, to_id, _ring)
1084        cdef poly *res = p_Copy(res_id.m[0], _ring)
1085
1086        id_Delete(&to_id, _ring)
1087        id_Delete(&from_id, _ring)
1088        id_Delete(&res_id, _ring)
1089        return new_MP(parent, res)
1090           
1091    def __richcmp__(left, right, int op):
1092        return (<Element>left)._richcmp(right, op)
1093
1094    cdef int _cmp_c_impl(left, Element right) except -2:
1095        """
1096        Compare left and right and return -1, 0, and 1 for <,==, and > respectively.
1097
1098        EXAMPLES:
1099            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1100            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3, order='degrevlex')
1101            sage: x == x
1102            True
1103
1104            sage: x > y
1105            True
1106            sage: y^2 > x
1107            True
1108
1109            sage: (2/3*x^2 + 1/2*y + 3) > (2/3*x^2 + 1/4*y + 10)
1110            True
1111
1112        TESTS:
1113            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1114            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3, order='degrevlex')
1115            sage: x > P(0)
1116            True
1117
1118            sage: P(0) == P(0)
1119            True
1120
1121            sage: P(0) < P(1)
1122            True
1123
1124            sage: x > P(1)
1125            True
1126           
1127            sage: 1/2*x < 3/4*x
1128            True
1129
1130            sage: (x+1) > x
1131            True
1132
1133            sage: f = 3/4*x^2*y + 1/2*x + 2/7
1134            sage: f > f
1135            False
1136            sage: f < f
1137            False
1138            sage: f == f
1139            True
1140
1141            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1142            sage: P.<x,y,z> = MPolynomialRing_libsingular(GF(127),3, order='degrevlex')
1143            sage: (66*x^2 + 23) > (66*x^2 + 2)
1144            True
1145
1146       
1147        """
1148        cdef ring *r
1149        cdef poly *p, *q
1150        cdef number *h
1151        cdef int ret = 0
1152       
1153        r = (<MPolynomialRing_libsingular>left._parent)._ring
1154        if(r != currRing): rChangeCurrRing(r)
1155        p = (<MPolynomial_libsingular>left)._poly
1156        q = (<MPolynomial_libsingular>right)._poly
1157
1158        # handle special cases first (slight slowdown, as special
1159        # cases are - well - special
1160        if p==NULL:
1161            if q==NULL:
1162                # compare 0, 0
1163                return 0
1164            elif p_IsConstant(q,r):
1165                # compare 0, const
1166                return 1-2*n_GreaterZero(p_GetCoeff(q,r), r) # -1: <, 1: > #
1167        elif q==NULL:
1168            if p_IsConstant(p,r):
1169                # compare const, 0
1170                return -1+2*n_GreaterZero(p_GetCoeff(p,r), r) # -1: <, 1: >
1171        #else
1172       
1173        while ret==0 and p!=NULL and q!=NULL:
1174            ret = p_Cmp( p, q, r)
1175       
1176            if ret==0:
1177                h = n_Sub(p_GetCoeff(p, r),p_GetCoeff(q, r), r)
1178                # compare coeffs
1179                ret = -1+n_IsZero(h, r)+2*n_GreaterZero(h, r) # -1: <, 0:==, 1: >
1180                n_Delete(&h, r)
1181            p = pNext(p)
1182            q = pNext(q)
1183
1184        if ret==0:
1185            if p==NULL and q != NULL:
1186                ret = -1
1187            elif p!=NULL and q==NULL:
1188                ret = 1
1189
1190        return ret
1191
1192    cdef ModuleElement _add_c_impl( left, ModuleElement right):
1193        """
1194        Add left and right.
1195
1196        EXAMPLE:
1197            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1198            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
1199            sage: 3/2*x + 1/2*y + 1
1200            3/2*x + 1/2*y + 1
1201
1202        """
1203        cdef MPolynomial_libsingular res
1204       
1205        cdef poly *_l, *_r, *_p
1206        cdef ring *_ring
1207
1208        _ring = (<MPolynomialRing_libsingular>left._parent)._ring
1209
1210        _l = p_Copy(left._poly, _ring)
1211        _r = p_Copy((<MPolynomial_libsingular>right)._poly, _ring)
1212
1213        if(_ring != currRing): rChangeCurrRing(_ring)
1214        _p= p_Add_q(_l, _r, _ring)
1215
1216        return new_MP((<MPolynomialRing_libsingular>left._parent),_p)
1217
1218    cdef ModuleElement _sub_c_impl( left, ModuleElement right):
1219        """
1220        Subtract left and right.
1221
1222        EXAMPLE:
1223            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1224            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
1225            sage: 3/2*x - 1/2*y - 1
1226            3/2*x - 1/2*y - 1
1227
1228        """
1229        cdef MPolynomial_libsingular res
1230       
1231        cdef poly *_l, *_r, *_p
1232        cdef ring *_ring
1233
1234        _ring = (<MPolynomialRing_libsingular>left._parent)._ring
1235
1236        _l = p_Copy(left._poly, _ring)
1237        _r = p_Copy((<MPolynomial_libsingular>right)._poly, _ring)
1238
1239        if(_ring != currRing): rChangeCurrRing(_ring)
1240        _p= p_Add_q(_l, p_Neg(_r, _ring), _ring)
1241
1242        return new_MP((<MPolynomialRing_libsingular>left._parent),_p)
1243
1244
1245    cdef ModuleElement _rmul_c_impl(self, RingElement left):
1246        """
1247        Multiply self with a base ring element.
1248
1249        EXAMPLE:
1250            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1251            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
1252            sage: 3/2*x
1253            3/2*x
1254        """
1255
1256        cdef number *_n
1257        cdef ring *_ring
1258        cdef poly *_p
1259
1260        _ring = (<MPolynomialRing_libsingular>self._parent)._ring
1261       
1262        if(_ring != currRing): rChangeCurrRing(_ring)
1263       
1264        if PY_TYPE_CHECK((<MPolynomialRing_libsingular>self._parent)._base, FiniteField_prime_modn):
1265            _n = n_Init(int(left),_ring)
1266           
1267        elif PY_TYPE_CHECK((<MPolynomialRing_libsingular>self._parent)._base, RationalField):
1268            _n = co.sa2si_QQ(left,_ring)
1269
1270        _p = pp_Mult_nn(self._poly,_n,_ring)
1271        n_Delete(&_n, _ring)
1272        return new_MP((<MPolynomialRing_libsingular>self._parent),_p)
1273   
1274    cdef RingElement  _mul_c_impl(left, RingElement right):
1275        """
1276        Multiply left and right.
1277
1278        EXAMPLE:
1279            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1280            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
1281            sage: (3/2*x - 1/2*y - 1) * (3/2*x + 1/2*y + 1)
1282            9/4*x^2 - 1/4*y^2 - y - 1
1283        """
1284        cdef poly *_l, *_r, *_p
1285        cdef ring *_ring
1286
1287        _ring = (<MPolynomialRing_libsingular>left._parent)._ring
1288
1289        if(_ring != currRing): rChangeCurrRing(_ring)
1290        _p = pp_Mult_qq(left._poly, (<MPolynomial_libsingular>right)._poly, _ring)
1291
1292        return new_MP(left._parent,_p)
1293
1294    cdef RingElement  _div_c_impl(left, RingElement right):
1295        """
1296        Divide left by right
1297
1298        EXAMPLES:
1299            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1300            sage: R.<x,y>=MPolynomialRing_libsingular(QQ,2)
1301            sage: f = (x + y)/3
1302            sage: f.parent()
1303            Polynomial Ring in x, y over Rational Field
1304
1305        Note that / is still a constructor for elements of the
1306        fraction field in all cases as long as both arguments have the
1307        same parent.
1308
1309            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1310            sage: R.<x,y>=MPolynomialRing_libsingular(QQ,2)
1311            sage: f = x^3 + y
1312            sage: g = x
1313            sage: h = f/g; h
1314            (x^3 + y)/x
1315            sage: h.parent()
1316            Fraction Field of Polynomial Ring in x, y over Rational Field
1317
1318        TESTS:
1319            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1320            sage: R.<x,y>=MPolynomialRing_libsingular(QQ,2)
1321            sage: x/0
1322            Traceback (most recent call last):
1323            ...
1324            ZeroDivisionError
1325
1326        """
1327        cdef poly *p
1328        cdef ring *r
1329        cdef number *n
1330        if (<MPolynomial_libsingular>right).is_constant_c():
1331
1332            p = (<MPolynomial_libsingular>right)._poly
1333            if p == NULL:
1334                raise ZeroDivisionError
1335            r = (<MPolynomialRing_libsingular>(<MPolynomial_libsingular>left)._parent)._ring
1336            n = p_GetCoeff(p, r)
1337            n = nInvers(n)
1338            p = pp_Mult_nn(left._poly, n,  r)
1339            n_Delete(&n,r)
1340            return new_MP(left._parent, p)
1341        else:
1342            return (<MPolynomialRing_libsingular>left._parent).fraction_field()(left,right)
1343
1344    def __pow__(MPolynomial_libsingular self,int exp,ignored):
1345        """
1346        """
1347        cdef ring *_ring
1348        _ring = (<MPolynomialRing_libsingular>self._parent)._ring
1349
1350        cdef poly *_p
1351
1352
1353        if exp < 0:
1354            raise ArithmeticError, "Cannot comput negative powers of polynomials"
1355
1356        if(_ring != currRing): rChangeCurrRing(_ring)
1357        _p = pPower( p_Copy(self._poly,(<MPolynomialRing_libsingular>self._parent)._ring),exp)
1358        return new_MP((<MPolynomialRing_libsingular>self._parent),_p)
1359
1360
1361    def __neg__(self):
1362        cdef ring *_ring
1363        _ring = (<MPolynomialRing_libsingular>self._parent)._ring
1364        if(_ring != currRing): rChangeCurrRing(_ring)
1365
1366        return new_MP((<MPolynomialRing_libsingular>self._parent),\
1367                                           p_Neg(p_Copy(self._poly,_ring),_ring))
1368
1369    def _repr_(self):
1370        s =  self._repr_short_c()
1371        s = s.replace("+"," + ").replace("-"," - ")
1372        if s.startswith(" - "):
1373            return "-" + s[3:]
1374        else:
1375            return s
1376
1377    def _repr_short(self):
1378        """
1379        This is a faster but less pretty way to print polynomials. If available
1380        it uses the short SINGULAR notation.
1381        """
1382        cdef ring *_ring = (<MPolynomialRing_libsingular>self._parent)._ring
1383        if _ring.CanShortOut:
1384            _ring.ShortOut = 1
1385            s = self._repr_short_c()
1386            _ring.ShortOut = 0
1387        else:
1388            s = self._repr_short_c()
1389        return s
1390
1391    cdef _repr_short_c(self):
1392        """
1393        Raw SINGULAR printing.
1394        """
1395        s = p_String(self._poly, (<MPolynomialRing_libsingular>self._parent)._ring, (<MPolynomialRing_libsingular>self._parent)._ring)
1396        return s
1397
1398    def _latex(self):
1399        #do it for real, may be slow
1400        raise NotImplementedError
1401
1402    def _repr_with_changed_varnames(self, varnames):
1403        #plan: replace self._ring.names for the moment with varnames
1404        raise NotImplementedError
1405   
1406    def degree(self, MPolynomial_libsingular x=None):
1407        """
1408        Return the maximal degree of self in x, where x must be one of the
1409        generators for the parent of self.
1410
1411        INPUT:
1412            x -- multivariate polynmial (a generator of the parent of self)
1413                 If x is not specified (or is None), return the total degree,
1414                 which is the maximum degree of any monomial.
1415
1416        OUTPUT:
1417            integer
1418       
1419        EXAMPLE:
1420            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1421            sage: R.<x, y> = MPolynomialRing_libsingular(QQ, 2)
1422            sage: f = y^2 - x^9 - x
1423            sage: f.degree(x)
1424            9
1425            sage: f.degree(y)
1426            2
1427            sage: (y^10*x - 7*x^2*y^5 + 5*x^3).degree(x)
1428            3
1429            sage: (y^10*x - 7*x^2*y^5 + 5*x^3).degree(y)
1430            10
1431
1432        TESTS:
1433            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1434            sage: P.<x, y> = MPolynomialRing_libsingular(QQ, 2)
1435            sage: P(0).degree(x)
1436            0
1437            sage: P(1).degree(x)
1438            0
1439
1440        """
1441        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
1442        cdef poly *p = self._poly
1443        cdef int deg, _deg
1444
1445        deg = 0
1446       
1447        if not x:
1448            return self.total_degree()
1449
1450        # TODO: we can do this faster
1451        if not x in self._parent.gens():
1452            raise TypeError, "x must be one of the generators of the parent."
1453        for i from 1 <= i <= r.N:
1454            if p_GetExp(x._poly, i, r):
1455                break
1456        while p:
1457            _deg =  p_GetExp(p,i,r)
1458            if _deg > deg:
1459                deg = _deg
1460            p = pNext(p)
1461
1462        return deg
1463
1464    def newton_polytope(self):
1465        """
1466        Return the Newton polytope of this polynomial.
1467
1468        You should have the optional polymake package installed.
1469
1470        EXAMPLES:
1471            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1472            sage: R.<x,y> = MPolynomialRing_libsingular(QQ,2)
1473            sage: f = 1 + x*y + x^3 + y^3
1474            sage: P = f.newton_polytope()
1475            sage: P
1476            Convex hull of points [[1, 0, 0], [1, 0, 3], [1, 1, 1], [1, 3, 0]]
1477            sage: P.facets()
1478            [(0, 1, 0), (3, -1, -1), (0, 0, 1)]
1479            sage: P.is_simple()
1480            True
1481
1482        TESTS:
1483            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1484            sage: R.<x,y> = MPolynomialRing_libsingular(QQ,2)
1485            sage: R(0).newton_polytope()
1486            Convex hull of points []
1487            sage: R(1).newton_polytope()
1488            Convex hull of points [[1, 0, 0]]
1489           
1490        """
1491        from sage.geometry.all import polymake
1492        e = self.exponents()
1493        a = [[1] + list(v) for v in e]
1494        P = polymake.convex_hull(a)
1495        return P
1496
1497    def total_degree(self):
1498        """
1499        Return the total degree of self, which is the maximum degree
1500        of all monomials in self.
1501
1502        EXAMPLES:
1503            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1504            sage: R.<x,y,z> = MPolynomialRing_libsingular(QQ, 3)
1505            sage: f=2*x*y^3*z^2
1506            sage: f.total_degree()
1507            6
1508            sage: f=4*x^2*y^2*z^3
1509            sage: f.total_degree()
1510            7
1511            sage: f=99*x^6*y^3*z^9
1512            sage: f.total_degree()
1513            18
1514            sage: f=x*y^3*z^6+3*x^2
1515            sage: f.total_degree()
1516            10
1517            sage: f=z^3+8*x^4*y^5*z
1518            sage: f.total_degree()
1519            10
1520            sage: f=z^9+10*x^4+y^8*x^2
1521            sage: f.total_degree()
1522            10
1523
1524        TESTS:
1525            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1526            sage: R.<x,y,z> = MPolynomialRing_libsingular(QQ, 3)
1527            sage: R(0).total_degree()
1528            0
1529            sage: R(1).total_degree()
1530            0
1531        """
1532        cdef poly *p = self._poly
1533        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
1534        cdef int l
1535        if self._poly == NULL:
1536            return 0
1537        if(r != currRing): rChangeCurrRing(r)
1538        return pLDeg(p,&l,r)
1539
1540    def monomial_coefficient(self, MPolynomial_libsingular mon):
1541        """
1542        Return the coefficient of the monomial mon in self, where mon
1543        must have the same parent as self.
1544
1545        INPUT:
1546            mon -- a monomial
1547
1548        OUTPUT:
1549            ring element
1550       
1551        EXAMPLE:
1552            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
1553            sage: P.<x,y> = MPolynomialRing_libsingular(QQ, 2)
1554
1555        The coefficient returned is an element of the base ring of self; in
1556        this case, QQ.
1557       
1558            sage: f = 2 * x * y
1559            sage: c = f.monomial_coefficient(x*y); c
1560            2
1561            sage: c in QQ
1562            True
1563
1564            sage: f = y^2 - x^9 - 7*x + 5*x*y
1565            sage: f.monomial_coefficient(y^2)
1566            1
1567            sage: f.monomial_coefficient(x*y)
1568            5
1569            sage: f.monomial_coefficient(x^9)
1570            -1
1571            sage: f.monomial_coefficient(x^10)
1572            0
1573        """
1574        cdef poly *p = self._poly
1575        cdef poly *m = mon._poly
1576        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
1577
1578        if not mon._parent is self._parent:
1579            raise TypeError, "mon must have same parent as self"
1580       
1581        while(p):
1582            if p_ExpVectorEqual(p, m, r) == 1:
1583                return co.si2sa(p_GetCoeff(p, r), r, (<MPolynomialRing_libsingular>self._parent)._base)
1584            p = pNext(p)
1585
1586        return (<MPolynomialRing_libsingular>self._parent)._base(0)
1587
1588    def dict(self):
1589        cdef poly *p
1590        cdef ring *r
1591        cdef int n
1592        cdef int v
1593        r = (<MPolynomialRing_libsingular>self._parent)._ring
1594        base = (<MPolynomialRing_libsingular>self._parent)._base
1595        p = self._poly
1596        pd = dict()
1597        while p:
1598            d = dict()
1599            # assuming that SINGULAR stores monomial dense
1600            for v from 1 <= v <= r.N:
1601                n = p_GetExp(p,v,r)
1602                if n!=0:
1603                    d[v-1] = n
1604               
1605            pd[ETuple(d,r.N)] = co.si2sa(p_GetCoeff(p, r), r, base)
1606
1607            p = pNext(p)
1608        return pd
1609
1610    def __getitem__(self,x):
1611        """
1612        same as self.monomial_coefficent but for exponent vectors.
1613       
1614        INPUT:
1615            x -- a tuple or, in case of a single-variable MPolynomial
1616                 ring x can also be an integer.
1617       
1618        EXAMPLES:
1619            sage: R.<x, y> = PolynomialRing(QQ, 2)
1620            sage: f = -10*x^3*y + 17*x*y
1621            sage: f[3,1]
1622            -10
1623            sage: f[1,1]
1624            17
1625            sage: f[0,1]
1626            0
1627
1628            sage: R.<x> = PolynomialRing(GF(7),1); R
1629            Polynomial Ring in x over Finite Field of size 7
1630            sage: f = 5*x^2 + 3; f
1631            3 + 5*x^2
1632            sage: f[2]
1633            5
1634        """
1635
1636        cdef poly *m
1637        cdef poly *p = self._poly
1638        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
1639        cdef int i
1640
1641        if PY_TYPE_CHECK(x, MPolynomial_libsingular):
1642            return self.monomial_coefficient(x)
1643        if not PY_TYPE_CHECK(x, tuple):
1644            try:
1645                x = tuple(x)
1646            except TypeError:
1647                x = (x,)
1648
1649        if len(x) != (<MPolynomialRing_libsingular>self._parent).__ngens:
1650            raise TypeError, "x must have length self.ngens()"
1651
1652        m = p_ISet(1,r)
1653        i = 1
1654        for e in x:
1655            p_SetExp(m, i, int(e), r)
1656            i += 1
1657        p_Setm(m, r)
1658
1659        while(p):
1660            if p_ExpVectorEqual(p, m, r) == 1:
1661                p_Delete(&m,r)
1662                return co.si2sa(p_GetCoeff(p, r), r, (<MPolynomialRing_libsingular>self._parent)._base)
1663            p = pNext(p)
1664
1665        p_Delete(&m,r)
1666        return (<MPolynomialRing_libsingular>self._parent)._base(0)
1667
1668    def coefficient(self, mon):
1669        raise NotImplementedError
1670
1671    def exponents(self):
1672        """
1673        Return the exponents of the monomials appearing in self.
1674       
1675        EXAMPLES:
1676            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
1677            sage: R.<a,b,c> = MPolynomialRing_libsingular(QQ, 3)
1678            sage: f = a^3 + b + 2*b^2
1679            sage: f.exponents()
1680            [(3, 0, 0), (0, 2, 0), (0, 1, 0)]
1681        """
1682        cdef poly *p
1683        cdef ring *r
1684        cdef int v
1685        r = (<MPolynomialRing_libsingular>self._parent)._ring
1686
1687        p = self._poly
1688
1689        pl = list()
1690        while p:
1691            ml = list()
1692            for v from 1 <= v <= r.N:
1693                ml.append(p_GetExp(p,v,r))
1694            pl.append(tuple(ml))
1695
1696            p = pNext(p)
1697        return pl
1698
1699    def is_unit(self):
1700        """
1701        Return True if self is a unit.
1702
1703        EXAMPLES:
1704            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
1705            sage: R.<x,y> = MPolynomialRing_libsingular(QQ, 2)
1706            sage: (x+y).is_unit()
1707            False
1708            sage: R(0).is_unit()
1709            False
1710            sage: R(-1).is_unit()
1711            True
1712            sage: R(-1 + x).is_unit()
1713            False
1714            sage: R(2).is_unit()
1715            True
1716        """
1717        return bool(p_IsUnit(self._poly, (<MPolynomialRing_libsingular>self._parent)._ring))
1718
1719    def inverse_of_unit(self):
1720        """
1721        Return the inverse of self if self is a unit.
1722
1723        EXAMPLES:
1724            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
1725            sage: R.<x,y> = MPolynomialRing_libsingular(QQ, 2)
1726        """
1727        cdef ring *_ring = (<MPolynomialRing_libsingular>self._parent)._ring
1728        if(_ring != currRing): rChangeCurrRing(_ring)
1729
1730        if not p_IsUnit(self._poly, _ring):
1731            raise ArithmeticError, "is not a unit"
1732        else:
1733            return new_MP(self._parent,pInvers(0,self._poly,NULL))
1734
1735    def is_homogeneous(self):
1736        """
1737        Return True if self is a homogeneous polynomial.
1738
1739        EXAMPLES:
1740            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
1741            sage: P.<x,y> = MPolynomialRing_libsingular(RationalField(), 2)
1742            sage: (x+y).is_homogeneous()
1743            True
1744            sage: (x.parent()(0)).is_homogeneous()
1745            True
1746            sage: (x+y^2).is_homogeneous()
1747            False
1748            sage: (x^2 + y^2).is_homogeneous()
1749            True
1750            sage: (x^2 + y^2*x).is_homogeneous()
1751            False
1752            sage: (x^2*y + y^2*x).is_homogeneous()
1753            True
1754       
1755        """
1756        cdef ring *_ring = (<MPolynomialRing_libsingular>self._parent)._ring
1757        if(_ring != currRing): rChangeCurrRing(_ring)
1758        return bool(pIsHomogeneous(self._poly))
1759
1760    def homogenize(self, var='h'):
1761        raise NotImplementedError
1762
1763    def is_monomial(self):
1764        return not self._poly.next
1765
1766    def fix(self, fixed):
1767        """
1768        Fixes some given variables in a given multivariate polynomial and
1769        returns the changed multivariate polynomials. The polynomial
1770        itself is not affected.  The variable,value pairs for fixing are
1771        to be provided as dictionary of the form {variable:value}.
1772
1773        This is a special case of evaluating the polynomial with some of
1774        the variables constants and the others the original variables, but
1775        should be much faster if only few variables are to be fixed.
1776
1777        INPUT:
1778            fixed -- dict with variable:value pairs
1779
1780        OUTPUT:
1781            new MPolynomial
1782
1783        EXAMPLES:
1784            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
1785            sage: x, y = MPolynomialRing_libsingular(QQ,2,'xy').gens() 
1786            sage: f = x^2 + y + x^2*y^2 + 5
1787            sage: f(5,y)
1788            25*y^2 + y + 30
1789            sage: f.fix({x:5})
1790            25*y^2 + y + 30
1791           
1792        """
1793        cdef int mi, i
1794       
1795        cdef MPolynomialRing_libsingular parent = <MPolynomialRing_libsingular>self._parent
1796        cdef ring *_ring = parent._ring
1797
1798        if(_ring != currRing): rChangeCurrRing(_ring)
1799
1800        cdef poly *_p = p_Copy(self._poly, _ring)
1801       
1802        for m,v in fixed.iteritems():
1803            if PY_TYPE_CHECK(m,int) or PY_TYPE_CHECK(m,Integer):
1804                mi = m+1
1805            elif PY_TYPE_CHECK(m,MPolynomial_libsingular) and <MPolynomialRing_libsingular>m.parent() is parent:
1806                for i from 0 < i <= _ring.N:
1807                    mi = p_GetExp((<MPolynomial_libsingular>m)._poly,i, _ring)
1808                    if mi!=0:
1809                        break
1810            else:
1811                raise TypeError, "keys do not match self's parent"
1812                                 
1813            _p = pSubst(p_Copy(_p, _ring), mi, (<MPolynomial_libsingular>parent._coerce_c(v))._poly)
1814
1815        return new_MP(parent,_p)
1816
1817    def monomials(self):
1818        raise NotImplementedError
1819
1820    def constant_coefficent(self):
1821        raise NotImplementedError
1822
1823    def is_univariate(self):
1824        return bool(len(self._variable_indices_(sort=False))<2)
1825
1826    def univariate_polynomial(self, R=None):
1827        raise NotImplementedError
1828
1829    def _variable_indices_(self, sort=True):
1830        cdef poly *p
1831        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
1832        cdef int i
1833        s = set()
1834        p = self._poly
1835        while p:
1836            for i from 1 <= i <= r.N:
1837                if p_GetExp(p,i,r):
1838                    s.add(i)
1839            p = pNext(p)
1840        if sort:
1841            return sorted(s)
1842        else:
1843            return list(s)
1844
1845    def variables(self, sort=True):
1846        cdef poly *p, *v
1847        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
1848        cdef int i
1849        l = list()
1850        si = set()
1851        p = self._poly
1852        while p:
1853            for i from 1 <= i <= r.N:
1854                if i not in si and p_GetExp(p,i,r):
1855                    v = p_ISet(1,r)
1856                    p_SetExp(v, i, 1, r)
1857                    p_Setm(v, r)
1858                    l.append(new_MP(self._parent, v))
1859                    si.add(i)
1860            p = pNext(p)
1861        if sort:
1862            return sorted(l)
1863        else:
1864            return l
1865
1866    def variable(self, i=0):
1867        return self.variables()[i]
1868       
1869    def nvariables(self):
1870        return self._variable_indices_(sort=False)
1871
1872    def is_constant(self):
1873        return bool(p_IsConstant(self._poly, (<MPolynomialRing_libsingular>self._parent)._ring))
1874
1875    cdef int is_constant_c(self):
1876        return p_IsConstant(self._poly, (<MPolynomialRing_libsingular>self._parent)._ring)
1877
1878    def __hash__(self):
1879        s = p_String(self._poly, (<MPolynomialRing_libsingular>self._parent)._ring, (<MPolynomialRing_libsingular>self._parent)._ring)
1880        return hash(s)
1881       
1882
1883    def lm(MPolynomial_libsingular self):
1884        """
1885        Returns the lead monomial of self with respect to the term
1886        order of self.parent(). In SAGE a monomial is a product of
1887        variables in some power without a coefficient.
1888
1889        EXAMPLES:
1890            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
1891
1892            sage: R.<x,y,z>=MPolynomialRing_libsingular(GF(7),3,order='lex')
1893            sage: f = x^1*y^2 + y^3*z^4
1894            sage: f.lm()
1895            x*y^2
1896            sage: f = x^3*y^2*z^4 + x^3*y^2*z^1
1897            sage: f.lm()
1898            x^3*y^2*z^4
1899
1900            sage: R.<x,y,z>=MPolynomialRing_libsingular(QQ,3,order='deglex')
1901            sage: f = x^1*y^2*z^3 + x^3*y^2*z^0
1902            sage: f.lm()
1903            x*y^2*z^3
1904            sage: f = x^1*y^2*z^4 + x^1*y^1*z^5
1905            sage: f.lm()
1906            x*y^2*z^4
1907
1908            sage: R.<x,y,z>=PolynomialRing(GF(127),3,order='degrevlex')
1909            sage: f = x^1*y^5*z^2 + x^4*y^1*z^3
1910            sage: f.lm()
1911            x*y^5*z^2
1912            sage: f = x^4*y^7*z^1 + x^4*y^2*z^3
1913            sage: f.lm()
1914            x^4*y^7*z
1915
1916        """
1917        cdef poly *_p
1918        cdef ring *_ring
1919        _ring = (<MPolynomialRing_libsingular>self._parent)._ring
1920        _p = p_Head(self._poly, _ring)
1921        p_SetCoeff(_p, n_Init(int(1),_ring), _ring)
1922        p_Setm(_p,_ring)
1923        return new_MP((<MPolynomialRing_libsingular>self._parent), _p)
1924
1925
1926    def lc(MPolynomial_libsingular self):
1927        """
1928        Leading coefficient of self. See self.lm() for details.
1929        """
1930
1931        cdef poly *_p
1932        cdef ring *_ring
1933        cdef number *_n
1934        _ring = (<MPolynomialRing_libsingular>self._parent)._ring
1935
1936        if(_ring != currRing): rChangeCurrRing(_ring)
1937       
1938        _p = p_Head(self._poly, _ring)
1939        _n = p_GetCoeff(_p, _ring)
1940
1941        return co.si2sa(_n, _ring, (<MPolynomialRing_libsingular>self._parent)._base)
1942
1943    def lt(MPolynomial_libsingular self):
1944        """
1945        Leading term of self. In SAGE a term is a product of variables
1946        in some power AND a coefficient.
1947
1948        See self.lm() for details
1949        """
1950        return new_MP((<MPolynomialRing_libsingular>self._parent),
1951                                           p_Head(self._poly,(<MPolynomialRing_libsingular>self._parent)._ring))
1952
1953    def is_zero(self):
1954        if self._poly is NULL:
1955            return True
1956        else:
1957            return False
1958
1959    def __nonzero__(self):
1960        if self._poly:
1961            return True
1962        else:
1963            return False
1964
1965    def __floordiv__(self,right):
1966        raise NotImplementedError
1967
1968    def factor(self, param=0):
1969        """
1970        """
1971        cdef ring *_ring
1972        cdef intvec *iv
1973        cdef int *ivv
1974        cdef ideal *I
1975        cdef MPolynomialRing_libsingular parent
1976        cdef int i
1977       
1978        parent = self._parent
1979        _ring = parent._ring
1980
1981        if(_ring != currRing): rChangeCurrRing(_ring)
1982
1983        iv = NULL
1984        I = singclap_factorize ( self._poly, &iv , int(param)) #delete iv at some point
1985
1986        if param!=1:
1987            ivv = iv.ivGetVec()
1988            v = [(new_MP(parent, p_Copy(I.m[i],_ring)) , ivv[i])   for i in range(I.ncols)]
1989        else:
1990            v = [(new_MP(parent, p_Copy(I.m[i],_ring)) , 1)   for i in range(I.ncols)]
1991
1992        # TODO: peel of 1
1993       
1994        F = Factorization(v)
1995        F.sort()
1996
1997        # delete intvec
1998        # delete ideal
1999       
2000        return F
2001
2002    def lift(self, I):
2003        #m = idLift(Ideal(self), I, NULL, FALSE, TRUE );
2004        #matrix_to_list(m)
2005        raise NotImplementedError
2006
2007    def gcd(self, right):
2008        """
2009        Return the greates common divisor of self and right.
2010
2011        INPUT:
2012            right -- polynomial
2013
2014        EXAMPLES:
2015            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
2016            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
2017            sage: f = (x*y*z)^6 - 1
2018            sage: g = (x*y*z)^4 - 1
2019            sage: f.gcd(g)
2020            x^2*y^2*z^2 - 1
2021            sage: GCD([x^3 - 3*x + 2, x^4 - 1, x^6 -1])
2022            x - 1
2023
2024        TESTS:
2025            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
2026            sage: Q.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
2027            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3)
2028            sage: P(0).gcd(Q(0))
2029            0
2030            sage: x.gcd(1)
2031            1
2032
2033        """
2034        cdef MPolynomial_libsingular _right
2035        cdef poly *_res
2036        cdef ring *_ring
2037
2038        _ring = (<MPolynomialRing_libsingular>self._parent)._ring
2039
2040        if(_ring != currRing): rChangeCurrRing(_ring)
2041       
2042        if not PY_TYPE_CHECK(right, MPolynomial_libsingular):
2043            _right = (<MPolynomialRing_libsingular>self._parent)._coerce_c(right)
2044        else:
2045            _right = (<MPolynomial_libsingular>right)
2046
2047        _res = singclap_gcd(p_Copy(self._poly, _ring), p_Copy(_right._poly, _ring))
2048
2049        return new_MP((<MPolynomialRing_libsingular>self._parent), _res)
2050
2051    def lcm(self, g):
2052        #poly *singclap_alglcm ( poly *f, poly *g );
2053        raise NotImplementedError
2054
2055    def is_square_free(self):
2056        #BOOLEAN singclap_isSqrFree(poly f);
2057        raise NotImplementedError
2058
2059    def quo_rem(self, MPolynomial_libsingular right):
2060        if self._parent is not right._parent:
2061            right = self._parent._coerce_c(right)
2062        raise NotImplementedError
2063
2064    def _magma_(self):
2065        raise NotImplementedError
2066
2067    def _singular_(self, singular=singular_default, have_ring=False):
2068        """
2069        Return a SINGULAR (as in the CAS) element for this
2070        element. The result is cached.
2071
2072        INPUT:
2073            singular -- interpreter (default: singular_default)
2074            have_ring -- should the correct ring not be set in SINGULAR first (default:False)
2075
2076        EXAMPLES:
2077            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
2078            sage: P.<x,y,z> = MPolynomialRing_libsingular(GF(127),3)
2079            sage: x._singular_()
2080            x
2081            sage: f =(x^2 + 35*y + 128); f
2082            x^2 + 35*y + 1
2083            sage: x._singular_().name() == x._singular_().name()
2084            True
2085           
2086
2087        TESTS:
2088            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
2089            sage: P.<x,y,z> = MPolynomialRing_libsingular(GF(127),3)
2090            sage: P(0)._singular_()
2091            0
2092
2093        """
2094        if not have_ring:
2095            self.parent()._singular_(singular).set_ring() #this is expensive
2096
2097        try:
2098            if self.__singular is None:
2099                return self._singular_init_c(singular, True)
2100           
2101            self.__singular._check_valid()
2102
2103            if self.__singular.parent() is singular:
2104                return self.__singular
2105           
2106        except (AttributeError, ValueError):
2107            pass
2108
2109        return self._singular_init_c(singular, True)       
2110
2111    def _singular_init_(self,singular=singular_default, have_ring=False):
2112        """
2113        Return a new SINGULAR (as in the CAS) element for this element.
2114
2115        INPUT:
2116            singular -- interpreter (default: singular_default)
2117            have_ring -- should the correct ring not be set in SINGULAR first (default:False)
2118
2119        EXAMPLES:
2120            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
2121            sage: P.<x,y,z> = MPolynomialRing_libsingular(GF(127),3)
2122            sage: x._singular_init_()
2123            x
2124            sage: (x^2+37*y+128)._singular_init_()
2125            x^2+37*y+1
2126            sage: x._singular_init_().name() == x._singular_init_().name()
2127            False
2128
2129        TESTS:
2130            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
2131            sage: P(0)._singular_init_()
2132            0
2133        """
2134        return self._singular_init_c(singular, have_ring)
2135
2136    cdef _singular_init_c(self,singular, have_ring):
2137        """
2138        See MPolynomial_libsingular._singular_init_
2139       
2140        """
2141        if not have_ring:
2142            self.parent()._singular_(singular).set_ring() #this is expensive
2143
2144        self.__singular = singular(str(self))
2145        return self.__singular
2146   
2147    def sub_m_mul_q(self, MPolynomial_libsingular m, MPolynomial_libsingular q):
2148        """
2149        Return self - m*q, where m must be a monomial and q a
2150        polynomial.
2151
2152        INPUT:
2153            m -- a monomial
2154            q -- a polynomial
2155
2156        EXAMPLE:
2157            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
2158            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
2159            sage: x.sub_m_mul_q(y,z)
2160            -y*z + x
2161
2162        TESTS:
2163            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
2164            sage: Q.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
2165            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
2166            sage: P(0).sub_m_mul_q(P(0),P(1))
2167            0
2168            sage: x.sub_m_mul_q(Q.gen(1),Q.gen(2))
2169            -y*z + x
2170
2171         """
2172        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
2173
2174        if not self._parent is m._parent:
2175            m = self._parent._coerce_c(m)
2176        if not self._parent is q._parent:
2177            q = self._parent._coerce_c(q)
2178
2179        if m._poly and m._poly.next:
2180            raise ArithmeticError, "m must be a monomial"
2181        elif not m._poly:
2182            return self
2183
2184        return new_MP(self._parent, p_Minus_mm_Mult_qq(p_Copy(self._poly, r), m._poly, q._poly, r))
2185
2186    def add_m_mul_q(self, MPolynomial_libsingular m, MPolynomial_libsingular q):
2187        """
2188        Return self + m*q, where m must be a monomial and q a
2189        polynomial.
2190
2191       INPUT:
2192            m -- a monomial
2193            q -- a polynomial
2194
2195        EXAMPLE:
2196            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
2197            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
2198            sage: x.add_m_mul_q(y,z)
2199            y*z + x
2200
2201        TESTS:
2202            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular       
2203            sage: R.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
2204            sage: P.<x,y,z>=MPolynomialRing_libsingular(QQ,3)
2205            sage: P(0).add_m_mul_q(P(0),P(1))
2206            0
2207            sage: x.add_m_mul_q(R.gen(),R.gen(1))
2208            x*y + x
2209         """
2210
2211        cdef ring *r = (<MPolynomialRing_libsingular>self._parent)._ring
2212
2213        if not self._parent is m._parent:
2214            m = self._parent._coerce_c(m)
2215        if not self._parent is q._parent:
2216            q = self._parent._coerce_c(q)
2217
2218        if m._poly and m._poly.next:
2219            raise ArithmeticError, "m must be a monomial"
2220        elif not m._poly:
2221            return self
2222
2223        return new_MP(self._parent, p_Plus_mm_Mult_qq(p_Copy(self._poly, r), m._poly, q._poly, r))
2224       
2225
2226    def __reduce__(self):
2227        """
2228
2229        Serialize self.
2230
2231        EXAMPLES:
2232            sage: from sage.rings.multi_polynomial_libsingular import MPolynomialRing_libsingular
2233            sage: P.<x,y,z> = MPolynomialRing_libsingular(QQ,3, order='degrevlex')
2234            sage: f = 27/113 * x^2 + y*z + 1/2
2235            sage: f == loads(dumps(f))
2236            True
2237           
2238            sage: P = MPolynomialRing_libsingular(GF(127),3,names='abc')
2239            sage: a,b,c = P.gens()
2240            sage: f = 57 * a^2*b + 43 * c + 1
2241            sage: f == loads(dumps(f))
2242            True
2243
2244        """
2245        return sage.rings.multi_polynomial_libsingular.unpickle_MPolynomial_libsingular, ( self._parent, self.dict() )
2246
2247def unpickle_MPolynomial_libsingular(MPolynomialRing_libsingular R, d):
2248    """
2249    Deserialize a MPolynomial_libsingular object
2250
2251    INPUT:
2252        R -- the base ring
2253        d -- a Python dictionary as returned by MPolynomial_libsingular.dict
2254   
2255    """
2256    cdef ring *r = R._ring
2257    cdef poly *m, *p
2258    cdef int _i, _e
2259    p = p_ISet(0,r)
2260    for mon,c in d.iteritems():
2261        m = p_Init(r)
2262        for i,e in mon.sparse_iter():
2263            _i = i
2264            if _i >= r.N:
2265                p_Delete(&p,r)
2266                p_Delete(&m,r)
2267                raise TypeError, "variable index too big"
2268            _e = e
2269            if _e <= 0:
2270                p_Delete(&p,r)
2271                p_Delete(&m,r)
2272                raise TypeError, "exponent too small"
2273            p_SetExp(m, _i+1,_e, r)
2274        p_SetCoeff(m, co.sa2si(c, r), r)
2275        p_Setm(m,r)
2276        p = p_Add_q(p,m,r)
2277    return new_MP(R,p)
Note: See TracBrowser for help on using the repository browser.