The coercion model

The coercion model manages how elements of one parent get related to elements of another. For example, the integer 2 can canonically be viewed as an element of the rational numbers. (The parent of a non-element is its Python type.)

sage: ZZ(2).parent()
Integer Ring
sage: QQ(2).parent()
Rational Field
>>> from sage.all import *
>>> ZZ(Integer(2)).parent()
Integer Ring
>>> QQ(Integer(2)).parent()
Rational Field
ZZ(2).parent()
QQ(2).parent()

The most prominent role of the coercion model is to make sense of binary operations between elements that have distinct parents. It does this by finding a parent where both elements make sense, and doing the operation there. For example:

sage: a = 1/2; a.parent()
Rational Field
sage: b = ZZ['x'].gen(); b.parent()
Univariate Polynomial Ring in x over Integer Ring
sage: a + b
x + 1/2
sage: (a + b).parent()
Univariate Polynomial Ring in x over Rational Field
>>> from sage.all import *
>>> a = Integer(1)/Integer(2); a.parent()
Rational Field
>>> b = ZZ['x'].gen(); b.parent()
Univariate Polynomial Ring in x over Integer Ring
>>> a + b
x + 1/2
>>> (a + b).parent()
Univariate Polynomial Ring in x over Rational Field
a = 1/2; a.parent()
b = ZZ['x'].gen(); b.parent()
a + b
(a + b).parent()

If there is a coercion (see below) from one of the parents to the other, the operation is always performed in the codomain of that coercion. Otherwise a reasonable attempt to create a new parent with coercion maps from both original parents is made. The results of these discoveries are cached. On failure, a TypeError is always raised.

Some arithmetic operations (such as multiplication) can indicate an action rather than arithmetic in a common parent. For example:

sage: E = EllipticCurve('37a')                                                      # needs sage.schemes
sage: P = E(0,0)                                                                    # needs sage.schemes
sage: 5*P                                                                           # needs sage.schemes
(1/4 : -5/8 : 1)
>>> from sage.all import *
>>> E = EllipticCurve('37a')                                                      # needs sage.schemes
>>> P = E(Integer(0),Integer(0))                                                                    # needs sage.schemes
>>> Integer(5)*P                                                                           # needs sage.schemes
(1/4 : -5/8 : 1)
E = EllipticCurve('37a')                                                      # needs sage.schemes
P = E(0,0)                                                                    # needs sage.schemes
5*P                                                                           # needs sage.schemes

where there is action of \(\ZZ\) on the points of \(E\) given by the additive group law. Parents can specify how they act on or are acted upon by other parents.

There are two kinds of ways to get from one parent to another, coercions and conversions.

Coercions are canonical (possibly modulo a finite number of deterministic choices) morphisms, and the set of all coercions between all parents forms a commuting diagram (modulo possibly rounding issues). \(\ZZ \rightarrow \QQ\) is an example of a coercion. These are invoked implicitly by the coercion model.

Conversions try to construct an element out of their input if at all possible. Examples include sections of coercions, creating an element from a string or list, etc. and may fail on some inputs of a given type while succeeding on others (i.e. they may not be defined on the whole domain). Conversions are always explicitly invoked, and never used by the coercion model to resolve binary operations.

For more information on how to specify coercions, conversions, and actions, see the documentation for Parent.

class sage.structure.coerce.CoercionModel[source]

Bases: object

See also sage.categories.pushout

EXAMPLES:

sage: f = ZZ['t', 'x'].0 + QQ['x'].0 + CyclotomicField(13).gen(); f             # needs sage.rings.number_field
t + x + zeta13
sage: f.parent()                                                                # needs sage.rings.number_field
Multivariate Polynomial Ring in t, x
 over Cyclotomic Field of order 13 and degree 12
sage: ZZ['x','y'].0 + ~Frac(QQ['y']).0
(x*y + 1)/y
sage: MatrixSpace(ZZ['x'], 2, 2)(2) + ~Frac(QQ['x']).0                          # needs sage.modules
[(2*x + 1)/x           0]
[          0 (2*x + 1)/x]
sage: f = ZZ['x,y,z'].0 + QQ['w,x,z,a'].0; f
w + x
sage: f.parent()
Multivariate Polynomial Ring in w, x, y, z, a over Rational Field
sage: ZZ['x,y,z'].0 + ZZ['w,x,z,a'].1
2*x
>>> from sage.all import *
>>> f = ZZ['t', 'x'].gen(0) + QQ['x'].gen(0) + CyclotomicField(Integer(13)).gen(); f             # needs sage.rings.number_field
t + x + zeta13
>>> f.parent()                                                                # needs sage.rings.number_field
Multivariate Polynomial Ring in t, x
 over Cyclotomic Field of order 13 and degree 12
>>> ZZ['x','y'].gen(0) + ~Frac(QQ['y']).gen(0)
(x*y + 1)/y
>>> MatrixSpace(ZZ['x'], Integer(2), Integer(2))(Integer(2)) + ~Frac(QQ['x']).gen(0)                          # needs sage.modules
[(2*x + 1)/x           0]
[          0 (2*x + 1)/x]
>>> f = ZZ['x,y,z'].gen(0) + QQ['w,x,z,a'].gen(0); f
w + x
>>> f.parent()
Multivariate Polynomial Ring in w, x, y, z, a over Rational Field
>>> ZZ['x,y,z'].gen(0) + ZZ['w,x,z,a'].gen(1)
2*x
f = ZZ['t', 'x'].0 + QQ['x'].0 + CyclotomicField(13).gen(); f             # needs sage.rings.number_field
f.parent()                                                                # needs sage.rings.number_field
ZZ['x','y'].0 + ~Frac(QQ['y']).0
MatrixSpace(ZZ['x'], 2, 2)(2) + ~Frac(QQ['x']).0                          # needs sage.modules
f = ZZ['x,y,z'].0 + QQ['w,x,z,a'].0; f
f.parent()
ZZ['x,y,z'].0 + ZZ['w,x,z,a'].1

AUTHOR:

  • Robert Bradshaw

analyse(xp, yp, op='mul')[source]

Emulate the process of doing arithmetic between xp and yp, returning a list of steps and the parent that the result will live in.

The explain() method is easier to use, but if one wants access to the actual morphism and action objects (rather than their string representations), then this is the function to use.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: GF7 = GF(7)
sage: steps, res = cm.analyse(GF7, ZZ)
sage: steps
['Coercion on right operand via',
 Natural morphism:
  From: Integer Ring
  To:   Finite Field of size 7,
 'Arithmetic performed after coercions.']
sage: res
Finite Field of size 7
sage: f = steps[1]; type(f)
<class 'sage.rings.finite_rings.integer_mod.Integer_to_IntegerMod'>
sage: f(100)
2
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> GF7 = GF(Integer(7))
>>> steps, res = cm.analyse(GF7, ZZ)
>>> steps
['Coercion on right operand via',
 Natural morphism:
  From: Integer Ring
  To:   Finite Field of size 7,
 'Arithmetic performed after coercions.']
>>> res
Finite Field of size 7
>>> f = steps[Integer(1)]; type(f)
<class 'sage.rings.finite_rings.integer_mod.Integer_to_IntegerMod'>
>>> f(Integer(100))
2
cm = sage.structure.element.get_coercion_model()
GF7 = GF(7)
steps, res = cm.analyse(GF7, ZZ)
steps
res
f = steps[1]; type(f)
f(100)
bin_op(x, y, op)[source]

Execute the operation op on \(x\) and \(y\).

It first looks for an action corresponding to op, and failing that, it tries to coerce \(x\) and \(y\) into a common parent and calls op on them.

If it cannot make sense of the operation, a TypeError is raised.

INPUT:

  • x – the left operand

  • y – the right operand

  • op – a python function taking 2 arguments

    Note

    op is often an arithmetic operation, but need not be so.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: cm.bin_op(1/2, 5, operator.mul)
5/2
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.bin_op(Integer(1)/Integer(2), Integer(5), operator.mul)
5/2
cm = sage.structure.element.get_coercion_model()
cm.bin_op(1/2, 5, operator.mul)

The operator can be any callable:

sage: R.<x> = ZZ['x']
sage: cm.bin_op(x^2 - 1, x + 1, gcd)
x + 1
>>> from sage.all import *
>>> R = ZZ['x']; (x,) = R._first_ngens(1)
>>> cm.bin_op(x**Integer(2) - Integer(1), x + Integer(1), gcd)
x + 1
R.<x> = ZZ['x']
cm.bin_op(x^2 - 1, x + 1, gcd)

Actions are detected and performed:

sage: M = matrix(ZZ, 2, 2, range(4))                                        # needs sage.modules
sage: V = vector(ZZ, [5,7])                                                 # needs sage.modules
sage: cm.bin_op(M, V, operator.mul)                                         # needs sage.modules
(7, 31)
>>> from sage.all import *
>>> M = matrix(ZZ, Integer(2), Integer(2), range(Integer(4)))                                        # needs sage.modules
>>> V = vector(ZZ, [Integer(5),Integer(7)])                                                 # needs sage.modules
>>> cm.bin_op(M, V, operator.mul)                                         # needs sage.modules
(7, 31)
M = matrix(ZZ, 2, 2, range(4))                                        # needs sage.modules
V = vector(ZZ, [5,7])                                                 # needs sage.modules
cm.bin_op(M, V, operator.mul)                                         # needs sage.modules
canonical_coercion(x, y)[source]

Given two elements \(x\) and \(y\), with parents \(S\) and \(R\) respectively, find a common parent \(Z\) such that there are coercions \(f: S \to Z\) and \(g: R \to Z\) and return \(f(x), g(y)\), which will have the same parent.

Raises a type error if no such \(Z\) can be found.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: cm.canonical_coercion(mod(2, 10), 17)
(2, 7)

