Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/libs/singular/groebner_strategy.pyx
4085 views
1
"""
2
Singular's Groebner Strategy Objects
3
4
AUTHORS:
5
6
- Martin Albrecht (2009-07): initial implementation
7
- Michael Brickenstein (2009-07): initial implementation
8
- Hans Schoenemann (2009-07): initial implementation
9
"""
10
11
#*****************************************************************************
12
# Copyright (C) 2009 Martin Albrecht <[email protected]>
13
# Copyright (C) 2009 Michael Brickenstein <[email protected]>
14
# Copyright (C) 2009 Hans Schoenemann <[email protected]>
15
#
16
# Distributed under the terms of the GNU General Public License (GPL)
17
# http://www.gnu.org/licenses/
18
#*****************************************************************************
19
20
cdef extern from "":
21
int unlikely(int)
22
int likely(int)
23
24
from sage.libs.singular.decl cimport ideal, ring, poly, currRing
25
from sage.libs.singular.decl cimport rChangeCurrRing
26
from sage.libs.singular.decl cimport new_skStrategy, delete_skStrategy, idRankFreeModule
27
from sage.libs.singular.decl cimport initEcartBBA, enterSBba, initBuchMoraCrit, initS, pNorm, id_Delete, kTest
28
from sage.libs.singular.decl cimport omfree, redNF, p_Copy, redtailBba
29
30
from sage.libs.singular.ring cimport singular_ring_reference, singular_ring_delete
31
32
from sage.rings.polynomial.multi_polynomial_ideal import MPolynomialIdeal, NCPolynomialIdeal
33
from sage.rings.polynomial.multi_polynomial_ideal_libsingular cimport sage_ideal_to_singular_ideal
34
from sage.rings.polynomial.multi_polynomial_libsingular cimport MPolynomial_libsingular, MPolynomialRing_libsingular, new_MP
35
from sage.rings.polynomial.plural cimport new_NCP
36
37
cdef class GroebnerStrategy(SageObject):
38
"""
39
A Wrapper for Singular's Groebner Strategy Object.
40
41
This object provides functions for normal form computations and
42
other functions for Groebner basis computation.
43
44
ALGORITHM:
45
46
Uses Singular via libSINGULAR
47
"""
48
def __cinit__(self, L):
49
"""
50
Create a new :class:`GroebnerStrategy` object for the
51
generators of the ideal ``L``.
52
53
INPUT:
54
55
- ``L`` - a multivariate polynomial ideal
56
57
EXAMPLES::
58
59
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
60
sage: P.<x,y,z> = PolynomialRing(QQ)
61
sage: I = Ideal([x+z,y+z+1])
62
sage: strat = GroebnerStrategy(I); strat
63
Groebner Strategy for ideal generated by 2 elements
64
over Multivariate Polynomial Ring in x, y, z over Rational Field
65
66
TESTS::
67
68
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
69
sage: strat = GroebnerStrategy(None)
70
Traceback (most recent call last):
71
...
72
TypeError: First parameter must be a multivariate polynomial ideal.
73
74
sage: P.<x,y,z> = PolynomialRing(QQ,order='neglex')
75
sage: I = Ideal([x+z,y+z+1])
76
sage: strat = GroebnerStrategy(I)
77
Traceback (most recent call last):
78
...
79
NotImplementedError: The local case is not implemented yet.
80
81
sage: P.<x,y,z> = PolynomialRing(CC,order='neglex')
82
sage: I = Ideal([x+z,y+z+1])
83
sage: strat = GroebnerStrategy(I)
84
Traceback (most recent call last):
85
...
86
TypeError: First parameter's ring must be multivariate polynomial ring via libsingular.
87
88
sage: P.<x,y,z> = PolynomialRing(ZZ)
89
sage: I = Ideal([x+z,y+z+1])
90
sage: strat = GroebnerStrategy(I)
91
Traceback (most recent call last):
92
...
93
NotImplementedError: Only coefficient fields are implemented so far.
94
95
"""
96
if not isinstance(L, MPolynomialIdeal):
97
raise TypeError("First parameter must be a multivariate polynomial ideal.")
98
99
if not isinstance(L.ring(), MPolynomialRing_libsingular):
100
raise TypeError("First parameter's ring must be multivariate polynomial ring via libsingular.")
101
102
self._ideal = L
103
104
cdef MPolynomialRing_libsingular R = <MPolynomialRing_libsingular>L.ring()
105
self._parent = R
106
self._parent_ring = singular_ring_reference(R._ring)
107
108
if not R.term_order().is_global():
109
raise NotImplementedError("The local case is not implemented yet.")
110
111
if not R.base_ring().is_field():
112
raise NotImplementedError("Only coefficient fields are implemented so far.")
113
114
if (R._ring != currRing):
115
rChangeCurrRing(R._ring)
116
117
cdef ideal *i = sage_ideal_to_singular_ideal(L)
118
self._strat = new_skStrategy()
119
120
self._strat.ak = idRankFreeModule(i, R._ring)
121
#- creating temp data structures
122
initBuchMoraCrit(self._strat)
123
self._strat.initEcart = initEcartBBA
124
self._strat.enterS = enterSBba
125
#- set S
126
self._strat.sl = -1
127
#- init local data struct
128
initS(i, NULL, self._strat)
129
130
cdef int j
131
cdef bint base_ring_is_field = R.base_ring().is_field()
132
if (R._ring != currRing): rChangeCurrRing(R._ring)
133
if base_ring_is_field:
134
for j in range(self._strat.sl+1)[::-1]:
135
pNorm(self._strat.S[j])
136
137
id_Delete(&i, R._ring)
138
kTest(self._strat)
139
140
def __dealloc__(self):
141
"""
142
TEST::
143
144
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
145
sage: P.<x,y,z> = PolynomialRing(GF(32003))
146
sage: I = Ideal([x + z, y + z])
147
sage: strat = GroebnerStrategy(I)
148
sage: del strat
149
"""
150
# WARNING: the Cython class self._parent is no longer accessible!
151
# see http://trac.sagemath.org/sage_trac/ticket/11339
152
cdef ring *oldRing = NULL
153
if self._strat:
154
omfree(self._strat.sevS)
155
omfree(self._strat.ecartS)
156
omfree(self._strat.T)
157
omfree(self._strat.sevT)
158
omfree(self._strat.R)
159
omfree(self._strat.S_2_R)
160
omfree(self._strat.L)
161
omfree(self._strat.B)
162
omfree(self._strat.fromQ)
163
id_Delete(&self._strat.Shdl, self._parent_ring)
164
165
if self._parent_ring != currRing:
166
oldRing = currRing
167
rChangeCurrRing(self._parent_ring)
168
delete_skStrategy(self._strat)
169
rChangeCurrRing(oldRing)
170
else:
171
delete_skStrategy(self._strat)
172
singular_ring_delete(self._parent_ring)
173
174
def _repr_(self):
175
"""
176
TEST::
177
178
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
179
sage: P.<x,y,z> = PolynomialRing(GF(32003))
180
sage: I = Ideal([x + z, y + z])
181
sage: strat = GroebnerStrategy(I)
182
sage: strat # indirect doctest
183
Groebner Strategy for ideal generated by 2 elements over
184
Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003
185
"""
186
return "Groebner Strategy for ideal generated by %d elements over %s"%(self._ideal.ngens(),self._parent)
187
188
def ideal(self):
189
"""
190
Return the ideal this strategy object is defined for.
191
192
EXAMPLE::
193
194
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
195
sage: P.<x,y,z> = PolynomialRing(GF(32003))
196
sage: I = Ideal([x + z, y + z])
197
sage: strat = GroebnerStrategy(I)
198
sage: strat.ideal()
199
Ideal (x + z, y + z) of Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003
200
"""
201
return self._ideal
202
203
def ring(self):
204
"""
205
Return the ring this strategy object is defined over.
206
207
EXAMPLE::
208
209
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
210
sage: P.<x,y,z> = PolynomialRing(GF(32003))
211
sage: I = Ideal([x + z, y + z])
212
sage: strat = GroebnerStrategy(I)
213
sage: strat.ring()
214
Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003
215
"""
216
return self._parent
217
218
def __cmp__(self, other):
219
"""
220
EXAMPLE::
221
222
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
223
sage: P.<x,y,z> = PolynomialRing(GF(19))
224
sage: I = Ideal([P(0)])
225
sage: strat = GroebnerStrategy(I)
226
sage: strat == GroebnerStrategy(I)
227
True
228
sage: I = Ideal([x+1,y+z])
229
sage: strat == GroebnerStrategy(I)
230
False
231
"""
232
if not isinstance(other, GroebnerStrategy):
233
return cmp(type(self),other(type))
234
else:
235
return cmp(self._ideal.gens(),(<GroebnerStrategy>other)._ideal.gens())
236
237
def __reduce__(self):
238
"""
239
EXAMPLE::
240
241
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
242
sage: P.<x,y,z> = PolynomialRing(GF(32003))
243
sage: I = Ideal([x + z, y + z])
244
sage: strat = GroebnerStrategy(I)
245
sage: loads(dumps(strat)) == strat
246
True
247
"""
248
return unpickle_GroebnerStrategy0, (self._ideal,)
249
250
cpdef MPolynomial_libsingular normal_form(self, MPolynomial_libsingular p):
251
"""
252
Compute the normal form of ``p`` with respect to the
253
generators of this object.
254
255
EXAMPLE::
256
257
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
258
sage: P.<x,y,z> = PolynomialRing(QQ)
259
sage: I = Ideal([x + z, y + z])
260
sage: strat = GroebnerStrategy(I)
261
sage: strat.normal_form(x*y) # indirect doctest
262
z^2
263
sage: strat.normal_form(x + 1)
264
-z + 1
265
266
TESTS::
267
268
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
269
sage: P.<x,y,z> = PolynomialRing(QQ)
270
sage: I = Ideal([P(0)])
271
sage: strat = GroebnerStrategy(I)
272
sage: strat.normal_form(x)
273
x
274
sage: strat.normal_form(P(0))
275
0
276
"""
277
if unlikely(p._parent is not self._parent):
278
raise TypeError("parent(p) must be the same as this object's parent.")
279
if unlikely(self._parent._ring != currRing):
280
rChangeCurrRing(self._parent._ring)
281
282
cdef int max_ind
283
cdef poly *_p = redNF(p_Copy(p._poly, self._parent._ring), max_ind, 0, self._strat)
284
if likely(_p!=NULL):
285
_p = redtailBba(_p, max_ind, self._strat)
286
return new_MP(self._parent, _p)
287
288
cdef class NCGroebnerStrategy(SageObject):
289
"""
290
A Wrapper for Singular's Groebner Strategy Object.
291
292
This object provides functions for normal form computations and
293
other functions for Groebner basis computation.
294
295
ALGORITHM:
296
297
Uses Singular via libSINGULAR
298
"""
299
def __init__(self, L):
300
"""
301
Create a new :class:`GroebnerStrategy` object for the
302
generators of the ideal ``L``.
303
304
INPUT:
305
306
- ``L`` - an ideal in a g-algebra
307
308
EXAMPLES::
309
310
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
311
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
312
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
313
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
314
sage: NCGroebnerStrategy(I)
315
Groebner Strategy for ideal generated by 3 elements over Noncommutative Multivariate Polynomial Ring in x, y, z over Rational Field, nc-relations: {y*x: x*y - z, z*y: y*z - 2*y, z*x: x*z + 2*x}
316
317
TESTS::
318
319
sage: strat = NCGroebnerStrategy(None)
320
Traceback (most recent call last):
321
...
322
TypeError: First parameter must be an ideal in a g-algebra.
323
324
sage: P.<x,y,z> = PolynomialRing(CC,order='neglex')
325
sage: I = Ideal([x+z,y+z+1])
326
sage: strat = NCGroebnerStrategy(I)
327
Traceback (most recent call last):
328
...
329
TypeError: First parameter must be an ideal in a g-algebra.
330
331
"""
332
if not isinstance(L, NCPolynomialIdeal):
333
raise TypeError("First parameter must be an ideal in a g-algebra.")
334
335
if not isinstance(L.ring(), NCPolynomialRing_plural):
336
raise TypeError("First parameter's ring must be a g-algebra.")
337
338
self._ideal = L
339
340
cdef NCPolynomialRing_plural R = <NCPolynomialRing_plural>L.ring()
341
self._parent = R
342
343
if not R.term_order().is_global():
344
raise NotImplementedError("The local case is not implemented yet.")
345
346
if not R.base_ring().is_field():
347
raise NotImplementedError("Only coefficient fields are implemented so far.")
348
349
if (R._ring != currRing):
350
rChangeCurrRing(R._ring)
351
352
cdef ideal *i = sage_ideal_to_singular_ideal(L)
353
self._strat = new_skStrategy()
354
355
self._strat.ak = idRankFreeModule(i, R._ring)
356
#- creating temp data structures
357
initBuchMoraCrit(self._strat)
358
self._strat.initEcart = initEcartBBA
359
self._strat.enterS = enterSBba
360
#- set S
361
self._strat.sl = -1
362
#- init local data struct
363
initS(i, NULL, self._strat)
364
365
cdef int j
366
if R.base_ring().is_field():
367
for j in range(self._strat.sl+1)[::-1]:
368
pNorm(self._strat.S[j])
369
370
id_Delete(&i, R._ring)
371
kTest(self._strat)
372
373
def __dealloc__(self):
374
"""
375
TEST::
376
377
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
378
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
379
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
380
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
381
sage: strat = NCGroebnerStrategy(I)
382
sage: del strat # indirect doctest
383
"""
384
cdef ring *oldRing = NULL
385
if self._strat:
386
omfree(self._strat.sevS)
387
omfree(self._strat.ecartS)
388
omfree(self._strat.T)
389
omfree(self._strat.sevT)
390
omfree(self._strat.R)
391
omfree(self._strat.S_2_R)
392
omfree(self._strat.L)
393
omfree(self._strat.B)
394
omfree(self._strat.fromQ)
395
id_Delete(&self._strat.Shdl, self._parent._ring)
396
397
if self._parent._ring != currRing:
398
oldRing = currRing
399
rChangeCurrRing(self._parent._ring)
400
delete_skStrategy(self._strat)
401
rChangeCurrRing(oldRing)
402
else:
403
delete_skStrategy(self._strat)
404
405
def _repr_(self):
406
"""
407
TEST::
408
409
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
410
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
411
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
412
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
413
sage: strat = NCGroebnerStrategy(I)
414
sage: strat # indirect doctest
415
Groebner Strategy for ideal generated by 3 elements over Noncommutative Multivariate Polynomial Ring in x, y, z over Rational Field, nc-relations: {y*x: x*y - z, z*y: y*z - 2*y, z*x: x*z + 2*x}
416
"""
417
return "Groebner Strategy for ideal generated by %d elements over %s"%(self._ideal.ngens(),self._parent)
418
419
def ideal(self):
420
"""
421
Return the ideal this strategy object is defined for.
422
423
EXAMPLE::
424
425
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
426
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
427
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
428
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
429
sage: strat = NCGroebnerStrategy(I)
430
sage: strat.ideal() == I
431
True
432
433
"""
434
return self._ideal
435
436
def ring(self):
437
"""
438
Return the ring this strategy object is defined over.
439
440
EXAMPLE::
441
442
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
443
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
444
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
445
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
446
sage: strat = NCGroebnerStrategy(I)
447
sage: strat.ring() is H
448
True
449
"""
450
return self._parent
451
452
def __cmp__(self, other):
453
"""
454
EXAMPLE::
455
456
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
457
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
458
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
459
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
460
sage: strat = NCGroebnerStrategy(I)
461
sage: strat == NCGroebnerStrategy(I)
462
True
463
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()], side='twosided')
464
sage: strat == NCGroebnerStrategy(I)
465
False
466
"""
467
if not isinstance(other, NCGroebnerStrategy):
468
return cmp(type(self),other(type))
469
else:
470
return cmp((self._ideal.gens(),self._ideal.side()),
471
((<NCGroebnerStrategy>other)._ideal.gens(),
472
(<NCGroebnerStrategy>other)._ideal.side()))
473
474
def __reduce__(self):
475
"""
476
EXAMPLE::
477
478
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
479
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
480
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
481
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
482
sage: strat = NCGroebnerStrategy(I)
483
sage: loads(dumps(strat)) == strat
484
True
485
"""
486
return unpickle_NCGroebnerStrategy0, (self._ideal,)
487
488
cpdef NCPolynomial_plural normal_form(self, NCPolynomial_plural p):
489
"""
490
Compute the normal form of ``p`` with respect to the
491
generators of this object.
492
493
EXAMPLE::
494
495
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
496
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
497
sage: JL = H.ideal([x^3, y^3, z^3 - 4*z])
498
sage: JT = H.ideal([x^3, y^3, z^3 - 4*z], side='twosided')
499
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
500
sage: SL = NCGroebnerStrategy(JL.std())
501
sage: ST = NCGroebnerStrategy(JT.std())
502
sage: SL.normal_form(x*y^2)
503
x*y^2
504
sage: ST.normal_form(x*y^2)
505
y*z
506
507
"""
508
if unlikely(p._parent is not self._parent):
509
raise TypeError("parent(p) must be the same as this object's parent.")
510
if unlikely(self._parent._ring != currRing):
511
rChangeCurrRing(self._parent._ring)
512
513
cdef int max_ind
514
cdef poly *_p = redNF(p_Copy(p._poly, self._parent._ring), max_ind, 0, self._strat)
515
if likely(_p!=NULL):
516
_p = redtailBba(_p, max_ind, self._strat)
517
return new_NCP(self._parent, _p)
518
519
def unpickle_NCGroebnerStrategy0(I):
520
"""
521
EXAMPLE::
522
523
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
524
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
525
sage: H.<x,y,z> = A.g_algebra({y*x:x*y-z, z*x:x*z+2*x, z*y:y*z-2*y})
526
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
527
sage: strat = NCGroebnerStrategy(I)
528
sage: loads(dumps(strat)) == strat # indirect doctest
529
True
530
"""
531
return NCGroebnerStrategy(I)
532
533
def unpickle_GroebnerStrategy0(I):
534
"""
535
EXAMPLE::
536
537
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
538
sage: P.<x,y,z> = PolynomialRing(GF(32003))
539
sage: I = Ideal([x + z, y + z])
540
sage: strat = GroebnerStrategy(I)
541
sage: loads(dumps(strat)) == strat # indirect doctest
542
True
543
"""
544
return GroebnerStrategy(I)
545
546