Tropical Varieties

A tropical variety is a piecewise-linear geometric object derived from a classical algebraic variety by using tropical mathematics, where the tropical semiring replaces the usual arithmetic operations.

AUTHORS:

  • Verrel Rievaldo Wijaya (2024-06): initial version

REFERENCES:

class sage.rings.semirings.tropical_variety.TropicalCurve(poly)[source]

Bases: TropicalVariety

A tropical curve in \(\RR^2\).

The tropical curve consists of line segments and half-lines, which we call edges. These edges are connected in such a way that they form a piecewise linear graph embedded in the plane. These edges meet at a vertices, where the balancing condition is satisfied. This balancing condition ensures that the sum of the outgoing slopes at each vertex is zero, reflecting the equilibrium.

EXAMPLES:

sage: T = TropicalSemiring(QQ, use_min=False)
sage: R.<x,y> = PolynomialRing(T)
sage: p1 = x + y + R(0)
sage: tv = p1.tropical_variety(); tv
Tropical curve of 0*x + 0*y + 0
sage: tv.components()
[[(t1, t1), [t1 >= 0], 1], [(0, t1), [t1 <= 0], 1], [(t1, 0), [t1 <= 0], 1]]
>>> from sage.all import *
>>> T = TropicalSemiring(QQ, use_min=False)
>>> R = PolynomialRing(T, names=('x', 'y',)); (x, y,) = R._first_ngens(2)
>>> p1 = x + y + R(Integer(0))
>>> tv = p1.tropical_variety(); tv
Tropical curve of 0*x + 0*y + 0
>>> tv.components()
[[(t1, t1), [t1 >= 0], 1], [(0, t1), [t1 <= 0], 1], [(t1, 0), [t1 <= 0], 1]]
plot()[source]

Return the plot of self.

Generates a visual representation of the tropical curve in cartesian coordinates. The plot shows piecewise-linear segments representing each components. The axes are centered around the vertices.

OUTPUT:

A Graphics object. The weight of the component will be written if it is greater or equal than 2. The weight is written near the vertex.

EXAMPLES:

A polynomial with only two terms will give one straight line:

sage: T = TropicalSemiring(QQ)
sage: R.<x,y> = PolynomialRing(T)
sage: (y+R(1)).tropical_variety().components()
[[(t1, 1), [-Infinity < t1, t1 < +Infinity], 1]]
sage: (y+R(1)).tropical_variety().plot()
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> T = TropicalSemiring(QQ)
>>> R = PolynomialRing(T, names=('x', 'y',)); (x, y,) = R._first_ngens(2)
>>> (y+R(Integer(1))).tropical_variety().components()
[[(t1, 1), [-Infinity < t1, t1 < +Infinity], 1]]
>>> (y+R(Integer(1))).tropical_variety().plot()
Graphics object consisting of 1 graphics primitive
../../../_images/tropical_variety-1.svg

An intriguing and fascinating tropical curve can be obtained with a more complex tropical polynomial:

sage: p1 = R(1) + R(2)*x + R(3)*y + R(6)*x*y + R(10)*x*y^2
sage: p1.tropical_variety().components()
[[(-1, t1), [-2 <= t1], 1],
[(t1, -2), [-1 <= t1], 1],
[(t1 + 1, t1), [-4 <= t1, t1 <= -2], 1],
[(t1, -4), [t1 <= -3], 2],
[(-t1 - 7, t1), [t1 <= -4], 1]]
sage: p1.tropical_variety().plot()
Graphics object consisting of 6 graphics primitives
>>> from sage.all import *
>>> p1 = R(Integer(1)) + R(Integer(2))*x + R(Integer(3))*y + R(Integer(6))*x*y + R(Integer(10))*x*y**Integer(2)
>>> p1.tropical_variety().components()
[[(-1, t1), [-2 <= t1], 1],
[(t1, -2), [-1 <= t1], 1],
[(t1 + 1, t1), [-4 <= t1, t1 <= -2], 1],
[(t1, -4), [t1 <= -3], 2],
[(-t1 - 7, t1), [t1 <= -4], 1]]
>>> p1.tropical_variety().plot()
Graphics object consisting of 6 graphics primitives
../../../_images/tropical_variety-2.svg

Another tropical polynomial with numerous components, resulting in a more intricate structure:

sage: p2 = (x^6 + R(4)*x^4*y^2 + R(2)*x^3*y^3 + R(3)*x^2*y^4 + x*y^5
....:       + R(7)*x^2 + R(5)*x*y + R(3)*y^2 + R(2)*x + y + R(10))
sage: p2.tropical_variety().plot()
Graphics object consisting of 11 graphics primitives
>>> from sage.all import *
>>> p2 = (x**Integer(6) + R(Integer(4))*x**Integer(4)*y**Integer(2) + R(Integer(2))*x**Integer(3)*y**Integer(3) + R(Integer(3))*x**Integer(2)*y**Integer(4) + x*y**Integer(5)
...       + R(Integer(7))*x**Integer(2) + R(Integer(5))*x*y + R(Integer(3))*y**Integer(2) + R(Integer(2))*x + y + R(Integer(10)))
>>> p2.tropical_variety().plot()
Graphics object consisting of 11 graphics primitives
../../../_images/tropical_variety-3.svg
vertices()[source]

Return all vertices of self, which is the point where three or more edges intersect.

OUTPUT: A set of \((x,y)\) points

EXAMPLES:

sage: T = TropicalSemiring(QQ)
sage: R.<x,y> = PolynomialRing(T)
sage: p1 = x + y
sage: p1.tropical_variety().vertices()
set()
sage: p2 = R(-2)*x^2 + R(-1)*x + R(1/2)*y + R(1/6)
sage: p2.tropical_variety().vertices()
{(1, -1/2), (7/6, -1/3)}
>>> from sage.all import *
>>> T = TropicalSemiring(QQ)
>>> R = PolynomialRing(T, names=('x', 'y',)); (x, y,) = R._first_ngens(2)
>>> p1 = x + y
>>> p1.tropical_variety().vertices()
set()
>>> p2 = R(-Integer(2))*x**Integer(2) + R(-Integer(1))*x + R(Integer(1)/Integer(2))*y + R(Integer(1)/Integer(6))
>>> p2.tropical_variety().vertices()
{(1, -1/2), (7/6, -1/3)}
class sage.rings.semirings.tropical_variety.TropicalSurface(poly)[source]

Bases: TropicalVariety

A tropical surface in \(\RR^3\).

The tropical surface consists of planar regions and facets, which we can call cells. These cells are connected in such a way that they form a piecewise linear structure embedded in three-dimensional space. These cells meet along edges, where the balancing condition is satisfied. This balancing condition ensures that the sum of the outgoing normal vectors at each edge is zero, reflecting the equilibrium.

EXAMPLES:

sage: T = TropicalSemiring(QQ, use_min=False)
sage: R.<x,y,z> = PolynomialRing(T)
sage: p1 = x + y + z + R(0)
sage: tv = p1.tropical_variety(); tv
Tropical surface of 0*x + 0*y + 0*z + 0
sage: tv.components()
[[(t1, t1, t2), [t2 <= t1, 0 <= t1], 1],
[(t1, t2, t1), [max(0, t2) <= t1], 1],
[(0, t1, t2), [t2 <= 0, t1 <= 0], 1],
[(t1, t2, t2), [max(0, t1) <= t2], 1],
[(t1, 0, t2), [t2 <= 0, t1 <= 0], 1],
[(t1, t2, 0), [t1 <= 0, t2 <= 0], 1]]
>>> from sage.all import *
>>> T = TropicalSemiring(QQ, use_min=False)
>>> R = PolynomialRing(T, names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3)
>>> p1 = x + y + z + R(Integer(0))
>>> tv = p1.tropical_variety(); tv
Tropical surface of 0*x + 0*y + 0*z + 0
>>> tv.components()
[[(t1, t1, t2), [t2 <= t1, 0 <= t1], 1],
[(t1, t2, t1), [max(0, t2) <= t1], 1],
[(0, t1, t2), [t2 <= 0, t1 <= 0], 1],
[(t1, t2, t2), [max(0, t1) <= t2], 1],
[(t1, 0, t2), [t2 <= 0, t1 <= 0], 1],
[(t1, t2, 0), [t1 <= 0, t2 <= 0], 1]]
plot(color='random')[source]

Return the plot of self by constructing a polyhedron from vertices in self.polygon_vertices().

INPUT:

  • color – string or tuple that represent a color (default: random); random means each polygon will be assigned a different color. If instead a specific color is provided, then all polygon will be given the same color.

OUTPUT: Graphics3d Object

EXAMPLES:

A tropical surface that consist of only one cell:

sage: T = TropicalSemiring(QQ)
sage: R.<x,y,z> = PolynomialRing(T)
sage: p1 = x + z
sage: tv = p1.tropical_variety()
sage: tv.plot()
Graphics3d Object
>>> from sage.all import *
>>> T = TropicalSemiring(QQ)
>>> R = PolynomialRing(T, names=('x', 'y', 'z',)); (x, y, z,) = R._first_ngens(3)
>>> p1 = x + z
>>> tv = p1.tropical_variety()
>>> tv.plot()
Graphics3d Object
../../../_images/tropical_variety-4.svg

A tropical surface with multiple cell that exhibit complex and intriguing geometric structures:

sage: p2 = x^2 + x + y + z + R(1)
sage: tv = p2.tropical_variety()
sage: tv.plot()
Graphics3d Object
>>> from sage.all import *
>>> p2 = x**Integer(2) + x + y + z + R(Integer(1))
>>> tv = p2.tropical_variety()
>>> tv.plot()
Graphics3d Object
../../../_images/tropical_variety-5.svg
class sage.rings.semirings.tropical_variety.TropicalVariety(poly)[source]

Bases: UniqueRepresentation, SageObject

A tropical variety in \(\RR^n\).

A tropical variety is defined as a corner locus of tropical polynomial function. This means it consist of all points in \(\RR^n\) for which the minimum (maximum) of the function is attained at least twice.

We represent the tropical variety as a list of lists, where the inner list consist of three parts. The first one is a parametric equations for tropical roots. The second one is the condition for parameters. The third one is the order of the corresponding component.

INPUT:

  • poly – a TropicalMPolynomial

ALGORITHM:

We need to determine a corner locus of this tropical polynomial function, which is all points \((x_1, x_2, \ldots, x_n)\) for which the maximum (minimum) is obtained at least twice. First, we convert each monomial to its corresponding linear function. Then for each two monomials of polynomial, we find the points where their values are equal. Since we attempt to solve the equality of two equations in \(n\) variables, the solution set will be described by \(n-1\) parameters.

Next, we need to check if the value of previous two monomials at the points in solution set is really the maximum (minimum) of function. We do this by solving the inequality of the previous monomial with all other monomials in the polynomial after substituting the parameter. This will give us the condition of parameters. Each of this condition is then combined by union operator. If this final condition is not an empty set, then it represent one component of tropical root. Then we calculate the weight of this particular component by the maximum of gcd of the numbers \(|i-k|\) and \(|j-l|\) for all pairs \((i,j)\) and \((k,l)\) such that the value of on this component is given by the corresponding monomials.

EXAMPLES:

We construct a tropical variety in \(\RR^2\), where it is called a tropical curve:

sage: T = TropicalSemiring(QQ, use_min=False)
sage: R.<x,y> = PolynomialRing(T)
sage: p1 = R(1)*x + x*y + R(0); p1
0*x*y + 1*x + 0
sage: tv = p1.tropical_variety(); tv
Tropical curve of 0*x*y + 1*x + 0
sage: tv.components()
[[(t1, 1), [t1 >= -1], 1], [(-1, t1), [t1 <= 1], 1], [(-t1, t1), [t1 >= 1], 1]]
sage: tv.vertices()
{(-1, 1)}
sage: tv.plot()
Graphics object consisting of 3 graphics primitives
>>> from sage.all import *
>>> T = TropicalSemiring(QQ, use_min=False)
>>> R = PolynomialRing(T, names=('x', 'y',)); (x, y,) = R._first_ngens(2)
>>> p1 = R(Integer(1))*x + x*y + R(Integer(0)); p1
0*x*y + 1*x + 0
>>> tv = p1.tropical_variety(); tv
Tropical curve of 0*x*y + 1*x + 0
>>> tv.components()
[[(t1, 1), [t1 >= -1], 1], [(-1, t1), [t1 <= 1], 1], [(-t1, t1), [t1 >= 1], 1]]
>>> tv.vertices()
{(-1, 1)}
>>> tv.plot()
Graphics object consisting of 3 graphics primitives
../../../_images/tropical_variety-6.svg

A slightly different result will be obtained if we use min-plus algebra for the base tropical semiring:

sage: T = TropicalSemiring(QQ, use_min=True)
sage: R.<x,y> = PolynomialRing(T)
sage: p1 = R(1)*x + x*y + R(0)
sage: tv = p1.tropical_variety(); tv
Tropical curve of 0*x*y + 1*x + 0
sage: tv.components()
[[(t1, 1), [t1 <= -1], 1], [(-1, t1), [t1 >= 1], 1], [(-t1, t1), [t1 <= 1], 1]]
sage: tv.plot()
Graphics object consisting of 3 graphics primitives
>>> from sage.all import *
>>> T = TropicalSemiring(QQ, use_min=True)
>>> R = PolynomialRing(T, names=('x', 'y',)); (x, y,) = R._first_ngens(2)
>>> p1 = R(Integer(1))*x + x*y + R(Integer(0))
>>> tv = p1.tropical_variety(); tv
Tropical curve of 0*x*y + 1*x + 0
>>> tv.components()
[[(t1, 1), [t1 <= -1], 1], [(-1, t1), [t1 >= 1], 1], [(-t1, t1), [t1 <= 1], 1]]
>>> tv.plot()
Graphics object consisting of 3 graphics primitives
../../../_images/tropical_variety-7.svg

Tropical variety can consist of multiple components with varying orders:

sage: T = TropicalSemiring(QQ, use_min=False)
sage: R.<x,y> = PolynomialRing(T)
sage: p1 = R(7) + T(4)*x + y + R(4)*x*y + R(3)*y^2 + R(-3)*x^2
sage: tv = p1.tropical_variety(); tv
Tropical curve of (-3)*x^2 + 4*x*y + 3*y^2 + 4*x + 0*y + 7
sage: tv.components()
[[(3, t1), [t1 <= 0], 1],
[(-t1 + 3, t1), [0 <= t1, t1 <= 2], 1],
[(t1, 2), [t1 <= 1], 2],
[(t1, 0), [3 <= t1, t1 <= 7], 1],
[(7, t1), [t1 <= 0], 1],
[(t1 - 1, t1), [2 <= t1], 1],
[(t1 + 7, t1), [0 <= t1], 1]]
sage: tv.plot()
Graphics object consisting of 8 graphics primitives
>>> from sage.all import *
>>> T = TropicalSemiring(QQ, use_min=False)
>>> R = PolynomialRing(T, names=('x', 'y',)); (x, y,) = R._first_ngens(2)
>>> p1 = R(Integer(7)) + T(Integer(4))*x + y + R(Integer(4))*x*y + R(Integer(3))*y**Integer(2) + R(-Integer(3))*x**Integer(2)
>>> tv = p1.tropical_variety(); tv
Tropical curve of (-3)*x^2 + 4*x*y + 3*y^2 + 4*x + 0*y + 7
>>> tv.components()
[[(3, t1), [t1 <= 0], 1],
[(-t1 + 3, t1), [0 <= t1, t1 <= 2], 1],
[(t1, 2), [t1 <= 1], 2],
[(t1, 0), [3 <= t1, t1 <= 7], 1],
[(7, t1), [t1 <= 0], 1],
[(t1 - 1, t1), [2 <= t1], 1],
[(t1 + 7, t1), [0 <= t1], 1]]
>>> tv.plot()
Graphics object consisting of 8 graphics primitives
../../../_images/tropical_variety-8.svg

If the tropical polynomial have \(n>2\) variables, then the result will be a tropical hypersurface embedded in a real space \(\RR^n\):