sage: # needs sage.modules
sage: x, y = cm.canonical_coercion(1/2, matrix(ZZ, 2, 2, range(4)))
sage: x
[1/2   0]
[  0 1/2]
sage: y
[0 1]
[2 3]
sage: parent(x) is parent(y)
True
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.canonical_coercion(mod(Integer(2), Integer(10)), Integer(17))
(2, 7)

>>> # needs sage.modules
>>> x, y = cm.canonical_coercion(Integer(1)/Integer(2), matrix(ZZ, Integer(2), Integer(2), range(Integer(4))))
>>> x
[1/2   0]
[  0 1/2]
>>> y
[0 1]
[2 3]
>>> parent(x) is parent(y)
True
cm = sage.structure.element.get_coercion_model()
cm.canonical_coercion(mod(2, 10), 17)
# needs sage.modules
x, y = cm.canonical_coercion(1/2, matrix(ZZ, 2, 2, range(4)))
x
y
parent(x) is parent(y)

There is some support for non-Sage datatypes as well:

sage: x, y = cm.canonical_coercion(int(5), 10)
sage: type(x), type(y)
(<class 'sage.rings.integer.Integer'>, <class 'sage.rings.integer.Integer'>)

sage: x, y = cm.canonical_coercion(int(5), complex(3))
sage: type(x), type(y)
(<class 'complex'>, <class 'complex'>)

sage: class MyClass:
....:     def _sage_(self):
....:         return 13
sage: a, b = cm.canonical_coercion(MyClass(), 1/3)
sage: a, b
(13, 1/3)
sage: type(a)
<class 'sage.rings.rational.Rational'>
>>> from sage.all import *
>>> x, y = cm.canonical_coercion(int(Integer(5)), Integer(10))
>>> type(x), type(y)
(<class 'sage.rings.integer.Integer'>, <class 'sage.rings.integer.Integer'>)

>>> x, y = cm.canonical_coercion(int(Integer(5)), complex(Integer(3)))
>>> type(x), type(y)
(<class 'complex'>, <class 'complex'>)

>>> class MyClass:
...     def _sage_(self):
...         return Integer(13)
>>> a, b = cm.canonical_coercion(MyClass(), Integer(1)/Integer(3))
>>> a, b
(13, 1/3)
>>> type(a)
<class 'sage.rings.rational.Rational'>
x, y = cm.canonical_coercion(int(5), 10)
type(x), type(y)
x, y = cm.canonical_coercion(int(5), complex(3))
type(x), type(y)
class MyClass:
    def _sage_(self):
        return 13
a, b = cm.canonical_coercion(MyClass(), 1/3)
a, b
type(a)

We also make an exception for 0, even if \(\ZZ\) does not map in:

sage: canonical_coercion(vector([1, 2, 3]), 0)                              # needs sage.modules
((1, 2, 3), (0, 0, 0))
sage: canonical_coercion(GF(5)(0), float(0))
(0, 0)
>>> from sage.all import *
>>> canonical_coercion(vector([Integer(1), Integer(2), Integer(3)]), Integer(0))                              # needs sage.modules
((1, 2, 3), (0, 0, 0))
>>> canonical_coercion(GF(Integer(5))(Integer(0)), float(Integer(0)))
(0, 0)
canonical_coercion(vector([1, 2, 3]), 0)                              # needs sage.modules
canonical_coercion(GF(5)(0), float(0))
coercion_maps(R, S)[source]

Give two parents \(R\) and \(S\), return a pair of coercion maps \(f: R \rightarrow Z\) and \(g: S \rightarrow Z\) , if such a \(Z\) can be found.

In the (common) case that \(R=Z\) or \(S=Z\) then None is returned for \(f\) or \(g\) respectively rather than constructing (and subsequently calling) the identity morphism.

If no suitable \(f, g\) can be found, a single None is returned. This result is cached.

Note

By Issue #14711, coerce maps should be copied when using them outside of the coercion system, because they may become defunct by garbage collection.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: f, g = cm.coercion_maps(ZZ, QQ)
sage: print(copy(f))
Natural morphism:
  From: Integer Ring
  To:   Rational Field
sage: print(g)
None

sage: ZZx = ZZ['x']
sage: f, g = cm.coercion_maps(ZZx, QQ)
sage: print(f)
(map internal to coercion system -- copy before use)
Ring morphism:
  From: Univariate Polynomial Ring in x over Integer Ring
  To:   Univariate Polynomial Ring in x over Rational Field
sage: print(g)
(map internal to coercion system -- copy before use)
Polynomial base injection morphism:
  From: Rational Field
  To:   Univariate Polynomial Ring in x over Rational Field

sage: K = GF(7)
sage: cm.coercion_maps(QQ, K) is None
True
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> f, g = cm.coercion_maps(ZZ, QQ)
>>> print(copy(f))
Natural morphism:
  From: Integer Ring
  To:   Rational Field
>>> print(g)
None

>>> ZZx = ZZ['x']
>>> f, g = cm.coercion_maps(ZZx, QQ)
>>> print(f)
(map internal to coercion system -- copy before use)
Ring morphism:
  From: Univariate Polynomial Ring in x over Integer Ring
  To:   Univariate Polynomial Ring in x over Rational Field
>>> print(g)
(map internal to coercion system -- copy before use)
Polynomial base injection morphism:
  From: Rational Field
  To:   Univariate Polynomial Ring in x over Rational Field

>>> K = GF(Integer(7))
>>> cm.coercion_maps(QQ, K) is None
True
cm = sage.structure.element.get_coercion_model()
f, g = cm.coercion_maps(ZZ, QQ)
print(copy(f))
print(g)
ZZx = ZZ['x']
f, g = cm.coercion_maps(ZZx, QQ)
print(f)
print(g)
K = GF(7)
cm.coercion_maps(QQ, K) is None

Note that to break symmetry, if there is a coercion map in both directions, the parent on the left is used:

sage: # needs sage.modules
sage: V = QQ^3
sage: W = V.__class__(QQ, 3)
sage: V == W
True
sage: V is W
False
sage: cm = sage.structure.element.get_coercion_model()
sage: cm.coercion_maps(V, W)
(None,
 (map internal to coercion system -- copy before use)
 Coercion map:
   From: Vector space of dimension 3 over Rational Field
   To:   Vector space of dimension 3 over Rational Field)
sage: cm.coercion_maps(W, V)
(None,
 (map internal to coercion system -- copy before use)
 Coercion map:
   From: Vector space of dimension 3 over Rational Field
   To:   Vector space of dimension 3 over Rational Field)
sage: v = V([1,2,3])
sage: w = W([1,2,3])
sage: parent(v + w) is V
True
sage: parent(w + v) is W
True
>>> from sage.all import *
>>> # needs sage.modules
>>> V = QQ**Integer(3)
>>> W = V.__class__(QQ, Integer(3))
>>> V == W
True
>>> V is W
False
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.coercion_maps(V, W)
(None,
 (map internal to coercion system -- copy before use)
 Coercion map:
   From: Vector space of dimension 3 over Rational Field
   To:   Vector space of dimension 3 over Rational Field)
>>> cm.coercion_maps(W, V)
(None,
 (map internal to coercion system -- copy before use)
 Coercion map:
   From: Vector space of dimension 3 over Rational Field
   To:   Vector space of dimension 3 over Rational Field)
>>> v = V([Integer(1),Integer(2),Integer(3)])
>>> w = W([Integer(1),Integer(2),Integer(3)])
>>> parent(v + w) is V
True
>>> parent(w + v) is W
True
# needs sage.modules
V = QQ^3
W = V.__class__(QQ, 3)
V == W
V is W
cm = sage.structure.element.get_coercion_model()
cm.coercion_maps(V, W)
cm.coercion_maps(W, V)
v = V([1,2,3])
w = W([1,2,3])
parent(v + w) is V
parent(w + v) is W
common_parent(*args)[source]

Compute a common parent for all the inputs.

It’s essentially an \(n\)-ary canonical coercion except it can operate on parents rather than just elements.

INPUT:

  • args – set of elements and/or parents

OUTPUT:

A Parent into which each input should coerce, or raises a TypeError if no such Parent can be found.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: cm.common_parent(ZZ, QQ)
Rational Field
sage: cm.common_parent(ZZ, QQ, RR)                                          # needs sage.rings.real_mpfr
Real Field with 53 bits of precision
sage: ZZT = ZZ[['T']]
sage: QQT = QQ['T']
sage: cm.common_parent(ZZT, QQT, RDF)
Power Series Ring in T over Real Double Field
sage: cm.common_parent(4r, 5r)
<class 'int'>
sage: cm.common_parent(int, float, ZZ)
<class 'float'>
sage: real_fields = [RealField(prec) for prec in [10,20..100]]              # needs sage.rings.real_mpfr
sage: cm.common_parent(*real_fields)                                        # needs sage.rings.real_mpfr
Real Field with 10 bits of precision
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.common_parent(ZZ, QQ)
Rational Field
>>> cm.common_parent(ZZ, QQ, RR)                                          # needs sage.rings.real_mpfr
Real Field with 53 bits of precision
>>> ZZT = ZZ[['T']]
>>> QQT = QQ['T']
>>> cm.common_parent(ZZT, QQT, RDF)
Power Series Ring in T over Real Double Field
>>> cm.common_parent(4, 5)
<class 'int'>
>>> cm.common_parent(int, float, ZZ)
<class 'float'>
>>> real_fields = [RealField(prec) for prec in (ellipsis_range(Integer(10),Integer(20),Ellipsis,Integer(100)))]              # needs sage.rings.real_mpfr
>>> cm.common_parent(*real_fields)                                        # needs sage.rings.real_mpfr
Real Field with 10 bits of precision
cm = sage.structure.element.get_coercion_model()
cm.common_parent(ZZ, QQ)
cm.common_parent(ZZ, QQ, RR)                                          # needs sage.rings.real_mpfr
ZZT = ZZ[['T']]
QQT = QQ['T']
cm.common_parent(ZZT, QQT, RDF)
cm.common_parent(4r, 5r)
cm.common_parent(int, float, ZZ)
real_fields = [RealField(prec) for prec in [10,20..100]]              # needs sage.rings.real_mpfr
cm.common_parent(*real_fields)                                        # needs sage.rings.real_mpfr

