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

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

381

382

383

384

385

386

387

388

389

390

391

392

393

394

395

396

397

398

399

400

401

402

403

404

405

406

407

408

409

410

411

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

455

456

457

458

459

460

461

462

463

464

465

466

467

468

469

470

471

472

473

474

475

476

477

478

479

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

521

522

523

524

525

526

527

528

529

530

531

532

533

534

535

536

537

538

539

540

541

542

543

544

545

546

547

548

549

550

551

""" 

Singular's Groebner Strategy Objects 

  

AUTHORS: 

  

- Martin Albrecht (2009-07): initial implementation 

- Michael Brickenstein (2009-07): initial implementation 

- Hans Schoenemann (2009-07): initial implementation 

""" 

  

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

# Copyright (C) 2009 Martin Albrecht <M.R.Albrecht@rhul.ac.uk> 

# Copyright (C) 2009 Michael Brickenstein <brickenstein@mfo.de> 

# Copyright (C) 2009 Hans Schoenemann <hannes@mathematik.uni-kl.de> 

# 

# Distributed under the terms of the GNU General Public License (GPL) 

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

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

  

cdef extern from *: # hack to get at cython macro 

int unlikely(int) 

int likely(int) 

  

from sage.structure.richcmp cimport richcmp 

from sage.libs.singular.decl cimport ideal, ring, poly, currRing 

from sage.libs.singular.decl cimport rChangeCurrRing 

from sage.libs.singular.decl cimport new_skStrategy, delete_skStrategy, id_RankFreeModule 

from sage.libs.singular.decl cimport initEcartBBA, enterSBba, initBuchMoraCrit, initS, pNorm, id_Delete, kTest 

from sage.libs.singular.decl cimport omfree, redNF, p_Copy, redtailBba 

  

from sage.libs.singular.ring cimport singular_ring_reference, singular_ring_delete 

  

from sage.rings.polynomial.multi_polynomial_ideal import MPolynomialIdeal, NCPolynomialIdeal 

from sage.rings.polynomial.multi_polynomial_ideal_libsingular cimport sage_ideal_to_singular_ideal 

from sage.rings.polynomial.multi_polynomial_libsingular cimport MPolynomial_libsingular, MPolynomialRing_libsingular, new_MP 

from sage.rings.polynomial.plural cimport new_NCP 

  

cdef class GroebnerStrategy(SageObject): 

""" 

A Wrapper for Singular's Groebner Strategy Object. 

  

This object provides functions for normal form computations and 

other functions for Groebner basis computation. 

  

ALGORITHM: 

  

Uses Singular via libSINGULAR 

""" 

def __cinit__(self, L): 

""" 

Create a new :class:`GroebnerStrategy` object for the 

generators of the ideal ``L``. 

  

INPUT: 

  

- ``L`` - a multivariate polynomial ideal 

  

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x+z,y+z+1]) 

sage: strat = GroebnerStrategy(I); strat 

Groebner Strategy for ideal generated by 2 elements 

over Multivariate Polynomial Ring in x, y, z over Rational Field 

  

TESTS:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

sage: strat = GroebnerStrategy(None) 

Traceback (most recent call last): 

... 

TypeError: First parameter must be a multivariate polynomial ideal. 

  

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

sage: I = Ideal([x+z,y+z+1]) 

sage: strat = GroebnerStrategy(I) 

Traceback (most recent call last): 

... 

NotImplementedError: The local case is not implemented yet. 

  

sage: P.<x,y,z> = PolynomialRing(CC,order='neglex') 

sage: I = Ideal([x+z,y+z+1]) 

sage: strat = GroebnerStrategy(I) 

Traceback (most recent call last): 

... 

TypeError: First parameter's ring must be multivariate polynomial ring via libsingular. 

  

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

sage: I = Ideal([x+z,y+z+1]) 

sage: strat = GroebnerStrategy(I) 

Traceback (most recent call last): 

... 

NotImplementedError: Only coefficient fields are implemented so far. 

  

""" 

if not isinstance(L, MPolynomialIdeal): 

raise TypeError("First parameter must be a multivariate polynomial ideal.") 

  

if not isinstance(L.ring(), MPolynomialRing_libsingular): 

raise TypeError("First parameter's ring must be multivariate polynomial ring via libsingular.") 

  

