Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

""" 

Direct low-level access to SINGULAR's Groebner basis engine via libSINGULAR. 

  

AUTHOR: 

  

- Martin Albrecht (2007-08-08): initial version 

  

EXAMPLES:: 

  

sage: x,y,z = QQ['x,y,z'].gens() 

sage: I = ideal(x^5 + y^4 + z^3 - 1, x^3 + y^3 + z^2 - 1) 

sage: I.groebner_basis('libsingular:std') 

[y^6 + x*y^4 + 2*y^3*z^2 + x*z^3 + z^4 - 2*y^3 - 2*z^2 - x + 1, 

x^2*y^3 - y^4 + x^2*z^2 - z^3 - x^2 + 1, x^3 + y^3 + z^2 - 1] 

  

We compute a Groebner basis for cyclic 6, which is a standard 

benchmark and test ideal:: 

  

sage: R.<x,y,z,t,u,v> = QQ['x,y,z,t,u,v'] 

sage: I = sage.rings.ideal.Cyclic(R,6) 

sage: B = I.groebner_basis('libsingular:std') 

sage: len(B) 

45 

  

Two examples from the Mathematica documentation (done in Sage): 

  

- We compute a Groebner basis:: 

  

sage: R.<x,y> = PolynomialRing(QQ, order='lex') 

sage: ideal(x^2 - 2*y^2, x*y - 3).groebner_basis('libsingular:slimgb') 

[x - 2/3*y^3, y^4 - 9/2] 

  

- We show that three polynomials have no common root:: 

  

sage: R.<x,y> = QQ[] 

sage: ideal(x+y, x^2 - 1, y^2 - 2*x).groebner_basis('libsingular:slimgb') 

[1] 

""" 

  

#***************************************************************************** 

# Copyright (C) 2007 Martin Albrecht <malb@informatik.uni-bremen.de> 

# Copyright (C) 2007 William Stein <wstein@gmail.com> 

# 

# This program is free software: you can redistribute it and/or modify 

# it under the terms of the GNU General Public License as published by 

# the Free Software Foundation, either version 2 of the License, or 

# (at your option) any later version. 

# http://www.gnu.org/licenses/ 

#***************************************************************************** 

  

from cysignals.signals cimport sig_on, sig_off 

  

from sage.libs.singular.decl cimport tHomog, number, IDELEMS, p_Copy, rChangeCurrRing 

from sage.libs.singular.decl cimport idInit, id_Delete, currRing, Sy_bit, OPT_REDSB 

from sage.libs.singular.decl cimport scKBase, poly, testHomog, idSkipZeroes, id_RankFreeModule, kStd 

from sage.libs.singular.decl cimport OPT_REDTAIL, singular_options, kInterRed, t_rep_gb, p_GetCoeff 

from sage.libs.singular.decl cimport pp_Mult_nn, p_Delete, n_Delete 

from sage.libs.singular.decl cimport rIsPluralRing 

from sage.libs.singular.decl cimport n_unknown, n_Zp, n_Q, n_R, n_GF, n_long_R, n_algExt,n_transExt,n_long_C, n_Z, n_Zn, n_Znm, n_Z2m, n_CF 

  

from sage.rings.polynomial.multi_polynomial_libsingular cimport new_MP 

from sage.rings.polynomial.plural cimport new_NCP 

  

from sage.rings.polynomial.multi_polynomial_ideal import MPolynomialIdeal 

from sage.rings.polynomial.multi_polynomial_libsingular cimport MPolynomial_libsingular 

from sage.rings.polynomial.multi_polynomial_libsingular cimport MPolynomialRing_libsingular 

from sage.structure.sequence import Sequence 

  

from sage.rings.polynomial.plural cimport NCPolynomialRing_plural, NCPolynomial_plural 

  

cdef object singular_ideal_to_sage_sequence(ideal *i, ring *r, object parent): 