There are some cases where the ordering does matter, but if a parent can be found it is always the same:

sage: QQxy = QQ['x,y']
sage: QQyz = QQ['y,z']
sage: cm.common_parent(QQxy, QQyz) == cm.common_parent(QQyz, QQxy)
True
sage: QQzt = QQ['z,t']
sage: cm.common_parent(QQxy, QQyz, QQzt)
Multivariate Polynomial Ring in x, y, z, t over Rational Field
sage: cm.common_parent(QQxy, QQzt, QQyz)
Traceback (most recent call last):
...
TypeError: no common canonical parent for objects with parents:
'Multivariate Polynomial Ring in x, y over Rational Field' and
'Multivariate Polynomial Ring in z, t over Rational Field'
>>> from sage.all import *
>>> QQxy = QQ['x,y']
>>> QQyz = QQ['y,z']
>>> cm.common_parent(QQxy, QQyz) == cm.common_parent(QQyz, QQxy)
True
>>> QQzt = QQ['z,t']
>>> cm.common_parent(QQxy, QQyz, QQzt)
Multivariate Polynomial Ring in x, y, z, t over Rational Field
>>> cm.common_parent(QQxy, QQzt, QQyz)
Traceback (most recent call last):
...
TypeError: no common canonical parent for objects with parents:
'Multivariate Polynomial Ring in x, y over Rational Field' and
'Multivariate Polynomial Ring in z, t over Rational Field'
QQxy = QQ['x,y']
QQyz = QQ['y,z']
cm.common_parent(QQxy, QQyz) == cm.common_parent(QQyz, QQxy)
QQzt = QQ['z,t']
cm.common_parent(QQxy, QQyz, QQzt)
cm.common_parent(QQxy, QQzt, QQyz)
discover_action(R, S, op, r=None, s=None)[source]

INPUT:

  • R – the left Parent (or type)

  • S – the right Parent (or type)

  • op – the operand, typically an element of the operator module

  • r – (optional) element of \(R\)

  • s – (optional) element of \(S\)

OUTPUT: an action \(A\) such that \(s\) op \(r\) is given by \(A(s,r)\)

The steps taken are illustrated below.

EXAMPLES:

sage: P.<x> = ZZ['x']
sage: P.get_action(ZZ)
Right scalar multiplication by Integer Ring on
 Univariate Polynomial Ring in x over Integer Ring
sage: ZZ.get_action(P) is None
True
sage: cm = sage.structure.element.get_coercion_model()
>>> from sage.all import *
>>> P = ZZ['x']; (x,) = P._first_ngens(1)
>>> P.get_action(ZZ)
Right scalar multiplication by Integer Ring on
 Univariate Polynomial Ring in x over Integer Ring
>>> ZZ.get_action(P) is None
True
>>> cm = sage.structure.element.get_coercion_model()
P.<x> = ZZ['x']
P.get_action(ZZ)
ZZ.get_action(P) is None
cm = sage.structure.element.get_coercion_model()

If \(R\) or \(S\) is a Parent, ask it for an action by/on \(R\):

sage: cm.discover_action(ZZ, P, operator.mul)
Left scalar multiplication by Integer Ring on
 Univariate Polynomial Ring in x over Integer Ring
>>> from sage.all import *
>>> cm.discover_action(ZZ, P, operator.mul)
Left scalar multiplication by Integer Ring on
 Univariate Polynomial Ring in x over Integer Ring
cm.discover_action(ZZ, P, operator.mul)

If \(R\) or \(S\) a type, recursively call get_action() with the Sage versions of \(R\) and/or \(S\):

sage: cm.discover_action(P, int, operator.mul)
Right scalar multiplication by Integer Ring on
 Univariate Polynomial Ring in x over Integer Ring
 with precomposition on right by Native morphism:
  From: Set of Python objects of class 'int'
  To:   Integer Ring
>>> from sage.all import *
>>> cm.discover_action(P, int, operator.mul)
Right scalar multiplication by Integer Ring on
 Univariate Polynomial Ring in x over Integer Ring
 with precomposition on right by Native morphism:
  From: Set of Python objects of class 'int'
  To:   Integer Ring
cm.discover_action(P, int, operator.mul)

If op is division, look for action on right by inverse:

sage: cm.discover_action(P, ZZ, operator.truediv)
Right inverse action by Rational Field on
 Univariate Polynomial Ring in x over Integer Ring
with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Rational Field
>>> from sage.all import *
>>> cm.discover_action(P, ZZ, operator.truediv)
Right inverse action by Rational Field on
 Univariate Polynomial Ring in x over Integer Ring
with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Rational Field
cm.discover_action(P, ZZ, operator.truediv)

Check that Issue #17740 is fixed:

sage: R = GF(5)['x']
sage: cm.discover_action(R, ZZ, operator.truediv)
Right inverse action by Finite Field of size 5
 on Univariate Polynomial Ring in x over Finite Field of size 5
with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Finite Field of size 5
sage: cm.bin_op(R.gen(), 7, operator.truediv).parent()
Univariate Polynomial Ring in x over Finite Field of size 5
>>> from sage.all import *
>>> R = GF(Integer(5))['x']
>>> cm.discover_action(R, ZZ, operator.truediv)
Right inverse action by Finite Field of size 5
 on Univariate Polynomial Ring in x over Finite Field of size 5
with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Finite Field of size 5
>>> cm.bin_op(R.gen(), Integer(7), operator.truediv).parent()
Univariate Polynomial Ring in x over Finite Field of size 5
R = GF(5)['x']
cm.discover_action(R, ZZ, operator.truediv)
cm.bin_op(R.gen(), 7, operator.truediv).parent()

Check that Issue #18221 is fixed:

sage: # needs sage.combinat sage.modules
sage: F.<x> = FreeAlgebra(QQ)
sage: x / 2
1/2*x
sage: cm.discover_action(F, ZZ, operator.truediv)
Right inverse action by Rational Field on
 Free Algebra on 1 generators (x,) over Rational Field
 with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Rational Field
>>> from sage.all import *
>>> # needs sage.combinat sage.modules
>>> F = FreeAlgebra(QQ, names=('x',)); (x,) = F._first_ngens(1)
>>> x / Integer(2)
1/2*x
>>> cm.discover_action(F, ZZ, operator.truediv)
Right inverse action by Rational Field on
 Free Algebra on 1 generators (x,) over Rational Field
 with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Rational Field
# needs sage.combinat sage.modules
F.<x> = FreeAlgebra(QQ)
x / 2
cm.discover_action(F, ZZ, operator.truediv)
discover_coercion(R, S)[source]

This actually implements the finding of coercion maps as described in the coercion_maps() method.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
cm = sage.structure.element.get_coercion_model()

If R is S, then two identity morphisms suffice:

sage: cm.discover_coercion(SR, SR)                                          # needs sage.symbolic
(None, None)
>>> from sage.all import *
>>> cm.discover_coercion(SR, SR)                                          # needs sage.symbolic
(None, None)
cm.discover_coercion(SR, SR)                                          # needs sage.symbolic

If there is a coercion map either direction, use that:

sage: cm.discover_coercion(ZZ, QQ)
((map internal to coercion system -- copy before use)
Natural morphism:
  From: Integer Ring
  To:   Rational Field, None)
sage: cm.discover_coercion(RR, QQ)                                          # needs sage.rings.real_mpfr
(None, (map internal to coercion system -- copy before use)
 Generic map:
  From: Rational Field
  To:   Real Field with 53 bits of precision)
>>> from sage.all import *
>>> cm.discover_coercion(ZZ, QQ)
((map internal to coercion system -- copy before use)
Natural morphism:
  From: Integer Ring
  To:   Rational Field, None)
>>> cm.discover_coercion(RR, QQ)                                          # needs sage.rings.real_mpfr
(None, (map internal to coercion system -- copy before use)
 Generic map:
  From: Rational Field
  To:   Real Field with 53 bits of precision)
cm.discover_coercion(ZZ, QQ)
cm.discover_coercion(RR, QQ)                                          # needs sage.rings.real_mpfr

Otherwise, try and compute an appropriate cover:

sage: ZZxy = ZZ['x,y']
sage: cm.discover_coercion(ZZxy, RDF)
((map internal to coercion system -- copy before use)
Call morphism:
  From: Multivariate Polynomial Ring in x, y over Integer Ring
  To:   Multivariate Polynomial Ring in x, y over Real Double Field,
 (map internal to coercion system -- copy before use)
 Polynomial base injection morphism:
  From: Real Double Field
  To:   Multivariate Polynomial Ring in x, y over Real Double Field)
>>> from sage.all import *
>>> ZZxy = ZZ['x,y']
>>> cm.discover_coercion(ZZxy, RDF)
((map internal to coercion system -- copy before use)
Call morphism:
  From: Multivariate Polynomial Ring in x, y over Integer Ring
  To:   Multivariate Polynomial Ring in x, y over Real Double Field,
 (map internal to coercion system -- copy before use)
 Polynomial base injection morphism:
  From: Real Double Field
  To:   Multivariate Polynomial Ring in x, y over Real Double Field)
ZZxy = ZZ['x,y']
cm.discover_coercion(ZZxy, RDF)

Sometimes there is a reasonable “cover,” but no canonical coercion:

sage: sage.categories.pushout.pushout(QQ, QQ^3)                             # needs sage.modules
Vector space of dimension 3 over Rational Field
sage: print(cm.discover_coercion(QQ, QQ^3))                                 # needs sage.modules
None
>>> from sage.all import *
>>> sage.categories.pushout.pushout(QQ, QQ**Integer(3))                             # needs sage.modules
Vector space of dimension 3 over Rational Field
>>> print(cm.discover_coercion(QQ, QQ**Integer(3)))                                 # needs sage.modules
None
sage.categories.pushout.pushout(QQ, QQ^3)                             # needs sage.modules
print(cm.discover_coercion(QQ, QQ^3))                                 # needs sage.modules
division_parent(P)[source]

Deduces where the result of division in P lies by calculating the inverse of P.one() or P.an_element().

The result is cached.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: cm.division_parent(ZZ)
Rational Field
sage: cm.division_parent(QQ)
Rational Field
sage: ZZx = ZZ['x']
sage: cm.division_parent(ZZx)
Fraction Field of Univariate Polynomial Ring in x over Integer Ring
sage: K = GF(41)
sage: cm.division_parent(K)
Finite Field of size 41
sage: Zmod100 = Integers(100)
sage: cm.division_parent(Zmod100)
Ring of integers modulo 100
sage: S5 = SymmetricGroup(5)                                                # needs sage.groups
sage: cm.division_parent(S5)                                                # needs sage.groups
Symmetric group of order 5! as a permutation group
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.division_parent(ZZ)
Rational Field
>>> cm.division_parent(QQ)
Rational Field
>>> ZZx = ZZ['x']
>>> cm.division_parent(ZZx)
Fraction Field of Univariate Polynomial Ring in x over Integer Ring
>>> K = GF(Integer(41))
>>> cm.division_parent(K)
Finite Field of size 41
>>> Zmod100 = Integers(Integer(100))
>>> cm.division_parent(Zmod100)
Ring of integers modulo 100
>>> S5 = SymmetricGroup(Integer(5))                                                # needs sage.groups
>>> cm.division_parent(S5)                                                # needs sage.groups
Symmetric group of order 5! as a permutation group
cm = sage.structure.element.get_coercion_model()
cm.division_parent(ZZ)
cm.division_parent(QQ)
ZZx = ZZ['x']
cm.division_parent(ZZx)
K = GF(41)
cm.division_parent(K)
Zmod100 = Integers(100)
cm.division_parent(Zmod100)
S5 = SymmetricGroup(5)                                                # needs sage.groups
cm.division_parent(S5)                                                # needs sage.groups
exception_stack()[source]

Return the list of exceptions that were caught in the course of executing the last binary operation. Useful for diagnosis when user-defined maps or actions raise exceptions that are caught in the course of coercion detection.

If all went well, this should be the empty list. If things aren’t happening as you expect, this is a good place to check. See also coercion_traceback().

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: cm.record_exceptions()
sage: 1/2 + 2
5/2
sage: cm.exception_stack()
[]
sage: 1/2 + GF(3)(2)
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for +:
'Rational Field' and 'Finite Field of size 3'
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.record_exceptions()
>>> Integer(1)/Integer(2) + Integer(2)
5/2
>>> cm.exception_stack()
[]
>>> Integer(1)/Integer(2) + GF(Integer(3))(Integer(2))
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for +:
'Rational Field' and 'Finite Field of size 3'
cm = sage.structure.element.get_coercion_model()
cm.record_exceptions()
1/2 + 2
cm.exception_stack()
1/2 + GF(3)(2)

Now see what the actual problem was:

sage: import traceback
sage: cm.exception_stack()
['Traceback (most recent call last):...', 'Traceback (most recent call last):...']
sage: print(cm.exception_stack()[-1])
Traceback (most recent call last):
...
TypeError: no common canonical parent for objects with parents:
'Rational Field' and 'Finite Field of size 3'
>>> from sage.all import *
>>> import traceback
>>> cm.exception_stack()
['Traceback (most recent call last):...', 'Traceback (most recent call last):...']
>>> print(cm.exception_stack()[-Integer(1)])
Traceback (most recent call last):
...
TypeError: no common canonical parent for objects with parents:
'Rational Field' and 'Finite Field of size 3'
import traceback
cm.exception_stack()
print(cm.exception_stack()[-1])

This is typically accessed via the coercion_traceback() function.

sage: coercion_traceback()
Traceback (most recent call last):
...
TypeError: no common canonical parent for objects with parents:
'Rational Field' and 'Finite Field of size 3'
>>> from sage.all import *
>>> coercion_traceback()
Traceback (most recent call last):
...
TypeError: no common canonical parent for objects with parents:
'Rational Field' and 'Finite Field of size 3'
coercion_traceback()
explain(xp, yp, op='mul', verbosity=2)[source]

This function can be used to understand what coercions will happen for an arithmetic operation between xp and yp (which may be either elements or parents). If the parent of the result can be determined then it will be returned.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()

sage: cm.explain(ZZ, ZZ)
Identical parents, arithmetic performed immediately.
Result lives in Integer Ring
Integer Ring

sage: cm.explain(QQ, int)
Coercion on right operand via
    Native morphism:
      From: Set of Python objects of class 'int'
      To:   Rational Field
Arithmetic performed after coercions.
Result lives in Rational Field
Rational Field

sage: R = ZZ['x']
sage: cm.explain(R, QQ)
Action discovered.
    Right scalar multiplication by Rational Field
     on Univariate Polynomial Ring in x over Integer Ring
Result lives in Univariate Polynomial Ring in x over Rational Field
Univariate Polynomial Ring in x over Rational Field

sage: cm.explain(ZZ['x'], QQ, operator.add)
Coercion on left operand via
    Ring morphism:
      From: Univariate Polynomial Ring in x over Integer Ring
      To:   Univariate Polynomial Ring in x over Rational Field
      Defn: Induced from base ring by
            Natural morphism:
              From: Integer Ring
              To:   Rational Field
Coercion on right operand via
    Polynomial base injection morphism:
      From: Rational Field
      To:   Univariate Polynomial Ring in x over Rational Field
Arithmetic performed after coercions.
Result lives in Univariate Polynomial Ring in x over Rational Field
Univariate Polynomial Ring in x over Rational Field
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()

>>> cm.explain(ZZ, ZZ)
Identical parents, arithmetic performed immediately.
Result lives in Integer Ring
Integer Ring

>>> cm.explain(QQ, int)
Coercion on right operand via
    Native morphism:
      From: Set of Python objects of class 'int'
      To:   Rational Field
Arithmetic performed after coercions.
Result lives in Rational Field
Rational Field

>>> R = ZZ['x']
>>> cm.explain(R, QQ)
Action discovered.
    Right scalar multiplication by Rational Field
     on Univariate Polynomial Ring in x over Integer Ring
Result lives in Univariate Polynomial Ring in x over Rational Field
Univariate Polynomial Ring in x over Rational Field

>>> cm.explain(ZZ['x'], QQ, operator.add)
Coercion on left operand via
    Ring morphism:
      From: Univariate Polynomial Ring in x over Integer Ring
      To:   Univariate Polynomial Ring in x over Rational Field
      Defn: Induced from base ring by
            Natural morphism:
              From: Integer Ring
              To:   Rational Field
Coercion on right operand via
    Polynomial base injection morphism:
      From: Rational Field
      To:   Univariate Polynomial Ring in x over Rational Field
Arithmetic performed after coercions.
Result lives in Univariate Polynomial Ring in x over Rational Field
Univariate Polynomial Ring in x over Rational Field
cm = sage.structure.element.get_coercion_model()
cm.explain(ZZ, ZZ)
cm.explain(QQ, int)
R = ZZ['x']
cm.explain(R, QQ)
cm.explain(ZZ['x'], QQ, operator.add)

Sometimes with non-sage types there is not enough information to deduce what will actually happen:

sage: R100 = RealField(100)                                                 # needs sage.rings.real_mpfr
sage: cm.explain(R100, float, operator.add)                                 # needs sage.rings.real_mpfr
Right operand is numeric, will attempt coercion in both directions.
Unknown result parent.
sage: parent(R100(1) + float(1))                                            # needs sage.rings.real_mpfr
<class 'float'>
sage: cm.explain(QQ, float, operator.add)
Right operand is numeric, will attempt coercion in both directions.
Unknown result parent.
sage: parent(QQ(1) + float(1))
<class 'float'>
>>> from sage.all import *
>>> R100 = RealField(Integer(100))                                                 # needs sage.rings.real_mpfr
>>> cm.explain(R100, float, operator.add)                                 # needs sage.rings.real_mpfr
Right operand is numeric, will attempt coercion in both directions.
Unknown result parent.
>>> parent(R100(Integer(1)) + float(Integer(1)))                                            # needs sage.rings.real_mpfr
<class 'float'>
>>> cm.explain(QQ, float, operator.add)
Right operand is numeric, will attempt coercion in both directions.
Unknown result parent.
>>> parent(QQ(Integer(1)) + float(Integer(1)))
<class 'float'>
R100 = RealField(100)                                                 # needs sage.rings.real_mpfr
cm.explain(R100, float, operator.add)                                 # needs sage.rings.real_mpfr
parent(R100(1) + float(1))                                            # needs sage.rings.real_mpfr
cm.explain(QQ, float, operator.add)
parent(QQ(1) + float(1))

Special care is taken to deal with division:

sage: cm.explain(ZZ, ZZ, operator.truediv)
Identical parents, arithmetic performed immediately.
Result lives in Rational Field
Rational Field