self._ideal = L 

  

cdef MPolynomialRing_libsingular R = <MPolynomialRing_libsingular>L.ring() 

self._parent = R 

self._parent_ring = singular_ring_reference(R._ring) 

  

if not R.term_order().is_global(): 

raise NotImplementedError("The local case is not implemented yet.") 

  

if not R.base_ring().is_field(): 

raise NotImplementedError("Only coefficient fields are implemented so far.") 

  

if (R._ring != currRing): 

rChangeCurrRing(R._ring) 

  

cdef ideal *i = sage_ideal_to_singular_ideal(L) 

self._strat = new_skStrategy() 

  

self._strat.ak = id_RankFreeModule(i, R._ring) 

#- creating temp data structures 

initBuchMoraCrit(self._strat) 

self._strat.initEcart = initEcartBBA 

self._strat.enterS = enterSBba 

#- set S 

self._strat.sl = -1 

#- init local data struct 

initS(i, NULL, self._strat) 

  

cdef int j 

cdef bint base_ring_is_field = R.base_ring().is_field() 

if (R._ring != currRing): 

rChangeCurrRing(R._ring) 

if base_ring_is_field: 

for j in range(self._strat.sl, -1, -1): 

pNorm(self._strat.S[j]) 

  

id_Delete(&i, R._ring) 

kTest(self._strat) 

  

def __dealloc__(self): 

""" 

TESTS:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x + z, y + z]) 

sage: strat = GroebnerStrategy(I) 

sage: del strat 

""" 

# WARNING: the Cython class self._parent is no longer accessible! 

# see http://trac.sagemath.org/sage_trac/ticket/11339 

cdef ring *oldRing = NULL 

if self._strat: 

omfree(self._strat.sevS) 

omfree(self._strat.ecartS) 

omfree(self._strat.T) 

omfree(self._strat.sevT) 

omfree(self._strat.R) 

omfree(self._strat.S_2_R) 

omfree(self._strat.L) 

omfree(self._strat.B) 

omfree(self._strat.fromQ) 

id_Delete(&self._strat.Shdl, self._parent_ring) 

  

if self._parent_ring != currRing: 

oldRing = currRing 

rChangeCurrRing(self._parent_ring) 

delete_skStrategy(self._strat) 

rChangeCurrRing(oldRing) 

else: 

delete_skStrategy(self._strat) 

if self._parent_ring: 

singular_ring_delete(self._parent_ring) 

  

def _repr_(self): 

""" 

TESTS:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x + z, y + z]) 

sage: strat = GroebnerStrategy(I) 

sage: strat # indirect doctest 

Groebner Strategy for ideal generated by 2 elements over 

Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003 

""" 

return "Groebner Strategy for ideal generated by %d elements over %s"%(self._ideal.ngens(),self._parent) 

  

def ideal(self): 

""" 

Return the ideal this strategy object is defined for. 

  

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x + z, y + z]) 

sage: strat = GroebnerStrategy(I) 

sage: strat.ideal() 

Ideal (x + z, y + z) of Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003 

""" 

return self._ideal 

  

def ring(self): 

""" 

Return the ring this strategy object is defined over. 

  

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x + z, y + z]) 

sage: strat = GroebnerStrategy(I) 

sage: strat.ring() 

Multivariate Polynomial Ring in x, y, z over Finite Field of size 32003 

""" 

return self._parent 

  

def __richcmp__(self, other, op): 

""" 

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([P(0)]) 

sage: strat = GroebnerStrategy(I) 

sage: strat == GroebnerStrategy(I) 

True 

sage: I = Ideal([x+1,y+z]) 

sage: strat == GroebnerStrategy(I) 

False 

""" 

try: 

lx = <GroebnerStrategy?>self 

rx = <GroebnerStrategy?>other 

except TypeError: 

return NotImplemented 

return richcmp(lx._ideal.gens(), 

rx._ideal.gens(), op) 

  

def __reduce__(self): 

""" 

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x + z, y + z]) 

sage: strat = GroebnerStrategy(I) 

sage: loads(dumps(strat)) == strat 

True 

""" 

return unpickle_GroebnerStrategy0, (self._ideal,) 

  