sage: T = TropicalSemiring(QQ)
sage: R.<a,x,y,z> = PolynomialRing(T)
sage: p1 = x*y + R(-1/2)*x*z + R(4)*z^2 + a*x
sage: tv = p1.tropical_variety(); tv
Tropical hypersurface of 0*a*x + 0*x*y + (-1/2)*x*z + 4*z^2
sage: tv.components()
[[(t1, t2, t3 - 1/2, t3), [t2 - 9/2 <= t3, t3 <= t1 + 1/2, t2 - 5 <= t1], 1],
[(t1, 2*t2 - t3 + 4, t3, t2), [t3 + 1/2 <= t2, t3 <= t1], 1],
[(t1, t2, t1, t3), [max(t1 + 1/2, 1/2*t1 + 1/2*t2 - 2) <= t3], 1],
[(t1, t2 + 9/2, t3, t2), [t2 <= min(t3 + 1/2, t1 + 1/2)], 1],
[(t1 - 1/2, t2, t3, t1), [t2 - 9/2 <= t1, t1 <= t3 + 1/2, t2 - 5 <= t3], 1],
[(2*t1 - t2 + 4, t2, t3, t1), [t1 <= min(1/2*t2 + 1/2*t3 - 2, t2 - 9/2)], 1]]
>>> from sage.all import *
>>> T = TropicalSemiring(QQ)
>>> R = PolynomialRing(T, names=('a', 'x', 'y', 'z',)); (a, x, y, z,) = R._first_ngens(4)
>>> p1 = x*y + R(-Integer(1)/Integer(2))*x*z + R(Integer(4))*z**Integer(2) + a*x
>>> tv = p1.tropical_variety(); tv
Tropical hypersurface of 0*a*x + 0*x*y + (-1/2)*x*z + 4*z^2
>>> tv.components()
[[(t1, t2, t3 - 1/2, t3), [t2 - 9/2 <= t3, t3 <= t1 + 1/2, t2 - 5 <= t1], 1],
[(t1, 2*t2 - t3 + 4, t3, t2), [t3 + 1/2 <= t2, t3 <= t1], 1],
[(t1, t2, t1, t3), [max(t1 + 1/2, 1/2*t1 + 1/2*t2 - 2) <= t3], 1],
[(t1, t2 + 9/2, t3, t2), [t2 <= min(t3 + 1/2, t1 + 1/2)], 1],
[(t1 - 1/2, t2, t3, t1), [t2 - 9/2 <= t1, t1 <= t3 + 1/2, t2 - 5 <= t3], 1],
[(2*t1 - t2 + 4, t2, t3, t1), [t1 <= min(1/2*t2 + 1/2*t3 - 2, t2 - 9/2)], 1]]
components()[source]

Return all components of self.

EXAMPLES:

sage: T = TropicalSemiring(QQ)
sage: R.<a,x,y,z> = PolynomialRing(T)
sage: tv = (a+x+y+z).tropical_variety()
sage: tv.components()
[[(t1, t1, t2, t3), [t1 <= min(t3, t2)], 1],
 [(t1, t2, t1, t3), [t1 <= t3, t1 <= t2], 1],
 [(t1, t2, t3, t1), [t1 <= min(t3, t2)], 1],
 [(t1, t2, t2, t3), [t2 <= t3, t2 <= t1], 1],
 [(t1, t2, t3, t2), [t2 <= min(t3, t1)], 1],
 [(t1, t2, t3, t3), [t3 <= min(t1, t2)], 1]]
>>> from sage.all import *
>>> T = TropicalSemiring(QQ)
>>> R = PolynomialRing(T, names=('a', 'x', 'y', 'z',)); (a, x, y, z,) = R._first_ngens(4)
>>> tv = (a+x+y+z).tropical_variety()
>>> tv.components()
[[(t1, t1, t2, t3), [t1 <= min(t3, t2)], 1],
 [(t1, t2, t1, t3), [t1 <= t3, t1 <= t2], 1],
 [(t1, t2, t3, t1), [t1 <= min(t3, t2)], 1],
 [(t1, t2, t2, t3), [t2 <= t3, t2 <= t1], 1],
 [(t1, t2, t3, t2), [t2 <= min(t3, t1)], 1],
 [(t1, t2, t3, t3), [t3 <= min(t1, t2)], 1]]
dimension()[source]

Return the dimension of self.

EXAMPLES:

sage: T = TropicalSemiring(QQ)
sage: R.<a,x,y,z> = PolynomialRing(T)
sage: p1 = x*y + R(-1)*x*z
sage: p1.tropical_variety().dimension()
4
>>> from sage.all import *
>>> T = TropicalSemiring(QQ)
>>> R = PolynomialRing(T, names=('a', 'x', 'y', 'z',)); (a, x, y, z,) = R._first_ngens(4)
>>> p1 = x*y + R(-Integer(1))*x*z
>>> p1.tropical_variety().dimension()
4
number_of_components()[source]

Return the number of components that make up self.

EXAMPLES:

sage: T = TropicalSemiring(QQ)
sage: R.<a,x,y,z> = PolynomialRing(T)
sage: p1 = x*y*a + x*z + y^2 + a*x + y + z
sage: p1.tropical_variety().number_of_components()
13
>>> from sage.all import *
>>> T = TropicalSemiring(QQ)
>>> R = PolynomialRing(T, names=('a', 'x', 'y', 'z',)); (a, x, y, z,) = R._first_ngens(4)
>>> p1 = x*y*a + x*z + y**Integer(2) + a*x + y + z
>>> p1.tropical_variety().number_of_components()
13