sage: ZZx = ZZ['x']
sage: QQx = QQ['x']
sage: cm.explain(ZZx, QQx, operator.truediv)
Coercion on left operand via
    Ring morphism:
      From: Univariate Polynomial Ring in x over Integer Ring
      To:   Univariate Polynomial Ring in x over Rational Field
      Defn: Induced from base ring by
            Natural morphism:
              From: Integer Ring
              To:   Rational Field
Arithmetic performed after coercions.
Result lives in Fraction Field of Univariate Polynomial Ring in x over Rational Field
Fraction Field of Univariate Polynomial Ring in x over Rational Field

sage: cm.explain(int, ZZ, operator.truediv)
Coercion on left operand via
    Native morphism:
      From: Set of Python objects of class 'int'
      To:   Integer Ring
Arithmetic performed after coercions.
Result lives in Rational Field
Rational Field

sage: cm.explain(ZZx, ZZ, operator.truediv)
Action discovered.
    Right inverse action by Rational Field
     on Univariate Polynomial Ring in x over Integer Ring
    with precomposition on right by Natural morphism:
      From: Integer Ring
      To:   Rational Field
Result lives in Univariate Polynomial Ring in x over Rational Field
Univariate Polynomial Ring in x over Rational Field
>>> from sage.all import *
>>> cm.explain(ZZ, ZZ, operator.truediv)
Identical parents, arithmetic performed immediately.
Result lives in Rational Field
Rational Field

>>> ZZx = ZZ['x']
>>> QQx = QQ['x']
>>> cm.explain(ZZx, QQx, operator.truediv)
Coercion on left operand via
    Ring morphism:
      From: Univariate Polynomial Ring in x over Integer Ring
      To:   Univariate Polynomial Ring in x over Rational Field
      Defn: Induced from base ring by
            Natural morphism:
              From: Integer Ring
              To:   Rational Field
Arithmetic performed after coercions.
Result lives in Fraction Field of Univariate Polynomial Ring in x over Rational Field
Fraction Field of Univariate Polynomial Ring in x over Rational Field

>>> cm.explain(int, ZZ, operator.truediv)
Coercion on left operand via
    Native morphism:
      From: Set of Python objects of class 'int'
      To:   Integer Ring
Arithmetic performed after coercions.
Result lives in Rational Field
Rational Field

>>> cm.explain(ZZx, ZZ, operator.truediv)
Action discovered.
    Right inverse action by Rational Field
     on Univariate Polynomial Ring in x over Integer Ring
    with precomposition on right by Natural morphism:
      From: Integer Ring
      To:   Rational Field
Result lives in Univariate Polynomial Ring in x over Rational Field
Univariate Polynomial Ring in x over Rational Field
cm.explain(ZZ, ZZ, operator.truediv)
ZZx = ZZ['x']
QQx = QQ['x']
cm.explain(ZZx, QQx, operator.truediv)
cm.explain(int, ZZ, operator.truediv)
cm.explain(ZZx, ZZ, operator.truediv)

Note

This function is accurate only in so far as analyse() is kept in sync with the bin_op() and canonical_coercion() which are kept separate for maximal efficiency.

get_action(R, S, op='mul', r=None, s=None)[source]

Get the action of R on S or S on R associated to the operation op.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: ZZx = ZZ['x']
sage: cm.get_action(ZZx, ZZ, operator.mul)
Right scalar multiplication by Integer Ring
 on Univariate Polynomial Ring in x over Integer Ring
sage: cm.get_action(ZZx, QQ, operator.mul)
Right scalar multiplication by Rational Field
 on Univariate Polynomial Ring in x over Integer Ring
sage: QQx = QQ['x']
sage: cm.get_action(QQx, int, operator.mul)
Right scalar multiplication by Integer Ring
 on Univariate Polynomial Ring in x over Rational Field
with precomposition on right by Native morphism:
  From: Set of Python objects of class 'int'
  To:   Integer Ring

sage: A = cm.get_action(QQx, ZZ, operator.truediv); A
Right inverse action by Rational Field
 on Univariate Polynomial Ring in x over Rational Field
with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Rational Field
sage: x = QQx.gen()
sage: A(x+10, 5)
1/5*x + 2
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> ZZx = ZZ['x']
>>> cm.get_action(ZZx, ZZ, operator.mul)
Right scalar multiplication by Integer Ring
 on Univariate Polynomial Ring in x over Integer Ring
>>> cm.get_action(ZZx, QQ, operator.mul)
Right scalar multiplication by Rational Field
 on Univariate Polynomial Ring in x over Integer Ring
>>> QQx = QQ['x']
>>> cm.get_action(QQx, int, operator.mul)
Right scalar multiplication by Integer Ring
 on Univariate Polynomial Ring in x over Rational Field
with precomposition on right by Native morphism:
  From: Set of Python objects of class 'int'
  To:   Integer Ring

>>> A = cm.get_action(QQx, ZZ, operator.truediv); A
Right inverse action by Rational Field
 on Univariate Polynomial Ring in x over Rational Field
with precomposition on right by Natural morphism:
  From: Integer Ring
  To:   Rational Field
>>> x = QQx.gen()
>>> A(x+Integer(10), Integer(5))
1/5*x + 2
cm = sage.structure.element.get_coercion_model()
ZZx = ZZ['x']
cm.get_action(ZZx, ZZ, operator.mul)
cm.get_action(ZZx, QQ, operator.mul)
QQx = QQ['x']
cm.get_action(QQx, int, operator.mul)
A = cm.get_action(QQx, ZZ, operator.truediv); A
x = QQx.gen()
A(x+10, 5)
get_cache()[source]

This returns the current cache of coercion maps and actions, primarily useful for debugging and introspection.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: cm.canonical_coercion(1, 2/3)
(1, 2/3)
sage: maps, actions = cm.get_cache()
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.canonical_coercion(Integer(1), Integer(2)/Integer(3))
(1, 2/3)
>>> maps, actions = cm.get_cache()
cm = sage.structure.element.get_coercion_model()
cm.canonical_coercion(1, 2/3)
maps, actions = cm.get_cache()

Now let us see what happens when we do a binary operations with an integer and a rational:

sage: left_morphism_ref, right_morphism_ref = maps[ZZ, QQ]
>>> from sage.all import *
>>> left_morphism_ref, right_morphism_ref = maps[ZZ, QQ]
left_morphism_ref, right_morphism_ref = maps[ZZ, QQ]

Note that by Issue #14058 the coercion model only stores a weak reference to the coercion maps in this case:

sage: left_morphism_ref
<weakref at ...; to 'sage.rings.rational.Z_to_Q' at ...>
>>> from sage.all import *
>>> left_morphism_ref
<weakref at ...; to 'sage.rings.rational.Z_to_Q' at ...>
left_morphism_ref

Moreover, the weakly referenced coercion map uses only a weak reference to the codomain:

sage: left_morphism_ref()
(map internal to coercion system -- copy before use)
Natural morphism:
  From: Integer Ring
  To:   Rational Field
>>> from sage.all import *
>>> left_morphism_ref()
(map internal to coercion system -- copy before use)
Natural morphism:
  From: Integer Ring
  To:   Rational Field
left_morphism_ref()

To get an actual valid map, we simply copy the weakly referenced coercion map:

sage: print(copy(left_morphism_ref()))
Natural morphism:
  From: Integer Ring
  To:   Rational Field
sage: print(right_morphism_ref)
None
>>> from sage.all import *
>>> print(copy(left_morphism_ref()))
Natural morphism:
  From: Integer Ring
  To:   Rational Field
>>> print(right_morphism_ref)
None
print(copy(left_morphism_ref()))
print(right_morphism_ref)

We can see that it coerces the left operand from an integer to a rational, and doesn’t do anything to the right.

Now for some actions:

sage: R.<x> = ZZ['x']
sage: 1/2 * x
1/2*x
sage: maps, actions = cm.get_cache()
sage: act = actions[QQ, R, operator.mul]; act
Left scalar multiplication by Rational Field
 on Univariate Polynomial Ring in x over Integer Ring
sage: act.actor()
Rational Field
sage: act.domain()
Univariate Polynomial Ring in x over Integer Ring
sage: act.codomain()
Univariate Polynomial Ring in x over Rational Field
sage: act(1/5, x+10)
1/5*x + 2
>>> from sage.all import *
>>> R = ZZ['x']; (x,) = R._first_ngens(1)
>>> Integer(1)/Integer(2) * x
1/2*x
>>> maps, actions = cm.get_cache()
>>> act = actions[QQ, R, operator.mul]; act
Left scalar multiplication by Rational Field
 on Univariate Polynomial Ring in x over Integer Ring
>>> act.actor()
Rational Field
>>> act.domain()
Univariate Polynomial Ring in x over Integer Ring
>>> act.codomain()
Univariate Polynomial Ring in x over Rational Field
>>> act(Integer(1)/Integer(5), x+Integer(10))
1/5*x + 2
R.<x> = ZZ['x']
1/2 * x
maps, actions = cm.get_cache()
act = actions[QQ, R, operator.mul]; act
act.actor()
act.domain()
act.codomain()
act(1/5, x+10)
record_exceptions(value=True)[source]

Enables (or disables) recording of the exceptions suppressed during arithmetic.

Each time that record_exceptions is called (either enabling or disabling the record), the exception_stack is cleared.

reset_cache()[source]

Clear the coercion cache.

This should have no impact on the result of arithmetic operations, as the exact same coercions and actions will be re-discovered when needed.

It may be useful for debugging, and may also free some memory.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: len(cm.get_cache()[0])    # random
42
sage: cm.reset_cache()
sage: cm.get_cache()
({}, {})
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> len(cm.get_cache()[Integer(0)])    # random
42
>>> cm.reset_cache()
>>> cm.get_cache()
({}, {})
cm = sage.structure.element.get_coercion_model()
len(cm.get_cache()[0])    # random
cm.reset_cache()
cm.get_cache()
richcmp(x, y, op)[source]