cpdef MPolynomial_libsingular normal_form(self, MPolynomial_libsingular p): 

""" 

Compute the normal form of ``p`` with respect to the 

generators of this object. 

  

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x + z, y + z]) 

sage: strat = GroebnerStrategy(I) 

sage: strat.normal_form(x*y) # indirect doctest 

z^2 

sage: strat.normal_form(x + 1) 

-z + 1 

  

TESTS:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([P(0)]) 

sage: strat = GroebnerStrategy(I) 

sage: strat.normal_form(x) 

x 

sage: strat.normal_form(P(0)) 

0 

""" 

if unlikely(p._parent is not self._parent): 

raise TypeError("parent(p) must be the same as this object's parent.") 

if unlikely(self._parent._ring != currRing): 

rChangeCurrRing(self._parent._ring) 

  

cdef int max_ind = 0 

cdef poly *_p = redNF(p_Copy(p._poly, self._parent._ring), max_ind, 0, self._strat) 

if likely(_p!=NULL): 

_p = redtailBba(_p, max_ind, self._strat) 

return new_MP(self._parent, _p) 

  

cdef class NCGroebnerStrategy(SageObject): 

""" 

A Wrapper for Singular's Groebner Strategy Object. 

  

This object provides functions for normal form computations and 

other functions for Groebner basis computation. 

  

ALGORITHM: 

  

Uses Singular via libSINGULAR 

""" 

def __init__(self, L): 

""" 

Create a new :class:`GroebnerStrategy` object for the 

generators of the ideal ``L``. 

  

INPUT: 

  

- ``L`` - an ideal in a g-algebra 

  

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: NCGroebnerStrategy(I) #random 

Groebner Strategy for ideal generated by 3 elements over Noncommutative Multivariate Polynomial Ring in x, y, z over Rational Field, nc-relations: {z*x: x*z + 2*x, z*y: y*z - 2*y, y*x: x*y - z} 

  

TESTS:: 

  

sage: strat = NCGroebnerStrategy(None) 

Traceback (most recent call last): 

... 

TypeError: First parameter must be an ideal in a g-algebra. 

  

sage: P.<x,y,z> = PolynomialRing(CC,order='neglex') 

sage: I = Ideal([x+z,y+z+1]) 

sage: strat = NCGroebnerStrategy(I) 

Traceback (most recent call last): 

... 

TypeError: First parameter must be an ideal in a g-algebra. 

  

""" 

if not isinstance(L, NCPolynomialIdeal): 

raise TypeError("First parameter must be an ideal in a g-algebra.") 

  

if not isinstance(L.ring(), NCPolynomialRing_plural): 

raise TypeError("First parameter's ring must be a g-algebra.") 

  

self._ideal = L 

  

cdef NCPolynomialRing_plural R = <NCPolynomialRing_plural>L.ring() 

self._parent = R 

  

if not R.term_order().is_global(): 

raise NotImplementedError("The local case is not implemented yet.") 

  

if not R.base_ring().is_field(): 

raise NotImplementedError("Only coefficient fields are implemented so far.") 

  

if (R._ring != currRing): 

rChangeCurrRing(R._ring) 

  

cdef ideal *i = sage_ideal_to_singular_ideal(L) 

self._strat = new_skStrategy() 

  

self._strat.ak = id_RankFreeModule(i, R._ring) 

#- creating temp data structures 

initBuchMoraCrit(self._strat) 

self._strat.initEcart = initEcartBBA 

self._strat.enterS = enterSBba 

#- set S 

self._strat.sl = -1 

#- init local data struct 

initS(i, NULL, self._strat) 

  

cdef int j 

if R.base_ring().is_field(): 

for j in range(self._strat.sl,-1,-1): 

pNorm(self._strat.S[j]) 

  

id_Delete(&i, R._ring) 

kTest(self._strat) 

  

def __dealloc__(self): 

""" 

TESTS:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: strat = NCGroebnerStrategy(I) 

sage: del strat # indirect doctest 

""" 

cdef ring *oldRing = NULL 

if self._strat: 

omfree(self._strat.sevS) 

omfree(self._strat.ecartS) 

omfree(self._strat.T) 

omfree(self._strat.sevT) 

omfree(self._strat.R) 

omfree(self._strat.S_2_R) 

