Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/libs/singular/groebner_strategy.pyx
8815 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 *: # hack to get at cython macro
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
if self._parent_ring:
173
singular_ring_delete(self._parent_ring)
174
175
def _repr_(self):
176
"""
177
TEST::
178
179
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
180
sage: P.<x,y,z> = PolynomialRing(GF(32003))
181
sage: I = Ideal([x + z, y + z])
182
sage: strat = GroebnerStrategy(I)
183
sage: strat # indirect doctest
184
Groebner Strategy for ideal generated by 2 elements over
185
Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003
186
"""
187
return "Groebner Strategy for ideal generated by %d elements over %s"%(self._ideal.ngens(),self._parent)
188
189
def ideal(self):
190
"""
191
Return the ideal this strategy object is defined for.
192
193
EXAMPLE::
194
195
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
196
sage: P.<x,y,z> = PolynomialRing(GF(32003))
197
sage: I = Ideal([x + z, y + z])
198
sage: strat = GroebnerStrategy(I)
199
sage: strat.ideal()
200
Ideal (x + z, y + z) of Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003
201
"""
202
return self._ideal
203
204
def ring(self):
205
"""
206
Return the ring this strategy object is defined over.
207
208
EXAMPLE::
209
210
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
211
sage: P.<x,y,z> = PolynomialRing(GF(32003))
212
sage: I = Ideal([x + z, y + z])
213
sage: strat = GroebnerStrategy(I)
214
sage: strat.ring()
215
Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003
216
"""
217
return self._parent
218
219
def __cmp__(self, other):
220
"""
221
EXAMPLE::
222
223
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
224
sage: P.<x,y,z> = PolynomialRing(GF(19))
225
sage: I = Ideal([P(0)])
226
sage: strat = GroebnerStrategy(I)
227
sage: strat == GroebnerStrategy(I)
228
True
229
sage: I = Ideal([x+1,y+z])
230
sage: strat == GroebnerStrategy(I)
231
False
232
"""
233
if not isinstance(other, GroebnerStrategy):
234
return cmp(type(self),other(type))
235
else:
236
return cmp(self._ideal.gens(),(<GroebnerStrategy>other)._ideal.gens())
237
238
def __reduce__(self):
239
"""
240
EXAMPLE::
241
242
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
243
sage: P.<x,y,z> = PolynomialRing(GF(32003))
244
sage: I = Ideal([x + z, y + z])
245
sage: strat = GroebnerStrategy(I)
246
sage: loads(dumps(strat)) == strat
247
True
248
"""
249
return unpickle_GroebnerStrategy0, (self._ideal,)
250
251
cpdef MPolynomial_libsingular normal_form(self, MPolynomial_libsingular p):
252
"""
253
Compute the normal form of ``p`` with respect to the
254
generators of this object.
255
256
EXAMPLE::
257
258
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
259
sage: P.<x,y,z> = PolynomialRing(QQ)
260
sage: I = Ideal([x + z, y + z])
261
sage: strat = GroebnerStrategy(I)
262
sage: strat.normal_form(x*y) # indirect doctest
263
z^2
264
sage: strat.normal_form(x + 1)
265
-z + 1
266
267
TESTS::
268
269
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
270
sage: P.<x,y,z> = PolynomialRing(QQ)
271
sage: I = Ideal([P(0)])
272
sage: strat = GroebnerStrategy(I)
273
sage: strat.normal_form(x)
274
x
275
sage: strat.normal_form(P(0))
276
0
277
"""
278
if unlikely(p._parent is not self._parent):
279
raise TypeError("parent(p) must be the same as this object's parent.")
280
if unlikely(self._parent._ring != currRing):
281
rChangeCurrRing(self._parent._ring)
282
283
cdef int max_ind
284
cdef poly *_p = redNF(p_Copy(p._poly, self._parent._ring), max_ind, 0, self._strat)
285
if likely(_p!=NULL):
286
_p = redtailBba(_p, max_ind, self._strat)
287
return new_MP(self._parent, _p)
288
289
cdef class NCGroebnerStrategy(SageObject):
290
"""
291
A Wrapper for Singular's Groebner Strategy Object.
292
293
This object provides functions for normal form computations and
294
other functions for Groebner basis computation.
295
296
ALGORITHM:
297
298
Uses Singular via libSINGULAR
299
"""
300
def __init__(self, L):
301
"""
302
Create a new :class:`GroebnerStrategy` object for the
303
generators of the ideal ``L``.
304
305
INPUT:
306
307
- ``L`` - an ideal in a g-algebra
308
309
EXAMPLES::
310
311
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
312
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
313
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})
314
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
315
sage: NCGroebnerStrategy(I)
316
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}
317
318
TESTS::
319
320
sage: strat = NCGroebnerStrategy(None)
321
Traceback (most recent call last):
322
...
323
TypeError: First parameter must be an ideal in a g-algebra.
324
325
sage: P.<x,y,z> = PolynomialRing(CC,order='neglex')
326
sage: I = Ideal([x+z,y+z+1])
327
sage: strat = NCGroebnerStrategy(I)
328
Traceback (most recent call last):
329
...
330
TypeError: First parameter must be an ideal in a g-algebra.
331
332
"""
333
if not isinstance(L, NCPolynomialIdeal):
334
raise TypeError("First parameter must be an ideal in a g-algebra.")
335
336
if not isinstance(L.ring(), NCPolynomialRing_plural):
337
raise TypeError("First parameter's ring must be a g-algebra.")
338
339
self._ideal = L
340
341
cdef NCPolynomialRing_plural R = <NCPolynomialRing_plural>L.ring()
342
self._parent = R
343
344
if not R.term_order().is_global():
345
raise NotImplementedError("The local case is not implemented yet.")
346
347
if not R.base_ring().is_field():
348
raise NotImplementedError("Only coefficient fields are implemented so far.")
349
350
if (R._ring != currRing):
351
rChangeCurrRing(R._ring)
352
353
cdef ideal *i = sage_ideal_to_singular_ideal(L)
354
self._strat = new_skStrategy()
355
356
self._strat.ak = idRankFreeModule(i, R._ring)
357
#- creating temp data structures
358
initBuchMoraCrit(self._strat)
359
self._strat.initEcart = initEcartBBA
360
self._strat.enterS = enterSBba
361
#- set S
362
self._strat.sl = -1
363
#- init local data struct
364
initS(i, NULL, self._strat)
365
366
cdef int j
367
if R.base_ring().is_field():
368
for j in range(self._strat.sl+1)[::-1]:
369
pNorm(self._strat.S[j])
370
371
id_Delete(&i, R._ring)
372
kTest(self._strat)
373
374
def __dealloc__(self):
375
"""
376
TEST::
377
378
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
379
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
380
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})
381
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
382
sage: strat = NCGroebnerStrategy(I)
383
sage: del strat # indirect doctest
384
"""
385
cdef ring *oldRing = NULL
386
if self._strat:
387
omfree(self._strat.sevS)
388
omfree(self._strat.ecartS)
389
omfree(self._strat.T)
390
omfree(self._strat.sevT)
391
omfree(self._strat.R)
392
omfree(self._strat.S_2_R)
393
omfree(self._strat.L)
394
omfree(self._strat.B)
395
omfree(self._strat.fromQ)
396
id_Delete(&self._strat.Shdl, self._parent._ring)
397
398
if self._parent._ring != currRing:
399
oldRing = currRing
400
rChangeCurrRing(self._parent._ring)
401
delete_skStrategy(self._strat)
402
rChangeCurrRing(oldRing)
403
else:
404
delete_skStrategy(self._strat)
405
406
def _repr_(self):
407
"""
408
TEST::
409
410
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
411
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
412
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})
413
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
414
sage: strat = NCGroebnerStrategy(I)
415
sage: strat # indirect doctest
416
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}
417
"""
418
return "Groebner Strategy for ideal generated by %d elements over %s"%(self._ideal.ngens(),self._parent)
419
420
def ideal(self):
421
"""
422
Return the ideal this strategy object is defined for.
423
424
EXAMPLE::
425
426
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
427
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
428
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})
429
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
430
sage: strat = NCGroebnerStrategy(I)
431
sage: strat.ideal() == I
432
True
433
434
"""
435
return self._ideal
436
437
def ring(self):
438
"""
439
Return the ring this strategy object is defined over.
440
441
EXAMPLE::
442
443
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
444
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
445
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})
446
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
447
sage: strat = NCGroebnerStrategy(I)
448
sage: strat.ring() is H
449
True
450
"""
451
return self._parent
452
453
def __cmp__(self, other):
454
"""
455
EXAMPLE::
456
457
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
458
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
459
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})
460
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
461
sage: strat = NCGroebnerStrategy(I)
462
sage: strat == NCGroebnerStrategy(I)
463
True
464
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()], side='twosided')
465
sage: strat == NCGroebnerStrategy(I)
466
False
467
"""
468
if not isinstance(other, NCGroebnerStrategy):
469
return cmp(type(self),other(type))
470
else:
471
return cmp((self._ideal.gens(),self._ideal.side()),
472
((<NCGroebnerStrategy>other)._ideal.gens(),
473
(<NCGroebnerStrategy>other)._ideal.side()))
474
475
def __reduce__(self):
476
"""
477
EXAMPLE::
478
479
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
480
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
481
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})
482
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
483
sage: strat = NCGroebnerStrategy(I)
484
sage: loads(dumps(strat)) == strat
485
True
486
"""
487
return unpickle_NCGroebnerStrategy0, (self._ideal,)
488
489
cpdef NCPolynomial_plural normal_form(self, NCPolynomial_plural p):
490
"""
491
Compute the normal form of ``p`` with respect to the
492
generators of this object.
493
494
EXAMPLE::
495
496
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
497
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})
498
sage: JL = H.ideal([x^3, y^3, z^3 - 4*z])
499
sage: JT = H.ideal([x^3, y^3, z^3 - 4*z], side='twosided')
500
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
501
sage: SL = NCGroebnerStrategy(JL.std())
502
sage: ST = NCGroebnerStrategy(JT.std())
503
sage: SL.normal_form(x*y^2)
504
x*y^2
505
sage: ST.normal_form(x*y^2)
506
y*z
507
508
"""
509
if unlikely(p._parent is not self._parent):
510
raise TypeError("parent(p) must be the same as this object's parent.")
511
if unlikely(self._parent._ring != currRing):
512
rChangeCurrRing(self._parent._ring)
513
514
cdef int max_ind
515
cdef poly *_p = redNF(p_Copy(p._poly, self._parent._ring), max_ind, 0, self._strat)
516
if likely(_p!=NULL):
517
_p = redtailBba(_p, max_ind, self._strat)
518
return new_NCP(self._parent, _p)
519
520
def unpickle_NCGroebnerStrategy0(I):
521
"""
522
EXAMPLE::
523
524
sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy
525
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
526
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})
527
sage: I = H.ideal([y^2, x^2, z^2-H.one_element()])
528
sage: strat = NCGroebnerStrategy(I)
529
sage: loads(dumps(strat)) == strat # indirect doctest
530
True
531
"""
532
return NCGroebnerStrategy(I)
533
534
def unpickle_GroebnerStrategy0(I):
535
"""
536
EXAMPLE::
537
538
sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy
539
sage: P.<x,y,z> = PolynomialRing(GF(32003))
540
sage: I = Ideal([x + z, y + z])
541
sage: strat = GroebnerStrategy(I)
542
sage: loads(dumps(strat)) == strat # indirect doctest
543
True
544
"""
545
return GroebnerStrategy(I)
546
547