""" 

convert a SINGULAR ideal to a Sage Sequence (the format Sage 

stores a Groebner basis in) 

  

INPUT: 

  

- ``i`` -- a SINGULAR ideal 

- ``r`` -- a SINGULAR ring 

- ``parent`` -- a Sage ring matching r 

""" 

cdef int j 

cdef MPolynomial_libsingular p 

cdef NCPolynomial_plural p_nc 

  

l = [] 

  

if rIsPluralRing(r): 

for j from 0 <= j < IDELEMS(i): 

p_nc = new_NCP(parent, p_Copy(i.m[j], r)) 

l.append( p_nc ) 

else: 

for j from 0 <= j < IDELEMS(i): 

p = new_MP(parent, p_Copy(i.m[j], r)) 

l.append( p ) 

  

return Sequence(l, check=False, immutable=True) 

  

cdef ideal *sage_ideal_to_singular_ideal(I) except NULL: 

""" 

convert a Sage ideal to a SINGULAR ideal 

  

INPUT: 

  

- ``I`` -- a Sage ideal in a ring of type 

:class:`~sage.rings.polynomial.multi_polynomial_libsingular.MPolynomialRing_libsingular` or a list of generators. 

  

TESTS: 

  

  

We test conversion:: 

  

sage: P.<x,y,z> = QQ[] 

sage: sage.libs.singular.function_factory.ff.std(Sequence([x,y,z])) 

[z, y, x] 

sage: sage.libs.singular.function_factory.ff.std(Ideal([x,y,z])) 

[z, y, x] 

""" 

R = I.ring() 

try: 

gens = I.gens() 

except AttributeError: 

gens = I 

cdef ideal *result 

cdef ring *r 

cdef ideal *i 

cdef int j = 0 

  

if isinstance(R, MPolynomialRing_libsingular): 

r = (<MPolynomialRing_libsingular>R)._ring 

elif isinstance(R, NCPolynomialRing_plural): 

r = (<NCPolynomialRing_plural>R)._ring 

else: 

raise TypeError("Ring must be of type 'MPolynomialRing_libsingular'") 

  

rChangeCurrRing(r); 

  

i = idInit(len(gens),1) 

for j,f in enumerate(gens): 

if isinstance(f, MPolynomial_libsingular): 

i.m[j] = p_Copy((<MPolynomial_libsingular>f)._poly, r) 

elif isinstance(f, NCPolynomial_plural): 

i.m[j] = p_Copy((<NCPolynomial_plural>f)._poly, r) 

else: 

id_Delete(&i, r) 

raise TypeError("All generators must be of type MPolynomial_libsingular.") 

return i 

  

def kbase_libsingular(I): 

""" 

SINGULAR's ``kbase()`` algorithm. 

  

INPUT: 

  

- ``I`` -- a groebner basis of an ideal 

  

OUTPUT: 

  

Computes a vector space basis (consisting of monomials) of the quotient 

ring by the ideal, resp. of a free module by the module, in case it is 

finite dimensional and if the input is a standard basis with respect to 

the ring ordering. If the input is not a standard basis, the leading terms 

of the input are used and the result may have no meaning. 

  

EXAMPLES:: 

  

sage: R.<x,y> = PolynomialRing(QQ, order='lex') 

sage: I = R.ideal(x^2-2*y^2, x*y-3) 

sage: I.normal_basis() 

[y^3, y^2, y, 1] 

""" 

  

global singular_options 

  

cdef ideal *i = sage_ideal_to_singular_ideal(I) 

cdef ring *r = currRing 

cdef ideal *q = currRing.qideal 

  

cdef ideal *result 

singular_options = singular_options | Sy_bit(OPT_REDSB) 

  

sig_on() 

result = scKBase(-1, i, q) 

sig_off() 

  

id_Delete(&i, r) 

res = singular_ideal_to_sage_sequence(result,r,I.ring()) 

  

id_Delete(&result, r) 

  

return res 

  

def std_libsingular(I): 