omfree(self._strat.L) 

omfree(self._strat.B) 

omfree(self._strat.fromQ) 

id_Delete(&self._strat.Shdl, self._parent._ring) 

  

if self._parent._ring != currRing: 

oldRing = currRing 

rChangeCurrRing(self._parent._ring) 

delete_skStrategy(self._strat) 

rChangeCurrRing(oldRing) 

else: 

delete_skStrategy(self._strat) 

  

def _repr_(self): 

""" 

TESTS:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: strat = NCGroebnerStrategy(I) 

sage: strat # indirect doctest #random 

Groebner Strategy for ideal generated by 3 elements over Noncommutative Multivariate Polynomial Ring in x, y, z over Rational Field, nc-relations: {z*x: x*z + 2*x, z*y: y*z - 2*y, y*x: x*y - z} 

""" 

return "Groebner Strategy for ideal generated by %d elements over %s"%(self._ideal.ngens(),self._parent) 

  

def ideal(self): 

""" 

Return the ideal this strategy object is defined for. 

  

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: strat = NCGroebnerStrategy(I) 

sage: strat.ideal() == I 

True 

  

""" 

return self._ideal 

  

def ring(self): 

""" 

Return the ring this strategy object is defined over. 

  

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: strat = NCGroebnerStrategy(I) 

sage: strat.ring() is H 

True 

""" 

return self._parent 

  

def __richcmp__(self, other, op): 

""" 

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: strat = NCGroebnerStrategy(I) 

sage: strat == NCGroebnerStrategy(I) 

True 

sage: I = H.ideal([y^2, x^2, z^2-H.one()], side='twosided') 

sage: strat == NCGroebnerStrategy(I) 

False 

""" 

try: 

lx = <NCGroebnerStrategy?>self 

rx = <NCGroebnerStrategy?>other 

except TypeError: 

return NotImplemented 

return richcmp((lx._ideal.gens(), lx._ideal.side()), 

(rx._ideal.gens(), rx._ideal.side()), op) 

  

def __reduce__(self): 

""" 

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: strat = NCGroebnerStrategy(I) 

sage: loads(dumps(strat)) == strat 

True 

""" 

return unpickle_NCGroebnerStrategy0, (self._ideal,) 

  

cpdef NCPolynomial_plural normal_form(self, NCPolynomial_plural p): 

""" 

Compute the normal form of ``p`` with respect to the 

generators of this object. 

  

EXAMPLES:: 

  

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: JL = H.ideal([x^3, y^3, z^3 - 4*z]) 

sage: JT = H.ideal([x^3, y^3, z^3 - 4*z], side='twosided') 

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: SL = NCGroebnerStrategy(JL.std()) 

sage: ST = NCGroebnerStrategy(JT.std()) 

sage: SL.normal_form(x*y^2) 

x*y^2 

sage: ST.normal_form(x*y^2) 

y*z 

  

""" 

if unlikely(p._parent is not self._parent): 

raise TypeError("parent(p) must be the same as this object's parent.") 

if unlikely(self._parent._ring != currRing): 

rChangeCurrRing(self._parent._ring) 

  

cdef int max_ind 

cdef poly *_p = redNF(p_Copy(p._poly, self._parent._ring), max_ind, 0, self._strat) 

if likely(_p!=NULL): 

_p = redtailBba(_p, max_ind, self._strat) 

return new_NCP(self._parent, _p) 

  

def unpickle_NCGroebnerStrategy0(I): 

""" 

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import NCGroebnerStrategy 

sage: A.<x,y,z> = FreeAlgebra(QQ, 3) 

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}) 

sage: I = H.ideal([y^2, x^2, z^2-H.one()]) 

sage: strat = NCGroebnerStrategy(I) 

sage: loads(dumps(strat)) == strat # indirect doctest 

True 

""" 

return NCGroebnerStrategy(I) 

  

def unpickle_GroebnerStrategy0(I): 

""" 

EXAMPLES:: 

  

sage: from sage.libs.singular.groebner_strategy import GroebnerStrategy 

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

sage: I = Ideal([x + z, y + z]) 

sage: strat = GroebnerStrategy(I) 

sage: loads(dumps(strat)) == strat # indirect doctest 

True 

""" 

return GroebnerStrategy(I)