libSingular: Functions
Sage implements a C wrapper around the Singular interpreter which allows to call any function directly from Sage without string parsing or interprocess communication overhead. Users who do not want to call Singular functions directly, usually do not have to worry about this interface, since it is handled by higher level functions in Sage.
AUTHORS:
EXAMPLES:
The direct approach for loading a Singular function is to call the function singular_function() with the function name as parameter:
sage: from sage.libs.singular.function import singular_function
sage: P.<a,b,c,d> = PolynomialRing(GF(7))
sage: std = singular_function('std')
sage: I = sage.rings.ideal.Cyclic(P)
sage: std(I)
[a + b + c + d,
b^2 + 2*b*d + d^2,
b*c^2 + c^2*d - b*d^2 - d^3,
b*c*d^2 + c^2*d^2 - b*d^3 + c*d^3 - d^4 - 1,
b*d^4 + d^5 - b - d,
c^3*d^2 + c^2*d^3 - c - d,
c^2*d^4 + b*c - b*d + c*d - 2*d^2]
If a Singular library needs to be loaded before a certain function is available, use the lib() function as shown below:
sage: from sage.libs.singular.function import singular_function, lib as singular_lib
sage: primdecSY = singular_function('primdecSY')
Traceback (most recent call last):
...
NameError: Function 'primdecSY' is not defined.
sage: singular_lib('primdec.lib')
sage: primdecSY = singular_function('primdecSY')
There is also a short-hand notation for the above:
sage: primdecSY = sage.libs.singular.ff.primdec__lib.primdecSY
The above line will load “primdec.lib” first and then load the function primdecSY.
TESTS:
sage: from sage.libs.singular.function import singular_function
sage: std = singular_function('std')
sage: loads(dumps(std)) == std
True
Bases: object
A call handler is an abstraction which hides the details of the implementation differences between kernel and library functions.
Bases: sage.structure.sage_object.SageObject
A Converter interfaces between Sage objects and Singular interpreter objects.
Return the ring in which the arguments of this list live.
EXAMPLE:
sage: from sage.libs.singular.function import Converter
sage: P.<a,b,c> = PolynomialRing(GF(127))
sage: Converter([a,b,c],ring=P).ring()
Multivariate Polynomial Ring in a, b, c over Finite Field of size 127
Bases: sage.libs.singular.function.BaseCallHandler
A call handler is an abstraction which hides the details of the implementation differences between kernel and library functions.
This class implements calling a kernel function.
Note
Do not construct this class directly, use singular_function() instead.
Bases: sage.libs.singular.function.BaseCallHandler
A call handler is an abstraction which hides the details of the implementation differences between kernel and library functions.
This class implements calling a library function.
Note
Do not construct this class directly, use singular_function() instead.
Bases: object
A simple wrapper around Singular’s resolutions.
Bases: object
A simple wrapper around Singular’s rings.
Get characteristic.
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: ringlist = singular_function("ringlist")
sage: l = ringlist(P)
sage: ring = singular_function("ring")
sage: ring(l, ring=P).characteristic()
0
Determine whether a given ring is commutative.
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: ringlist = singular_function("ringlist")
sage: l = ringlist(P)
sage: ring = singular_function("ring")
sage: ring(l, ring=P).is_commutative()
True
Get number of generators.
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: ringlist = singular_function("ringlist")
sage: l = ringlist(P)
sage: ring = singular_function("ring")
sage: ring(l, ring=P).ngens()
3
Get number of parameters.
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: ringlist = singular_function("ringlist")
sage: l = ringlist(P)
sage: ring = singular_function("ring")
sage: ring(l, ring=P).npars()
0
Get Singular string defining monomial ordering.
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: ringlist = singular_function("ringlist")
sage: l = ringlist(P)
sage: ring = singular_function("ring")
sage: ring(l, ring=P).ordering_string()
'dp(3),C'
Get parameter names.
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: ringlist = singular_function("ringlist")
sage: l = ringlist(P)
sage: ring = singular_function("ring")
sage: ring(l, ring=P).par_names()
[]
Get names of variables.
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: ringlist = singular_function("ringlist")
sage: l = ringlist(P)
sage: ring = singular_function("ring")
sage: ring(l, ring=P).var_names()
['x', 'y', 'z']
Bases: sage.structure.sage_object.SageObject
The base class for Singular functions either from the kernel or from the library.
Bases: sage.libs.singular.function.SingularFunction
EXAMPLES:
sage: from sage.libs.singular.function import SingularKernelFunction
sage: R.<x,y> = PolynomialRing(QQ, order='lex')
sage: I = R.ideal(x, x+1)
sage: f = SingularKernelFunction("std")
sage: f(I)
[1]
Bases: sage.libs.singular.function.SingularFunction
EXAMPLES:
sage: from sage.libs.singular.function import SingularLibraryFunction
sage: R.<x,y> = PolynomialRing(QQ, order='lex')
sage: I = R.ideal(x, x+1)
sage: f = SingularLibraryFunction("groebner")
sage: f(I)
[1]
Tests for a sequence s, whether it consists of singular polynomials.
EXAMPLE:
sage: from sage.libs.singular.function import all_singular_poly_wrapper
sage: P.<x,y,z> = QQ[]
sage: all_singular_poly_wrapper([x+1, y])
True
sage: all_singular_poly_wrapper([x+1, y, 1])
False
Checks if a sequence s consists of free module elements over a singular ring.
EXAMPLE:
sage: from sage.libs.singular.function import all_vectors
sage: P.<x,y,z> = QQ[]
sage: M = P**2
sage: all_vectors([x])
False
sage: all_vectors([(x,y)])
False
sage: all_vectors([M(0), M((x,y))])
True
sage: all_vectors([M(0), M((x,y)),(0,0)])
False
Check whether wrapped ring arises from Singular or Singular/Plural.
EXAMPLE:
sage: from sage.libs.singular.function import is_sage_wrapper_for_singular_ring
sage: P.<x,y,z> = QQ[]
sage: is_sage_wrapper_for_singular_ring(P)
True
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
sage: P = A.g_algebra(relations={y*x:-x*y}, order = 'lex')
sage: is_sage_wrapper_for_singular_ring(P)
True
Checks if p is some data type corresponding to some singular poly.
EXAMPLE:
sage: from sage.libs.singular.function import is_singular_poly_wrapper
sage: A.<x,y,z> = FreeAlgebra(QQ, 3)
sage: H.<x,y,z> = A.g_algebra({z*x:x*z+2*x, z*y:y*z-2*y})
sage: is_singular_poly_wrapper(x+y)
True
Load the Singular library name.
INPUT:
EXAMPLE:
sage: from sage.libs.singular.function import singular_function
sage: from sage.libs.singular.function import lib as singular_lib
sage: singular_lib('general.lib')
sage: primes = singular_function('primes')
sage: primes(2,10, ring=GF(127)['x,y,z'])
(2, 3, 5, 7)
Return a list of all function names currently available.
INPUT:
EXAMPLE:
sage: 'groebner' in sage.libs.singular.function.list_of_functions()
True
Construct a new libSingular function object for the given name.
This function works both for interpreter and built-in functions.
INPUT:
EXAMPLES:
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: f = 3*x*y + 2*z + 1
sage: g = 2*x + 1/2
sage: I = Ideal([f,g])
sage: from sage.libs.singular.function import singular_function
sage: std = singular_function("std")
sage: std(I)
[3*y - 8*z - 4, 4*x + 1]
sage: size = singular_function("size")
sage: size([2, 3, 3], ring=P)
3
sage: size("sage", ring=P)
4
sage: size(["hello", "sage"], ring=P)
2
sage: factorize = singular_function("factorize")
sage: factorize(f)
[[1, 3*x*y + 2*z + 1], (1, 1)]
sage: factorize(f, 1)
[3*x*y + 2*z + 1]
We give a wrong number of arguments:
sage: factorize(ring=P)
Traceback (most recent call last):
...
RuntimeError: Error in Singular function call 'factorize':
Wrong number of arguments
sage: factorize(f, 1, 2)
Traceback (most recent call last):
...
RuntimeError: Error in Singular function call 'factorize':
Wrong number of arguments
sage: factorize(f, 1, 2, 3)
Traceback (most recent call last):
...
RuntimeError: Error in Singular function call 'factorize':
Wrong number of arguments
The Singular function list can be called with any number of arguments:
sage: singular_list = singular_function("list")
sage: singular_list(2, 3, 6, ring=P)
[2, 3, 6]
sage: singular_list(ring=P)
[]
sage: singular_list(1, ring=P)
[1]
sage: singular_list(1, 2, 3, 4, 5, 6, 7, 8, 9, 10, ring=P)
[1, 2, 3, 4, 5, 6, 7, 8, 9, 10]
We try to define a non-existing function:
sage: number_foobar = singular_function('number_foobar');
Traceback (most recent call last):
...
NameError: Function 'number_foobar' is not defined.
sage: from sage.libs.singular.function import lib as singular_lib
sage: singular_lib('general.lib')
sage: number_e = singular_function('number_e')
sage: number_e(10r,ring=P)
67957045707/25000000000
sage: RR(number_e(10r,ring=P))
2.71828182828000
sage: singular_lib('primdec.lib')
sage: primdecGTZ = singular_function("primdecGTZ")
sage: primdecGTZ(I)
[[[y - 8/3*z - 4/3, x + 1/4], [y - 8/3*z - 4/3, x + 1/4]]]
sage: singular_list((1,2,3),3,[1,2,3], ring=P)
[(1, 2, 3), 3, [1, 2, 3]]
sage: ringlist=singular_function("ringlist")
sage: l = ringlist(P)
sage: l[3].__class__
<class 'sage.rings.polynomial.multi_polynomial_sequence.PolynomialSequence_generic'>
sage: l
[0, ['x', 'y', 'z'], [['dp', (1, 1, 1)], ['C', (0,)]], [0]]
sage: ring=singular_function("ring")
sage: ring(l, ring=P)
<RingWrap>
sage: matrix = Matrix(P,2,2)
sage: matrix.randomize(terms=1)
sage: det = singular_function("det")
sage: det(matrix)
-3/5*x*y*z
sage: coeffs = singular_function("coeffs")
sage: coeffs(x*y+y+1,y)
[ 1]
[x + 1]
sage: F.<x,y,z> = GF(3)[]
sage: intmat = Matrix(ZZ, 2,2, [100,2,3,4])
sage: det(intmat, ring=F)
394
sage: random = singular_function("random")
sage: A = random(10,2,3, ring =F); A.nrows(), max(A.list()) <= 10
(2, True)
sage: P.<x,y,z> = PolynomialRing(QQ)
sage: M=P**3
sage: leadcoef = singular_function("leadcoef")
sage: v=M((100*x,5*y,10*z*x*y))
sage: leadcoef(v)
10
sage: v = M([x+y,x*y+y**3,z])
sage: lead = singular_function("lead")
sage: lead(v)
(0, y^3)
sage: jet = singular_function("jet")
sage: jet(v, 2)
(x + y, x*y, z)
sage: syz = singular_function("syz")
sage: I = P.ideal([x+y,x*y-y, y*2,x**2+1])
sage: M = syz(I)
sage: M
[(-2*y, 2, y + 1, 0), (0, -2, x - 1, 0), (x*y - y, -y + 1, 1, -y), (x^2 + 1, -x - 1, -1, -x)]
sage: singular_lib("mprimdec.lib")
sage: syz(M)
[(-x - 1, y - 1, 2*x, -2*y)]
sage: GTZmod = singular_function("GTZmod")
sage: GTZmod(M)
[[[(-2*y, 2, y + 1, 0), (0, x + 1, 1, -y), (0, -2, x - 1, 0), (x*y - y, -y + 1, 1, -y), (x^2 + 1, 0, 0, -x - y)], [0]]]
sage: mres = singular_function("mres")
sage: resolution = mres(M, 0)
sage: resolution
<Resolution>
sage: singular_list(resolution)
[[(-2*y, 2, y + 1, 0), (0, -2, x - 1, 0), (x*y - y, -y + 1, 1, -y), (x^2 + 1, -x - 1, -1, -x)], [(-x - 1, y - 1, 2*x, -2*y)], [(0)]]
sage: A.<x,y> = FreeAlgebra(QQ, 2)
sage: P.<x,y> = A.g_algebra({y*x:-x*y})
sage: I= Sequence([x*y,x+y], check=False, immutable=True)
sage: twostd = singular_function("twostd")
sage: twostd(I)
[x + y, y^2]
sage: M=syz(I)
doctest...
sage: M
[(x + y, x*y)]
sage: syz(M, ring=P)
[(0)]
sage: mres(I, 0)
<Resolution>
sage: M=P**3
sage: v=M((100*x,5*y,10*y*x*y))
sage: leadcoef(v)
-10
sage: v = M([x+y,x*y+y**3,x])
sage: lead(v)
(0, y^3)
sage: jet(v, 2)
(x + y, x*y, x)
sage: l = ringlist(P)
sage: len(l)
6
sage: ring(l, ring=P)
<noncommutative RingWrap>
sage: I=twostd(I)
sage: l[3]=I
sage: ring(l, ring=P)
<noncommutative RingWrap>