|
| 1 | +# -*- coding: utf-8 -*- |
| 2 | +""" |
| 3 | +Wrappers for Giac functions |
| 4 | +
|
| 5 | +We provide a python function to compute and convert to sage a groebner |
| 6 | +basis using the giacpy module. |
| 7 | +
|
| 8 | +AUTHORS: |
| 9 | +
|
| 10 | +- Martin Albrecht (2015-07-01): initial version |
| 11 | +- Han Frederic (2015-07-01): initial version |
| 12 | +
|
| 13 | +EXAMPLES:: |
| 14 | +
|
| 15 | + sage: from sage.libs.giac import groebner_basis as gb_giac # optional - giacpy |
| 16 | + sage: P = PolynomialRing(QQ, 6, 'x') |
| 17 | + sage: I = sage.rings.ideal.Cyclic(P) |
| 18 | + sage: B = gb_giac(I.gens()) # optional - giacpy, random |
| 19 | + sage: B # optional - giacpy |
| 20 | + Polynomial Sequence with 45 Polynomials in 6 Variables |
| 21 | +""" |
| 22 | + |
| 23 | +# ***************************************************************************** |
| 24 | +# Copyright (C) 2013 Frederic Han <[email protected]> |
| 25 | +# |
| 26 | +# This program is free software: you can redistribute it and/or modify |
| 27 | +# it under the terms of the GNU General Public License as published by |
| 28 | +# the Free Software Foundation, either version 2 of the License, or |
| 29 | +# (at your option) any later version. |
| 30 | +# http://www.gnu.org/licenses/ |
| 31 | +# ***************************************************************************** |
| 32 | + |
| 33 | +from sage.structure.proof.all import polynomial as proof_polynomial |
| 34 | +from sage.rings.polynomial.multi_polynomial_sequence import PolynomialSequence |
| 35 | + |
| 36 | +# Remarks for doctests: |
| 37 | +# 1) The first time that the c++ library giac is loaded a message appears. |
| 38 | +# This message is version and arch dependant. |
| 39 | +# 2) When proba_epsilon is too bad (>1e-6?) setting it to a better value |
| 40 | +# will give an additional message like the following one: |
| 41 | +# Restoring proba epsilon to 1e-6 from 1e-12 |
| 42 | +# (it looks like in internal giac changes this also to not work with a too bad probability) |
| 43 | + |
| 44 | + |
| 45 | +class GiacSettingsDefaultContext: |
| 46 | + """ |
| 47 | + Context preserve libgiac settings. |
| 48 | + """ |
| 49 | + |
| 50 | + def __enter__(self): |
| 51 | + """ |
| 52 | + EXAMPLE:: |
| 53 | +
|
| 54 | + sage: from sage.libs.giac import GiacSettingsDefaultContext # optional - giacpy |
| 55 | + sage: from giacpy import giacsettings # optional - giacpy |
| 56 | + sage: giacsettings.proba_epsilon = 1e-16 # optional - giacpy |
| 57 | + sage: with GiacSettingsDefaultContext(): giacsettings.proba_epsilon = 1e-12 # optional - giacpy |
| 58 | + sage: giacsettings.proba_epsilon < 1e-14 # optional - giacpy |
| 59 | + True |
| 60 | +
|
| 61 | + """ |
| 62 | + try: |
| 63 | + from giacpy import giacsettings, libgiac |
| 64 | + except ImportError: |
| 65 | + raise ImportError("""One of the optional packages giac or giacpy is missing""") |
| 66 | + |
| 67 | + self.proba_epsilon = giacsettings.proba_epsilon |
| 68 | + self.threads = giacsettings.threads |
| 69 | + # Change the debug level at the end to not have messages at each modification |
| 70 | + self.debuginfolevel = libgiac('debug_infolevel()') |
| 71 | + |
| 72 | + def __exit__(self, typ, value, tb): |
| 73 | + """ |
| 74 | + EXAMPLE:: |
| 75 | +
|
| 76 | + sage: from sage.libs.giac import GiacSettingsDefaultContext # optional - giacpy |
| 77 | + sage: from giacpy import giacsettings # optional - giacpy |
| 78 | + sage: giacsettings.proba_epsilon = 1e-16 # optional - giacpy |
| 79 | + sage: with GiacSettingsDefaultContext(): giacsettings.proba_epsilon = 1e-30 # optional - giacpy |
| 80 | + sage: giacsettings.proba_epsilon < 1e-20 # optional - giacpy |
| 81 | + False |
| 82 | +
|
| 83 | + """ |
| 84 | + try: |
| 85 | + from giacpy import giacsettings, libgiac |
| 86 | + except ImportError: |
| 87 | + raise ImportError("""One of the optional packages giac or giacpy is missing""") |
| 88 | + |
| 89 | + # Restore the debug level first to not have messages at each modification |
| 90 | + libgiac('debug_infolevel')(self.debuginfolevel) |
| 91 | + # NB: giacsettings.epsilon has a different meaning that giacsettings.proba_epsilon. |
| 92 | + giacsettings.proba_epsilon = self.proba_epsilon |
| 93 | + giacsettings.threads = self.threads |
| 94 | + |
| 95 | +def local_giacsettings(func): |
| 96 | + """ |
| 97 | + Decorator to preserve Giac's proba_epsilon and threads settings. |
| 98 | +
|
| 99 | + EXAMPLE:: |
| 100 | +
|
| 101 | + sage: def testf(a,b): # optional - giacpy |
| 102 | + ....: giacsettings.proba_epsilon = a/100 |
| 103 | + ....: giacsettings.threads = b+2 |
| 104 | + ....: return (giacsettings.proba_epsilon, giacsettings.threads) |
| 105 | +
|
| 106 | + sage: from giacpy import giacsettings # optional - giacpy |
| 107 | + sage: from sage.libs.giac import local_giacsettings # optional - giacpy |
| 108 | + sage: gporig, gtorig = (giacsettings.proba_epsilon,giacsettings.threads) # optional - giacpy |
| 109 | + sage: gp, gt = local_giacsettings(testf)(giacsettings.proba_epsilon,giacsettings.threads) # optional - giacpy |
| 110 | + sage: gporig == giacsettings.proba_epsilon # optional - giacpy |
| 111 | + True |
| 112 | + sage: gtorig == giacsettings.threads # optional - giacpy |
| 113 | + True |
| 114 | + sage: gp<gporig, gt-gtorig # optional - giacpy |
| 115 | + (True, 2) |
| 116 | +
|
| 117 | + """ |
| 118 | + from sage.misc.decorators import sage_wraps |
| 119 | + |
| 120 | + @sage_wraps(func) |
| 121 | + def wrapper(*args, **kwds): |
| 122 | + """ |
| 123 | + Execute function in ``GiacSettingsDefaultContext``. |
| 124 | + """ |
| 125 | + with GiacSettingsDefaultContext(): |
| 126 | + return func(*args, **kwds) |
| 127 | + |
| 128 | + return wrapper |
| 129 | + |
| 130 | + |
| 131 | +@local_giacsettings |
| 132 | +def groebner_basis(gens, proba_epsilon=None, threads=None, prot=False, *args, **kwds): |
| 133 | + """ |
| 134 | + Computes a Groebner Basis of an ideal using giacpy. The result is |
| 135 | + automatically converted to sage. |
| 136 | +
|
| 137 | + INPUT: |
| 138 | +
|
| 139 | + - ``gens`` - an ideal (or a list) of polynomials over a prime field |
| 140 | + of characteristic 0 or p<2^31 |
| 141 | +
|
| 142 | + - ``proba_epsilon`` - (default: None) majoration of the probability |
| 143 | + of a wrong answer when probabilistic algorithms are allowed. |
| 144 | +
|
| 145 | + * if ``proba_epsilon`` is None, the value of |
| 146 | + ``sage.structure.proof.all.polynomial()`` is taken. If it is |
| 147 | + false then the global ``giacpy.giacsettings.proba_epsilon`` is |
| 148 | + used. |
| 149 | +
|
| 150 | + * if ``proba_epsilon`` is 0, probabilistic algorithms are |
| 151 | + disabled. |
| 152 | +
|
| 153 | + - ``threads`` - (default: None) Maximal number of threads allowed |
| 154 | + for giac. If None, the global ``giacpy.giacsettings.threads`` is |
| 155 | + considered. |
| 156 | +
|
| 157 | + - ``prot`` - (default: False) if True print detailled informations |
| 158 | +
|
| 159 | + OUTPUT: |
| 160 | +
|
| 161 | + Polynomial sequence of the reduced Groebner basis. |
| 162 | +
|
| 163 | + EXAMPLES:: |
| 164 | +
|
| 165 | + sage: from sage.libs.giac import groebner_basis as gb_giac # optional - giacpy |
| 166 | + sage: P = PolynomialRing(GF(previous_prime(2**31)), 6, 'x') # optional - giacpy |
| 167 | + sage: I = sage.rings.ideal.Cyclic(P) # optional - giacpy |
| 168 | + sage: B=gb_giac(I.gens());B # optional - giacpy |
| 169 | + Polynomial Sequence with 45 Polynomials in 6 Variables |
| 170 | + sage: B.is_groebner() # optional - giacpy |
| 171 | + True |
| 172 | +
|
| 173 | + Computations over QQ can benefit from |
| 174 | +
|
| 175 | + * a probabilistic lifting:: |
| 176 | +
|
| 177 | + sage: P = PolynomialRing(QQ,5, 'x') # optional - giacpy |
| 178 | + sage: I = ideal([P.random_element(3,7) for j in range(5)]) # optional - giacpy |
| 179 | + sage: B1 = gb_giac(I.gens(),1e-16) # optional - giacpy, long time (1s) |
| 180 | + Running a probabilistic check for the reconstructed Groebner basis. |
| 181 | + If successfull, error probability is less than 1e-16 ... |
| 182 | + sage: sage.structure.proof.all.polynomial(True) # optional - giacpy |
| 183 | + sage: B2 = gb_giac(I.gens()) # optional - giacpy, long time (4s) |
| 184 | + sage: B1==B2 # optional - giacpy, long time |
| 185 | + True |
| 186 | + sage: B1.is_groebner() # optional - giacpy, long time (20s) |
| 187 | + True |
| 188 | +
|
| 189 | + * multi threaded operations:: |
| 190 | +
|
| 191 | + sage: P = PolynomialRing(QQ, 8, 'x') # optional - giacpy |
| 192 | + sage: I=sage.rings.ideal.Cyclic(P) # optional - giacpy |
| 193 | + sage: time B = gb_giac(I.gens(),1e-6,threads=2) # doctest: +SKIP |
| 194 | + Running a probabilistic check for the reconstructed Groebner basis... |
| 195 | + Time: CPU 168.98 s, Wall: 94.13 s |
| 196 | +
|
| 197 | + You can get detailled information by setting ``prot=True`` |
| 198 | +
|
| 199 | + :: |
| 200 | +
|
| 201 | + sage: I=sage.rings.ideal.Katsura(P) # optional - giacpy |
| 202 | + sage: gb_giac(I,prot=True) # optional - giacpy, random, long time (3s) |
| 203 | + 9381383 begin computing basis modulo 535718473 |
| 204 | + 9381501 begin new iteration zmod, number of pairs: 8, base size: 8 |
| 205 | + ...end, basis size 74 prime number 1 |
| 206 | + G=Vector [0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18,... |
| 207 | + ...creating reconstruction #0 |
| 208 | + ... |
| 209 | + ++++++++basis size 74 |
| 210 | + checking pairs for i=0, j= |
| 211 | + checking pairs for i=1, j=2,6,12,17,19,24,29,34,39,42,43,48,56,61,64,69, |
| 212 | + ... |
| 213 | + checking pairs for i=72, j=73, |
| 214 | + checking pairs for i=73, j= |
| 215 | + Number of critical pairs to check 373 |
| 216 | + +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++... |
| 217 | + Successfull check of 373 critical pairs |
| 218 | + 12380865 end final check |
| 219 | + Polynomial Sequence with 74 Polynomials in 8 Variables |
| 220 | +
|
| 221 | +
|
| 222 | + TESTS:: |
| 223 | +
|
| 224 | + sage: from giacpy import libgiac # optional - giacpy |
| 225 | + sage: libgiac("x2:=22; x4:='whywouldyoudothis'") # optional - giacpy |
| 226 | + 22,whywouldyoudothis |
| 227 | + sage: gb_giac(I) # optional - giacpy |
| 228 | + Traceback (most recent call last): |
| 229 | + ... |
| 230 | + ValueError: Variables names ['x2', 'x4'] conflict in giac. Change them or purge them from in giac with libgiac.purge('x2') |
| 231 | + sage: libgiac.purge('x2'),libgiac.purge('x4') # optional - giacpy |
| 232 | + (22, whywouldyoudothis) |
| 233 | + sage: gb_giac(I) # optional - giacpy, long time (3s) |
| 234 | + Polynomial Sequence with 74 Polynomials in 8 Variables |
| 235 | +
|
| 236 | + """ |
| 237 | + try: |
| 238 | + from giacpy import libgiac, giacsettings |
| 239 | + except ImportError: |
| 240 | + raise ImportError("""One of the optional packages giac or giacpy is missing""") |
| 241 | + |
| 242 | + try: |
| 243 | + iter(gens) |
| 244 | + except TypeError: |
| 245 | + gens = gens.gens() |
| 246 | + |
| 247 | + # get the ring from gens |
| 248 | + P = iter(gens).next().parent() |
| 249 | + K = P.base_ring() |
| 250 | + p = K.characteristic() |
| 251 | + |
| 252 | + # check for name confusions |
| 253 | + blackgiacconstants = ['i', 'e'] # NB e^k is expanded to exp(k) |
| 254 | + blacklist = blackgiacconstants + [str(j) for j in libgiac.VARS()] |
| 255 | + problematicnames = list(set(P.gens_dict().keys()).intersection(blacklist)) |
| 256 | + |
| 257 | + if(len(problematicnames)>0): |
| 258 | + raise ValueError("Variables names %s conflict in giac. Change them or purge them from in giac with libgiac.purge(\'%s\')" |
| 259 | + %(problematicnames, problematicnames[0])) |
| 260 | + |
| 261 | + if K.is_prime_field() and p == 0: |
| 262 | + F = libgiac(gens) |
| 263 | + elif K.is_prime_field() and p < 2**31: |
| 264 | + F = (libgiac(gens) % p) |
| 265 | + else: |
| 266 | + raise NotImplementedError("Only prime fields of cardinal < 2^31 are implemented in Giac for Groebner bases.") |
| 267 | + |
| 268 | + if P.term_order() != "degrevlex": |
| 269 | + raise NotImplementedError("Only degrevlex term orderings are supported in Giac Groebner bases.") |
| 270 | + |
| 271 | + # proof or probabilistic reconstruction |
| 272 | + if proba_epsilon is None: |
| 273 | + if proof_polynomial(): |
| 274 | + giacsettings.proba_epsilon = 0 |
| 275 | + else: |
| 276 | + giacsettings.proba_epsilon = 1e-15 |
| 277 | + else: |
| 278 | + giacsettings.proba_epsilon = proba_epsilon |
| 279 | + |
| 280 | + # prot |
| 281 | + if prot: |
| 282 | + libgiac('debug_infolevel(2)') |
| 283 | + |
| 284 | + # threads |
| 285 | + if threads is not None: |
| 286 | + giacsettings.threads = threads |
| 287 | + |
| 288 | + # compute de groebner basis with giac |
| 289 | + gb_giac = F.gbasis([P.gens()], "revlex") |
| 290 | + |
| 291 | + return PolynomialSequence(gb_giac, P, immutable=True) |
0 commit comments