""" 

SINGULAR's ``std()`` algorithm. 

  

INPUT: 

  

- ``I`` -- a Sage ideal 

""" 

global singular_options 

  

cdef ideal *i = sage_ideal_to_singular_ideal(I) 

cdef ring *r = currRing 

cdef tHomog hom = testHomog 

cdef ideal *result 

  

singular_options = singular_options | Sy_bit(OPT_REDSB) 

  

sig_on() 

result =kStd(i,NULL,hom,NULL) 

sig_off() 

  

idSkipZeroes(result) 

  

  

id_Delete(&i,r) 

  

res = singular_ideal_to_sage_sequence(result,r,I.ring()) 

  

id_Delete(&result,r) 

return res 

  

  

def slimgb_libsingular(I): 

""" 

SINGULAR's ``slimgb()`` algorithm. 

  

INPUT: 

  

- ``I`` -- a Sage ideal 

""" 

global singular_options 

  

cdef tHomog hom=testHomog 

cdef ideal *i 

cdef ring *r 

cdef ideal *result 

  

i = sage_ideal_to_singular_ideal(I) 

r = currRing 

  

if r.OrdSgn!=1 : 

id_Delete(&i, r) 

raise TypeError("ordering must be global for slimgb") 

  

if i.rank < id_RankFreeModule(i, r): 

id_Delete(&i, r) 

raise TypeError 

  

singular_options = singular_options | Sy_bit(OPT_REDSB) 

  

sig_on() 

result = t_rep_gb(r, i, i.rank, 0) 

sig_off() 

  

id_Delete(&i,r) 

  

res = singular_ideal_to_sage_sequence(result,r,I.ring()) 

  

id_Delete(&result,r) 

return res 

  

def interred_libsingular(I): 

""" 

SINGULAR's ``interred()`` command. 

  

INPUT: 

  

- ``I`` -- a Sage ideal 

  

EXAMPLES:: 

  

sage: P.<x,y,z> = PolynomialRing(ZZ) 

sage: I = ideal( x^2 - 3*y, y^3 - x*y, z^3 - x, x^4 - y*z + 1 ) 

sage: I.interreduced_basis() 

[y*z^2 - 81*x*y - 9*y - z, z^3 - x, x^2 - 3*y, 9*y^2 - y*z + 1] 

  

sage: P.<x,y,z> = PolynomialRing(QQ) 

sage: I = ideal( x^2 - 3*y, y^3 - x*y, z^3 - x, x^4 - y*z + 1 ) 

sage: I.interreduced_basis() 

[y*z^2 - 81*x*y - 9*y - z, z^3 - x, x^2 - 3*y, y^2 - 1/9*y*z + 1/9] 

""" 

global singular_options 

  

cdef ideal *i 

cdef ideal *result 

cdef ring *r 

cdef number *n 

cdef poly *p 

cdef int j 

cdef int bck 

  

try: 

if len(I.gens()) == 0: 

return Sequence([], check=False, immutable=True) 

except AttributeError: 

pass 

  

i = sage_ideal_to_singular_ideal(I) 

r = currRing 

  

bck = singular_options 

singular_options = singular_options | Sy_bit(OPT_REDTAIL)|Sy_bit(OPT_REDSB) 

sig_on() 

result = kInterRed(i,NULL) 

sig_off() 

singular_options = bck 

  

  

# divide head by coefficients 

if r.cf.type != n_Z and r.cf.type != n_Znm and r.cf.type != n_Zn and r.cf.type != n_Z2m : 

for j from 0 <= j < IDELEMS(result): 

p = result.m[j] 

if p: 

n = p_GetCoeff(p,r) 

n = r.cf.cfInvers(n,r.cf) 

result.m[j] = pp_Mult_nn(p, n, r) 

p_Delete(&p,r) 

n_Delete(&n,r) 

  

id_Delete(&i,r) 

  

res = sorted(singular_ideal_to_sage_sequence(result,r,I.ring()),reverse=True) 

  

id_Delete(&result,r) 

return res