Given two arbitrary objects x and y, coerce them to a common parent and compare them using rich comparison operator op.

EXAMPLES:

sage: from sage.structure.element import get_coercion_model
sage: from sage.structure.richcmp import op_LT, op_LE, op_EQ, op_NE, op_GT, op_GE
sage: richcmp = get_coercion_model().richcmp
sage: richcmp(None, None, op_EQ)
True
sage: richcmp(None, 1, op_LT)
True
sage: richcmp("hello", None, op_LE)
False
sage: richcmp(-1, 1, op_GE)
False
sage: richcmp(int(1), float(2), op_GE)
False
>>> from sage.all import *
>>> from sage.structure.element import get_coercion_model
>>> from sage.structure.richcmp import op_LT, op_LE, op_EQ, op_NE, op_GT, op_GE
>>> richcmp = get_coercion_model().richcmp
>>> richcmp(None, None, op_EQ)
True
>>> richcmp(None, Integer(1), op_LT)
True
>>> richcmp("hello", None, op_LE)
False
>>> richcmp(-Integer(1), Integer(1), op_GE)
False
>>> richcmp(int(Integer(1)), float(Integer(2)), op_GE)
False
from sage.structure.element import get_coercion_model
from sage.structure.richcmp import op_LT, op_LE, op_EQ, op_NE, op_GT, op_GE
richcmp = get_coercion_model().richcmp
richcmp(None, None, op_EQ)
richcmp(None, 1, op_LT)
richcmp("hello", None, op_LE)
richcmp(-1, 1, op_GE)
richcmp(int(1), float(2), op_GE)

If there is no coercion, we only support == and !=:

sage: x = QQ.one(); y = GF(2).one()
sage: richcmp(x, y, op_EQ)
False
sage: richcmp(x, y, op_NE)
True
sage: richcmp(x, y, op_GT)
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for >:
'Rational Field' and 'Finite Field of size 2'
>>> from sage.all import *
>>> x = QQ.one(); y = GF(Integer(2)).one()
>>> richcmp(x, y, op_EQ)
False
>>> richcmp(x, y, op_NE)
True
>>> richcmp(x, y, op_GT)
Traceback (most recent call last):
...
TypeError: unsupported operand parent(s) for >:
'Rational Field' and 'Finite Field of size 2'
x = QQ.one(); y = GF(2).one()
richcmp(x, y, op_EQ)
richcmp(x, y, op_NE)
richcmp(x, y, op_GT)

We support non-Sage types with the usual Python convention:

sage: class AlwaysEqual():
....:     def __eq__(self, other):
....:         return True
sage: x = AlwaysEqual()
sage: x == 1
True
sage: 1 == x
True
>>> from sage.all import *
>>> class AlwaysEqual():
...     def __eq__(self, other):
...         return True
>>> x = AlwaysEqual()
>>> x == Integer(1)
True
>>> Integer(1) == x
True
class AlwaysEqual():
    def __eq__(self, other):
        return True
x = AlwaysEqual()
x == 1
1 == x
verify_action(action, R, S, op, fix=True)[source]

Verify that action takes an element of R on the left and S on the right, raising an error if not.

This is used for consistency checking in the coercion model.

EXAMPLES:

sage: R.<x> = ZZ['x']
sage: cm = sage.structure.element.get_coercion_model()
sage: cm.verify_action(R.get_action(QQ), R, QQ, operator.mul)
Right scalar multiplication by Rational Field
 on Univariate Polynomial Ring in x over Integer Ring
sage: cm.verify_action(R.get_action(QQ), RDF, R, operator.mul)
Traceback (most recent call last):
...
RuntimeError: There is a BUG in the coercion model:
    Action found for R <built-in function mul> S does not have the correct domains
    R = Real Double Field
    S = Univariate Polynomial Ring in x over Integer Ring
    (should be Univariate Polynomial Ring in x over Integer Ring, Rational Field)
    action = Right scalar multiplication by Rational Field
             on Univariate Polynomial Ring in x over Integer Ring
    (<class 'sage.structure.coerce_actions.RightModuleAction'>)
>>> from sage.all import *
>>> R = ZZ['x']; (x,) = R._first_ngens(1)
>>> cm = sage.structure.element.get_coercion_model()
>>> cm.verify_action(R.get_action(QQ), R, QQ, operator.mul)
Right scalar multiplication by Rational Field
 on Univariate Polynomial Ring in x over Integer Ring
>>> cm.verify_action(R.get_action(QQ), RDF, R, operator.mul)
Traceback (most recent call last):
...
RuntimeError: There is a BUG in the coercion model:
    Action found for R <built-in function mul> S does not have the correct domains
    R = Real Double Field
    S = Univariate Polynomial Ring in x over Integer Ring
    (should be Univariate Polynomial Ring in x over Integer Ring, Rational Field)
    action = Right scalar multiplication by Rational Field
             on Univariate Polynomial Ring in x over Integer Ring
    (<class 'sage.structure.coerce_actions.RightModuleAction'>)
R.<x> = ZZ['x']
cm = sage.structure.element.get_coercion_model()
cm.verify_action(R.get_action(QQ), R, QQ, operator.mul)
cm.verify_action(R.get_action(QQ), RDF, R, operator.mul)
verify_coercion_maps(R, S, homs, fix=False)[source]

Make sure this is a valid pair of homomorphisms from \(R\) and \(S\) to a common parent. This function is used to protect the user against buggy parents.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: homs = QQ.coerce_map_from(ZZ), None
sage: cm.verify_coercion_maps(ZZ, QQ, homs) == homs
True
sage: homs = QQ.coerce_map_from(ZZ), RR.coerce_map_from(QQ)
sage: cm.verify_coercion_maps(ZZ, QQ, homs) == homs
Traceback (most recent call last):
...
RuntimeError: ('BUG in coercion model, codomains must be identical',
Natural morphism:
  From: Integer Ring
  To:   Rational Field,
Generic map:
  From: Rational Field
  To:   Real Field with 53 bits of precision)
>>> from sage.all import *
>>> cm = sage.structure.element.get_coercion_model()
>>> homs = QQ.coerce_map_from(ZZ), None
>>> cm.verify_coercion_maps(ZZ, QQ, homs) == homs
True
>>> homs = QQ.coerce_map_from(ZZ), RR.coerce_map_from(QQ)
>>> cm.verify_coercion_maps(ZZ, QQ, homs) == homs
Traceback (most recent call last):
...
RuntimeError: ('BUG in coercion model, codomains must be identical',
Natural morphism:
  From: Integer Ring
  To:   Rational Field,
Generic map:
  From: Rational Field
  To:   Real Field with 53 bits of precision)
cm = sage.structure.element.get_coercion_model()
homs = QQ.coerce_map_from(ZZ), None
cm.verify_coercion_maps(ZZ, QQ, homs) == homs
homs = QQ.coerce_map_from(ZZ), RR.coerce_map_from(QQ)
cm.verify_coercion_maps(ZZ, QQ, homs) == homs
sage.structure.coerce.is_mpmath_type(t)[source]

Check whether the type t is a type whose name starts with either mpmath. or sage.libs.mpmath..

EXAMPLES:

sage: # needs mpmath
sage: from sage.structure.coerce import is_mpmath_type
sage: is_mpmath_type(int)
False
sage: import mpmath
sage: is_mpmath_type(mpmath.mpc(2))
False
sage: is_mpmath_type(type(mpmath.mpc(2)))
True
sage: is_mpmath_type(type(mpmath.mpf(2)))
True
>>> from sage.all import *
>>> # needs mpmath
>>> from sage.structure.coerce import is_mpmath_type
>>> is_mpmath_type(int)
False
>>> import mpmath
>>> is_mpmath_type(mpmath.mpc(Integer(2)))
False
>>> is_mpmath_type(type(mpmath.mpc(Integer(2))))
True
>>> is_mpmath_type(type(mpmath.mpf(Integer(2))))
True
# needs mpmath
from sage.structure.coerce import is_mpmath_type
is_mpmath_type(int)
import mpmath
is_mpmath_type(mpmath.mpc(2))
is_mpmath_type(type(mpmath.mpc(2)))
is_mpmath_type(type(mpmath.mpf(2)))
sage.structure.coerce.is_numpy_type(t)[source]

Return True if and only if \(t\) is a type whose name starts with numpy.

EXAMPLES:

sage: from sage.structure.coerce import is_numpy_type

sage: # needs numpy
sage: import numpy
sage: is_numpy_type(numpy.int16)
True
sage: is_numpy_type(numpy.floating)
True
sage: is_numpy_type(numpy.ndarray)
True
sage: is_numpy_type(numpy.matrix)
True

sage: is_numpy_type(int)
False
sage: is_numpy_type(Integer)
False
sage: is_numpy_type(Sudoku)                                                     # needs sage.combinat
False
sage: is_numpy_type(None)
False
>>> from sage.all import *
>>> from sage.structure.coerce import is_numpy_type

>>> # needs numpy
>>> import numpy
>>> is_numpy_type(numpy.int16)
True
>>> is_numpy_type(numpy.floating)
True
>>> is_numpy_type(numpy.ndarray)
True
>>> is_numpy_type(numpy.matrix)
True

>>> is_numpy_type(int)
False
>>> is_numpy_type(Integer)
False
>>> is_numpy_type(Sudoku)                                                     # needs sage.combinat
False
>>> is_numpy_type(None)
False
from sage.structure.coerce import is_numpy_type
# needs numpy
import numpy
is_numpy_type(numpy.int16)
is_numpy_type(numpy.floating)
is_numpy_type(numpy.ndarray)
is_numpy_type(numpy.matrix)
is_numpy_type(int)
is_numpy_type(Integer)
is_numpy_type(Sudoku)                                                     # needs sage.combinat
is_numpy_type(None)
sage.structure.coerce.parent_is_integers(P)[source]

