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

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

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')
sage: P = E(0,0)
sage: 5*P
(1/4 : -5/8 : 1)

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_cache_maps

Bases: sage.structure.element.CoercionModel

See also sage.categories.pushout

EXAMPLES:

sage: f = ZZ['t','x'].0 + QQ['x'].0 + CyclotomicField(13).gen(); f
t + x + (zeta13)
sage: f.parent()
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
[(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

AUTHOR:

  • Robert Bradshaw
analyse(xp, yp, op='mul')

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 function 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: steps, res = cm.analyse(GF(7), ZZ)
sage: print steps
['Coercion on right operand via', Natural morphism:
  From: Integer Ring
  To:   Finite Field of size 7, 'Arithmetic performed after coercions.']
sage: print res
Finite Field of size 7
sage: f = steps[1]; type(f)
<type 'sage.rings.finite_rings.integer_mod.Integer_to_IntegerMod'>
sage: f(100)
2
bin_op(x, y, op)

Execute the operation op on x and y. It first looks for an action corresponding to op, and failing that, it tries to coerces x and y into the 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

The operator can be any callable:

set Rational Field Integer Ring <function <lambda> at 0xc0b2270> None None
(Rational Field, Rational Field)
sage: R.<x> = ZZ['x']
sage: cm.bin_op(x^2-1, x+1, gcd)
x + 1

Actions are detected and performed:

sage: M = matrix(ZZ, 2, 2, range(4))
sage: V = vector(ZZ, [5,7])
sage: cm.bin_op(M, V, operator.mul)
(7, 31)

TESTS:

sage: class Foo:                  
...      def __rmul__(self, left):
...          return 'hello'
...
sage: H = Foo()
sage: print int(3)*H
hello
sage: print Integer(3)*H 
hello
sage: print H*3
...
TypeError: unsupported operand parent(s) for '*': '<type 'instance'>' and 'Integer Ring'

sage: class Nonsense:
...       def __init__(self, s):
...           self.s = s
...       def __repr__(self):
...           return self.s
...       def __mul__(self, x):
...           return Nonsense(self.s + chr(x%256))
...       __add__ = __mul__
...       def __rmul__(self, x):
...           return Nonsense(chr(x%256) + self.s)
...       __radd__ = __rmul__
...
sage: a = Nonsense('blahblah')
sage: a*80
blahblahP
sage: 80*a
Pblahblah
sage: a+80
blahblahP
sage: 80+a
Pblahblah
canonical_coercion(x, y)

Given two elements x and y, with parents S and R respectively, find a common parent Z such that there are coercions f: S \mapsto Z and g: R \mapsto 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: 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

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

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


sage: x, y = cm.canonical_coercion(int(5), complex(3))
sage: type(x), type(y)
(<type 'complex'>, <type '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)
<type 'sage.rings.rational.Rational'>

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

sage: canonical_coercion(vector([1, 2, 3]), 0)
((1, 2, 3), (0, 0, 0))
coercion_maps(R, S)

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.

EXAMPLES:

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

sage: f, g = cm.coercion_maps(ZZ['x'], QQ)
sage: print f
Conversion map:
  From: Univariate Polynomial Ring in x over Integer Ring
  To:   Univariate Polynomial Ring in x over Rational Field
sage: print g
Polynomial base injection morphism:
  From: Rational Field
  To:   Univariate Polynomial Ring in x over Rational Field
  
sage: cm.coercion_maps(QQ, GF(7)) == None
True

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

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,
 Call morphism:
  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,
 Call morphism:
  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
discover_action(R, S, op)

INPUT

  • R - the left Parent (or type)
  • S - the right Parent (or type)
  • op - the operand, typically an element of the operator module.

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

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

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 type 'int'
  To:   Integer Ring

If op in an inplace operation, look for the non-inplace action:

sage: cm.discover_action(P, ZZ, operator.imul)
Right scalar multiplication by Integer Ring on Univariate Polynomial Ring in x over Integer Ring

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

sage: cm.discover_action(P, ZZ, operator.div)
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
discover_coercion(R, S)

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

EXAMPLES:

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

If R is S, then two identity morphisms suffice:

sage: cm.discover_coercion(SR, SR)
(None, None)

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

sage: cm.discover_coercion(ZZ, QQ)
(Natural morphism:
  From: Integer Ring
  To:   Rational Field, None)
sage: cm.discover_coercion(RR, QQ)
(None,
 Generic map:
  From: Rational Field
  To:   Real Field with 53 bits of precision)

Otherwise, try and compute an appropriate cover:

sage: cm.discover_coercion(ZZ['x,y'], RDF)
(Call morphism:
  From: Multivariate Polynomial Ring in x, y over Integer Ring
  To:   Multivariate Polynomial Ring in x, y over Real Double Field,
 Call morphism:
  From: Real Double Field
  To:   Multivariate Polynomial Ring in x, y over Real Double Field)

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

sage: sage.categories.pushout.pushout(QQ, QQ^3)
Vector space of dimension 3 over Rational Field
sage: print cm.discover_coercion(QQ, QQ^3)
None
division_parent(parent)

Deduces where the result of division in parent lies by calculating the inverse of parent.one_element() or parent.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: cm.division_parent(ZZ['x'])
Fraction Field of Univariate Polynomial Ring in x over Integer Ring
sage: cm.division_parent(GF(41))
Finite Field of size 41
sage: cm.division_parent(Integers(100))
Ring of integers modulo 100
sage: cm.division_parent(SymmetricGroup(5))
Symmetric group of order 5! as a permutation group
exception_stack()

Returns 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: 1/2 + 2
5/2
sage: cm.exception_stack()
[]
sage: 1/2 + GF(3)(2)
...
TypeError: unsupported operand parent(s) for '+': 'Rational Field' and 'Finite Field of size 3'

Now see what the actual problem was:

sage: import traceback
sage: cm.exception_stack()
[(<class 'sage.structure.coerce_exceptions.CoercionException'>, CoercionException("BUG: the base_extend method must be defined for 'Monoid of ideals of Integer Ring' (class '<class 'sage.rings.ideal_monoid.IdealMonoid_c'>')",), <traceback object at ...>), (<type 'exceptions.TypeError'>,  TypeError("no common canonical parent for objects with parents: 'Rational Field' and 'Finite Field of size 3'",),  <traceback object at ...>)]
sage: print ''.join(sum([traceback.format_exception(*info) for info in cm.exception_stack()], []))
...
TypeError: no common canonical parent for objects with parents: 'Rational Field' and 'Finite Field of size 3'

This is typically accessed via the coercion_traceback() function.

sage: coercion_traceback()
...
TypeError: no common canonical parent for objects with parents: 'Rational Field' and 'Finite Field of size 3'
explain(xp, yp, op='operator.mul', verbosity=2)

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 type 'int'
      To:   Rational Field
Arithmetic performed after coercions.
Result lives in Rational Field
Rational Field

sage: cm.explain(ZZ['x'], 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
    Conversion map:
      From: Univariate Polynomial Ring in x over Integer Ring
      To:   Univariate Polynomial Ring in x over 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

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

sage: cm.explain(RealField(100), float, operator.add)
Right operand is numeric, will attempt coercion in both directions.
Unknown result parent.
sage: parent(RealField(100)(1) + float(1))
<type '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))
<type 'float'>

Special care is taken to deal with division:

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

sage: cm.explain(ZZ['x'], QQ['x'], operator.div)
Coercion on left operand via
    Conversion map:
      From: Univariate Polynomial Ring in x over Integer Ring
      To:   Univariate Polynomial Ring in x over 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.div)
Coercion on left operand via
    Native morphism:
      From: Set of Python objects of type 'int'
      To:   Integer Ring
Arithmetic performed after coercions.
Result lives in Rational Field
Rational Field

sage: cm.explain(ZZ['x'], ZZ, operator.div)
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

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)

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: cm.get_action(ZZ['x'], ZZ, operator.mul)
Right scalar multiplication by Integer Ring on Univariate Polynomial Ring in x over Integer Ring
sage: cm.get_action(ZZ['x'], ZZ, operator.imul)
Right scalar multiplication by Integer Ring on Univariate Polynomial Ring in x over Integer Ring
sage: cm.get_action(ZZ['x'], QQ, operator.mul)
Right scalar multiplication by Rational Field on Univariate Polynomial Ring in x over Integer Ring
sage: cm.get_action(QQ['x'], 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 type 'int'
  To:   Integer Ring

sage: R.<x> = QQ['x']
sage: A = cm.get_action(R, ZZ, operator.div); 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: A(x+10, 5)
1/5*x + 2
get_cache()

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

EXAMPLES:

sage: 1 + 1/2
3/2
sage: cm = sage.structure.element.get_coercion_model()
sage: maps, actions = cm.get_cache()

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

sage: left_morphism, right_morphism = maps[ZZ, QQ]
sage: print left_morphism
Natural morphism:
  From: Integer Ring
  To:   Rational Field
sage: print right_morphism
None

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
get_stats()

This returns the state of the cache of coercion maps and actions, primarily useful for debugging and introspection. If a class is not part of the coercion system, we should call the __rmul__ method when it makes sense.

The coercion maps are stored in a specialized TripleDict hashtable, and the stats returned are (min, avg, max) of the number of items per bucket. The lower the better.

EXAMPLES:

sage: cm = sage.structure.element.get_coercion_model()
sage: cm.get_stats()                # random
((0, 0.16058394160583941, 2), (0, 0.13138686131386862, 3))
reset_cache(lookup_dict_size=127, lookup_dict_threshold=0.75)

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: cm.get_stats()                # random
((0, 0.3307086614173229, 3), (0, 0.1496062992125984, 2))
sage: cm.reset_cache()
sage: cm.get_stats()
((0, 0.0, 0), (0, 0.0, 0))
verify_action(action, R, S, op, fix=True)

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)
...
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 (<type 'sage.structure.coerce_actions.RightModuleAction'>)
verify_coercion_maps(R, S, homs, fix=False)

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
...
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)
sage.structure.coerce.parent(x)
sage.structure.coerce.py_scalar_parent(py_type)

Returns 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(long)
Integer Ring
sage: py_scalar_parent(float)
Real Double Field
sage: py_scalar_parent(complex)
Complex Double Field
sage: py_scalar_parent(dict),
(None,)

Previous topic

Coercion

Next topic

Coerce actions

This Page