Check whether the type or parent represents the ring of integers.

EXAMPLES:

sage: from sage.structure.coerce import parent_is_integers
sage: parent_is_integers(int)
True
sage: parent_is_integers(float)
False
sage: parent_is_integers(bool)
True
sage: parent_is_integers(dict)
False

sage: import numpy                                                              # needs numpy
sage: parent_is_integers(numpy.int16)                                           # needs numpy
True
sage: parent_is_integers(numpy.uint64)                                          # needs numpy
True
sage: parent_is_integers(float)
False

sage: import gmpy2
sage: parent_is_integers(gmpy2.mpz)
True
sage: parent_is_integers(gmpy2.mpq)
False
>>> from sage.all import *
>>> from sage.structure.coerce import parent_is_integers
>>> parent_is_integers(int)
True
>>> parent_is_integers(float)
False
>>> parent_is_integers(bool)
True
>>> parent_is_integers(dict)
False

>>> import numpy                                                              # needs numpy
>>> parent_is_integers(numpy.int16)                                           # needs numpy
True
>>> parent_is_integers(numpy.uint64)                                          # needs numpy
True
>>> parent_is_integers(float)
False

>>> import gmpy2
>>> parent_is_integers(gmpy2.mpz)
True
>>> parent_is_integers(gmpy2.mpq)
False
from sage.structure.coerce import parent_is_integers
parent_is_integers(int)
parent_is_integers(float)
parent_is_integers(bool)
parent_is_integers(dict)
import numpy                                                              # needs numpy
parent_is_integers(numpy.int16)                                           # needs numpy
parent_is_integers(numpy.uint64)                                          # needs numpy
parent_is_integers(float)
import gmpy2
parent_is_integers(gmpy2.mpz)
parent_is_integers(gmpy2.mpq)

Ensure (Issue #27893) is fixed:

sage: K.<f> = QQ[]
sage: gmpy2.mpz(2) * f
2*f
>>> from sage.all import *
>>> K = QQ['f']; (f,) = K._first_ngens(1)
>>> gmpy2.mpz(Integer(2)) * f
2*f
K.<f> = QQ[]
gmpy2.mpz(2) * f
sage.structure.coerce.parent_is_numerical(P)[source]

Test if elements of the parent or type P can be numerically evaluated as complex numbers (in a canonical way).

EXAMPLES:

sage: from sage.structure.coerce import parent_is_numerical
sage: import gmpy2
sage: [parent_is_numerical(R) for R in [QQ, int, complex, gmpy2.mpc]]           # needs sage.rings.complex_double
[True, True, True, True]
sage: [parent_is_numerical(R) for R in [RR, CC]]                                # needs sage.rings.real_mpfr
[True, True]
sage: parent_is_numerical(QuadraticField(-1))                                   # needs sage.rings.number_field
True
sage: import numpy; parent_is_numerical(numpy.complexfloating)                  # needs numpy
True
sage: parent_is_numerical(SR)                                                   # needs sage.symbolic
False
sage: [parent_is_numerical(R) for R in [QQ['x'], QQ[['x']], str]]
[False, False, False]
sage: [parent_is_numerical(R) for R in [RIF, RBF, CIF, CBF]]                    # needs sage.libs.flint
[False, False, False, False]
>>> from sage.all import *
>>> from sage.structure.coerce import parent_is_numerical
>>> import gmpy2
>>> [parent_is_numerical(R) for R in [QQ, int, complex, gmpy2.mpc]]           # needs sage.rings.complex_double
[True, True, True, True]
>>> [parent_is_numerical(R) for R in [RR, CC]]                                # needs sage.rings.real_mpfr
[True, True]
>>> parent_is_numerical(QuadraticField(-Integer(1)))                                   # needs sage.rings.number_field
True
>>> import numpy; parent_is_numerical(numpy.complexfloating)                  # needs numpy
True
>>> parent_is_numerical(SR)                                                   # needs sage.symbolic
False
>>> [parent_is_numerical(R) for R in [QQ['x'], QQ[['x']], str]]
[False, False, False]
>>> [parent_is_numerical(R) for R in [RIF, RBF, CIF, CBF]]                    # needs sage.libs.flint
[False, False, False, False]
from sage.structure.coerce import parent_is_numerical
import gmpy2
[parent_is_numerical(R) for R in [QQ, int, complex, gmpy2.mpc]]           # needs sage.rings.complex_double
[parent_is_numerical(R) for R in [RR, CC]]                                # needs sage.rings.real_mpfr
parent_is_numerical(QuadraticField(-1))                                   # needs sage.rings.number_field
import numpy; parent_is_numerical(numpy.complexfloating)                  # needs numpy
parent_is_numerical(SR)                                                   # needs sage.symbolic
[parent_is_numerical(R) for R in [QQ['x'], QQ[['x']], str]]
[parent_is_numerical(R) for R in [RIF, RBF, CIF, CBF]]                    # needs sage.libs.flint
sage.structure.coerce.parent_is_real_numerical(P)[source]

Test if elements of the parent or type P can be numerically evaluated as real numbers (in a canonical way).

EXAMPLES:

sage: from sage.structure.coerce import parent_is_real_numerical
sage: import gmpy2
sage: [parent_is_real_numerical(R) for R in [RR, QQ, ZZ, RLF, int, float, gmpy2.mpq]]
[True, True, True, True, True, True, True]
sage: parent_is_real_numerical(QuadraticField(2))                               # needs sage.rings.number_field
True
sage: import numpy; parent_is_real_numerical(numpy.integer)                     # needs numpy
True
sage: parent_is_real_numerical(QuadraticField(-1))                              # needs sage.rings.number_field
False
sage: [parent_is_real_numerical(R)                                              # needs numpy
....:  for R in [CC, complex, gmpy2.mpc, numpy.complexfloating]]
[False, False, False, False]
sage: [parent_is_real_numerical(R) for R in [QQ['x'], QQ[['x']], str]]
[False, False, False]
sage: parent_is_real_numerical(SR)                                              # needs sage.symbolic
False
sage: [parent_is_real_numerical(R) for R in [RIF, RBF, CIF, CBF]]               # needs sage.libs.flint
[False, False, False, False]
>>> from sage.all import *
>>> from sage.structure.coerce import parent_is_real_numerical
>>> import gmpy2
>>> [parent_is_real_numerical(R) for R in [RR, QQ, ZZ, RLF, int, float, gmpy2.mpq]]
[True, True, True, True, True, True, True]
>>> parent_is_real_numerical(QuadraticField(Integer(2)))                               # needs sage.rings.number_field
True
>>> import numpy; parent_is_real_numerical(numpy.integer)                     # needs numpy
True
>>> parent_is_real_numerical(QuadraticField(-Integer(1)))                              # needs sage.rings.number_field
False
>>> [parent_is_real_numerical(R)                                              # needs numpy
...  for R in [CC, complex, gmpy2.mpc, numpy.complexfloating]]
[False, False, False, False]
>>> [parent_is_real_numerical(R) for R in [QQ['x'], QQ[['x']], str]]
[False, False, False]
>>> parent_is_real_numerical(SR)                                              # needs sage.symbolic
False
>>> [parent_is_real_numerical(R) for R in [RIF, RBF, CIF, CBF]]               # needs sage.libs.flint
[False, False, False, False]
from sage.structure.coerce import parent_is_real_numerical
import gmpy2
[parent_is_real_numerical(R) for R in [RR, QQ, ZZ, RLF, int, float, gmpy2.mpq]]
parent_is_real_numerical(QuadraticField(2))                               # needs sage.rings.number_field
import numpy; parent_is_real_numerical(numpy.integer)                     # needs numpy
parent_is_real_numerical(QuadraticField(-1))                              # needs sage.rings.number_field
[parent_is_real_numerical(R)                                              # needs numpy
 for R in [CC, complex, gmpy2.mpc, numpy.complexfloating]]
[parent_is_real_numerical(R) for R in [QQ['x'], QQ[['x']], str]]
parent_is_real_numerical(SR)                                              # needs sage.symbolic
[parent_is_real_numerical(R) for R in [RIF, RBF, CIF, CBF]]               # needs sage.libs.flint
sage.structure.coerce.py_scalar_parent(py_type)[source]

Return the Sage equivalent of the given python type, if one exists. If there is no equivalent, return None.

EXAMPLES:

sage: from sage.structure.coerce import py_scalar_parent
sage: py_scalar_parent(int)
Integer Ring
sage: py_scalar_parent(float)
Real Double Field
sage: py_scalar_parent(complex)                                                 # needs sage.rings.complex_double
Complex Double Field
sage: py_scalar_parent(bool)
Integer Ring
sage: py_scalar_parent(dict),
(None,)

sage: import fractions
sage: py_scalar_parent(fractions.Fraction)
Rational Field

sage: # needs numpy
sage: import numpy
sage: py_scalar_parent(numpy.int16)
Integer Ring
sage: py_scalar_parent(numpy.int32)
Integer Ring
sage: py_scalar_parent(numpy.uint64)
Integer Ring
sage: py_scalar_parent(numpy.double)
Real Double Field

sage: import gmpy2
sage: py_scalar_parent(gmpy2.mpz)
Integer Ring
sage: py_scalar_parent(gmpy2.mpq)
Rational Field
sage: py_scalar_parent(gmpy2.mpfr)
Real Double Field
sage: py_scalar_parent(gmpy2.mpc)                                               # needs sage.rings.complex_double
Complex Double Field
>>> from sage.all import *
>>> from sage.structure.coerce import py_scalar_parent
>>> py_scalar_parent(int)
Integer Ring
>>> py_scalar_parent(float)
Real Double Field
>>> py_scalar_parent(complex)                                                 # needs sage.rings.complex_double
Complex Double Field
>>> py_scalar_parent(bool)
Integer Ring
>>> py_scalar_parent(dict),
(None,)

>>> import fractions
>>> py_scalar_parent(fractions.Fraction)
Rational Field

>>> # needs numpy
>>> import numpy
>>> py_scalar_parent(numpy.int16)
Integer Ring
>>> py_scalar_parent(numpy.int32)
Integer Ring
>>> py_scalar_parent(numpy.uint64)
Integer Ring
>>> py_scalar_parent(numpy.double)
Real Double Field

>>> import gmpy2
>>> py_scalar_parent(gmpy2.mpz)
Integer Ring
>>> py_scalar_parent(gmpy2.mpq)
Rational Field
>>> py_scalar_parent(gmpy2.mpfr)
Real Double Field
>>> py_scalar_parent(gmpy2.mpc)                                               # needs sage.rings.complex_double
Complex Double Field
from sage.structure.coerce import py_scalar_parent
py_scalar_parent(int)
py_scalar_parent(float)
py_scalar_parent(complex)                                                 # needs sage.rings.complex_double
py_scalar_parent(bool)
py_scalar_parent(dict),
import fractions
py_scalar_parent(fractions.Fraction)
# needs numpy
import numpy
py_scalar_parent(numpy.int16)
py_scalar_parent(numpy.int32)
py_scalar_parent(numpy.uint64)
py_scalar_parent(numpy.double)
import gmpy2
py_scalar_parent(gmpy2.mpz)
py_scalar_parent(gmpy2.mpq)
py_scalar_parent(gmpy2.mpfr)
py_scalar_parent(gmpy2.mpc)                                               # needs sage.rings.complex_double
sage.structure.coerce.py_scalar_to_element(x)[source]

Convert x to a Sage Element if possible.

If x was already an Element or if there is no obvious conversion possible, just return x itself.

EXAMPLES:

sage: from sage.structure.coerce import py_scalar_to_element
sage: x = py_scalar_to_element(42)
sage: x, parent(x)
(42, Integer Ring)
sage: x = py_scalar_to_element(int(42))
sage: x, parent(x)
(42, Integer Ring)
sage: x = py_scalar_to_element(float(42))
sage: x, parent(x)
(42.0, Real Double Field)
sage: x = py_scalar_to_element(complex(42))                                     # needs sage.rings.complex_double
sage: x, parent(x)                                                              # needs sage.rings.complex_double
(42.0, Complex Double Field)
sage: py_scalar_to_element('hello')
'hello'

sage: from fractions import Fraction
sage: f = Fraction(int(2^100), int(3^100))
sage: py_scalar_to_element(f)
1267650600228229401496703205376/515377520732011331036461129765621272702107522001
>>> from sage.all import *
>>> from sage.structure.coerce import py_scalar_to_element
>>> x = py_scalar_to_element(Integer(42))
>>> x, parent(x)
(42, Integer Ring)
>>> x = py_scalar_to_element(int(Integer(42)))
>>> x, parent(x)
(42, Integer Ring)
>>> x = py_scalar_to_element(float(Integer(42)))
>>> x, parent(x)
(42.0, Real Double Field)
>>> x = py_scalar_to_element(complex(Integer(42)))                                     # needs sage.rings.complex_double
>>> x, parent(x)                                                              # needs sage.rings.complex_double
(42.0, Complex Double Field)
>>> py_scalar_to_element('hello')
'hello'

>>> from fractions import Fraction
>>> f = Fraction(int(Integer(2)**Integer(100)), int(Integer(3)**Integer(100)))
>>> py_scalar_to_element(f)
1267650600228229401496703205376/515377520732011331036461129765621272702107522001
from sage.structure.coerce import py_scalar_to_element
x = py_scalar_to_element(42)
x, parent(x)
x = py_scalar_to_element(int(42))
x, parent(x)
x = py_scalar_to_element(float(42))
x, parent(x)
x = py_scalar_to_element(complex(42))                                     # needs sage.rings.complex_double
x, parent(x)                                                              # needs sage.rings.complex_double
py_scalar_to_element('hello')
from fractions import Fraction
f = Fraction(int(2^100), int(3^100))
py_scalar_to_element(f)

Note that bools are converted to 0 or 1:

sage: py_scalar_to_element(False), py_scalar_to_element(True)
(0, 1)
>>> from sage.all import *
>>> py_scalar_to_element(False), py_scalar_to_element(True)
(0, 1)
py_scalar_to_element(False), py_scalar_to_element(True)

Test gmpy2’s types:

sage: import gmpy2
sage: x = py_scalar_to_element(gmpy2.mpz(42))
sage: x, parent(x)
(42, Integer Ring)
sage: x = py_scalar_to_element(gmpy2.mpq('3/4'))
sage: x, parent(x)
(3/4, Rational Field)
sage: x = py_scalar_to_element(gmpy2.mpfr(42.57))
sage: x, parent(x)
(42.57, Real Double Field)
sage: x = py_scalar_to_element(gmpy2.mpc(int(42), int(42)))                     # needs sage.rings.complex_double
sage: x, parent(x)                                                              # needs sage.rings.complex_double
(42.0 + 42.0*I, Complex Double Field)
>>> from sage.all import *
>>> import gmpy2
>>> x = py_scalar_to_element(gmpy2.mpz(Integer(42)))
>>> x, parent(x)
(42, Integer Ring)
>>> x = py_scalar_to_element(gmpy2.mpq('3/4'))
>>> x, parent(x)
(3/4, Rational Field)
>>> x = py_scalar_to_element(gmpy2.mpfr(RealNumber('42.57')))
>>> x, parent(x)
(42.57, Real Double Field)
>>> x = py_scalar_to_element(gmpy2.mpc(int(Integer(42)), int(Integer(42))))                     # needs sage.rings.complex_double
>>> x, parent(x)                                                              # needs sage.rings.complex_double
(42.0 + 42.0*I, Complex Double Field)
import gmpy2
x = py_scalar_to_element(gmpy2.mpz(42))
x, parent(x)
x = py_scalar_to_element(gmpy2.mpq('3/4'))
x, parent(x)
x = py_scalar_to_element(gmpy2.mpfr(42.57))
x, parent(x)
x = py_scalar_to_element(gmpy2.mpc(int(42), int(42)))                     # needs sage.rings.complex_double
x, parent(x)                                                              # needs sage.rings.complex_double

Test compatibility with py_scalar_parent():

sage: from sage.structure.coerce import py_scalar_parent
sage: elt = [True, int(42), float(42), complex(42)]
sage: for x in elt:                                                             # needs sage.rings.complex_double
....:     assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()

sage: import numpy                                                              # needs numpy
sage: elt = [numpy.int8('-12'),  numpy.uint8('143'),                            # needs numpy
....:        numpy.int16('-33'), numpy.uint16('122'),
....:        numpy.int32('-19'), numpy.uint32('44'),
....:        numpy.int64('-3'),  numpy.uint64('552'),
....:        numpy.float16('-1.23'), numpy.float32('-2.22'),
....:        numpy.float64('-3.412'), numpy.complex64(1.2+I),
....:        numpy.complex128(-2+I)]
sage: for x in elt:                                                             # needs numpy
....:     assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()

sage: elt = [gmpy2.mpz(42), gmpy2.mpq('3/4'),
....:        gmpy2.mpfr(42.57), gmpy2.mpc(int(42), int(42))]
sage: for x in elt:                                                             # needs sage.rings.complex_double
....:     assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()
>>> from sage.all import *
>>> from sage.structure.coerce import py_scalar_parent
>>> elt = [True, int(Integer(42)), float(Integer(42)), complex(Integer(42))]
>>> for x in elt:                                                             # needs sage.rings.complex_double
...     assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()

>>> import numpy                                                              # needs numpy
>>> elt = [numpy.int8('-12'),  numpy.uint8('143'),                            # needs numpy
...        numpy.int16('-33'), numpy.uint16('122'),
...        numpy.int32('-19'), numpy.uint32('44'),
...        numpy.int64('-3'),  numpy.uint64('552'),
...        numpy.float16('-1.23'), numpy.float32('-2.22'),
...        numpy.float64('-3.412'), numpy.complex64(RealNumber('1.2')+I),
...        numpy.complex128(-Integer(2)+I)]
>>> for x in elt:                                                             # needs numpy
...     assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()

>>> elt = [gmpy2.mpz(Integer(42)), gmpy2.mpq('3/4'),
...        gmpy2.mpfr(RealNumber('42.57')), gmpy2.mpc(int(Integer(42)), int(Integer(42)))]
>>> for x in elt:                                                             # needs sage.rings.complex_double
...     assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()
from sage.structure.coerce import py_scalar_parent
elt = [True, int(42), float(42), complex(42)]
for x in elt:                                                             # needs sage.rings.complex_double
    assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()
import numpy                                                              # needs numpy
elt = [numpy.int8('-12'),  numpy.uint8('143'),                            # needs numpy
       numpy.int16('-33'), numpy.uint16('122'),
       numpy.int32('-19'), numpy.uint32('44'),
       numpy.int64('-3'),  numpy.uint64('552'),
       numpy.float16('-1.23'), numpy.float32('-2.22'),
       numpy.float64('-3.412'), numpy.complex64(1.2+I),
       numpy.complex128(-2+I)]
for x in elt:                                                             # needs numpy
    assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()
elt = [gmpy2.mpz(42), gmpy2.mpq('3/4'),
       gmpy2.mpfr(42.57), gmpy2.mpc(int(42), int(42))]
for x in elt:                                                             # needs sage.rings.complex_double
    assert py_scalar_parent(type(x)) == py_scalar_to_element(x).parent()