application: tropical
This application concentrates on tropical hypersurfaces and tropical cones. It provides the functionality for the computation of basic properties. Visualization and various constructions are possible.
uses: fan, group, ideal, matroid, polytope, topaz
Objects
A tropical cone is the tropical convex hull of finitely many points in tropical projective space. It should always be defined via POINTS instead of VERTICES, as those define the combinatorics of the induced subdivision.
Type Parameters
Addition Scalar Rational by default. The underlying type of ordered group.Properties of Cone
These properties capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.
- CONE_COVECTOR_DECOMPOSITION: CovectorLattice
This is a sublattice of TORUS_COVECTOR_DECOMPOSITION, containing only the cells that belong to the tropical span of POINTS.
- CONE_MAXIMAL_COVECTORS: common::Array<IncidenceMatrix<NonSymmetric>>
The covectors of the maximal cells of the cone subdivision. Entries correspond to rows of CONE_MAXIMAL_COVECTOR_CELLS.
- CONE_MAXIMAL_COVECTOR_CELLS: common::IncidenceMatrix<NonSymmetric>
This is a description of the tropical cone as a polyhedral complex. Each row is a maximal cell of the covector subdivision of the tropical cone. Indices refer to PSEUDOVERTICES.
- MAXIMAL_COVECTORS: common::Array<IncidenceMatrix<NonSymmetric>>
The covectors of the maximal cells of the torus subdivision. Entries correspond to rows of MAXIMAL_COVECTOR_CELLS.
- MAXIMAL_COVECTOR_CELLS: common::IncidenceMatrix<NonSymmetric>
These are the maximal cells of the covector decomposition of the tropical torus with respect to POINTS. Each row corresponds to a maximal cell, each column to an element of PSEUDOVERTICES.
- PSEUDOVERTEX_COARSE_COVECTORS: common::Matrix<Int, NonSymmetric>
Coarse types of PSEUDOVERTICES relative to POINTS. Each row corresponds to a row of PSEUDOVERTICES and encodes at position i, how many POINTS contain that pseudovertex in the i-th sector.
- PSEUDOVERTEX_COVECTORS: common::Array<IncidenceMatrix<NonSymmetric>>
Types of PSEUDOVERTICES relative to POINTS. Each type is encoded as an Incidence matrix, where rows correspond to coordinates and columns to POINTS. If the i-th row is a set S, that means that this pseudovertex is in the i-th sector of all points indexed by S. For bounded vertices, the type is computed as usual. For unbounded rays (i.e. starting with a 0), the type is computed as follows. Let g be a generator, with infinite entries at positions J and let the ray be e_J = sum_{j in J} +- e_j (the sign being the orientation of the addition). If J is contained in K, the ray is "contained" in all sectors of g. Otherwise, the ray is "contained" in the sectors indexed by g. NOTE: The latter is an artificial definition in the sense that it is not compatible with intersecting faces of the covector lattice. However, it is correct in the sense that faces spanned by a list of pseudovertices have as covector the intersection of the respective covectors.
- TORUS_COVECTOR_DECOMPOSITION: CovectorLattice
This is the face lattice of the polyhedral complex, whose vertices are PSEUDOVERTICES and whose cells are the cells of the covector decomposition. For each face in this lattice, we save the following information: 1) What PSEUDOVERTICES make up this face, i.e. a Set<int> 2) What is the covector of this face, i.e. an IncidenceMatrix (whose rows correspond to coordinates and whose columns to POINTS). NOTE: This lattice does not contain any far faces of the polyhedral cells, as they do not have well-defined covectors.
These properties capture geometric information of the object. Geometric properties depend on geometric information of the object, like, e.g., vertices or facets.
- DOME: polytope::Polytope
This is the dome of the tropical hyperplane arrangement defined by the POINTS. I.e. we take as function the (tropical) product of the tropical linear polynomials defined in the following manner: For each point (p_0,...,p_d) we get the linear polynomial sum_{i=1}^d (1/p_i) * x_i, where sum is the DUAL tropical addition and * and / is regular addition and subtraction, respectively.
- ENVELOPE: polytope::Polytope
Tropical cones have a natural description as the complex of certain faces of their envelopes. This envelope depends on the choice of the POINTS that generate the tropical polytope.
- POINTS: common::Matrix
Input points in tropical homogeneous coordinates. This is the fixed system of generators with respect to which many combinatorial properties are expressed.
- PROJECTIVE_AMBIENT_DIM: common::Int
Dimension of the tropical projective space which contains the tropical polytope.
User Methods of Cone
These methods capture geometric information of the object. Geometric properties depend on geometric information of the object, like, e.g., vertices or facets.
- cone_subdivision_as_complex (chart) → fan::PolyhedralComplex
This returns the subdivision of the cone induced by POINTS as a polyhedral complex on a chosen affine chart.
Parameters
Int chart Which coordinate to normalize to 0. This is 0 by default.Returns
fan::PolyhedralComplex - torus_subdivision_as_complex (chart) → fan::PolyhedralComplex
This returns the subdivision of the tropical torus induced by POINTS as a polyhedral complex on a chosen affine chart
Parameters
Int chart Which coordinate to normalize to 0. This is 0 by default.Returns
fan::PolyhedralComplex
These methods are for visualization.
- VISUAL_HYPERPLANE_ARRANGEMENT () → fan::Visual::PolyhedralFan
Visualize the arrangement of hyperplanes with apices in the POINTS of the tropical polytope.
- VISUAL_SUBDIVISION () → fan::Visual::PolyhedralFan
Visualize the subdivision of the torus induced by POINTS.
Permutations of Cone
- derived from: graph::FaceLattice
This is an augmented version of FaceLattice, which for each cell in a covector decomposition also saves the corresponding type.
Properties of CovectorLattice
These properties capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.
- COVECTORS: common::NodeMap<Directed, IncidenceMatrix<NonSymmetric>>
Each node in the face lattice is a cell of a covector decomposition (of either the tropical torus or the tropical span of some points). This property maps each cell to the corresponding covector. A covector is encoded as an IncidenceMatrix, where rows correspond to coordinates and columns to POINTS.
User Methods of CovectorLattice
These methods are for visualization.
- VISUAL ()
Visualizes the covector lattice. This works the same as a visualization of a Hasse diagram except that by default, covectors are displayed. This can be turned off by the option Covectors=>"hidden"
Options
option list: Visual::CovectorLattice::decorations
- derived from: fan::PolyhedralComplex<Rational>
A tropical cycle is a weighted, balanced, pure polyhedral complex. It is given as a polyhedral complex in tropical projective coordinates. To be precise: Each row of VERTICES and LINEALITY_SPACE has a leading 1 or 0, depending on whether it is a vertex or a ray. The remaining n coordinates are interpreted as an element of Rn modulo (1,..,1). IMPORTANT NOTE: VERTICES are assumed to be normalized such that the first coordinate (i.e. column index 1) is 0. If your input is not of that form, use PROJECTIVE_VERTICES. Note that there is a convenience method thomog, which converts affine coordinates into projective coordinates.
Type Parameters
Addition Properties of Cycle
These properties deal with affine and projective coordinates, conversion between those and properties like dimension that change in projective space.
- PROJECTIVE_AMBIENT_DIM: common::Int
This is the ambient projective dimension, i.e. it is FAN_AMBIENT_DIM-2.
- PROJECTIVE_CODIMENSION: common::Int
Codimension of the cycle. Same as PROJECTIVE_AMBIENT_DIM - PROJECTIVE_DIM
These properties capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.
- CODIMENSION_ONE_POLYTOPES: common::IncidenceMatrix<NonSymmetric>
Non-redundant list of all codimension one faces. Indices refer to VERTICES. Does not include any far faces.
- FACET_NORMALS_BY_PAIRS: common::Map<Pair<Int, Int>, Int>
This maps an index pair (i,j), where i corresponds to the i-th row of CODIMENSION_ONE_POLYTOPES and j to the j-th row of MAXIMAL_POLYTOPES, to the matching row number in FACET_NORMALS.
- MAXIMAL_AT_CODIM_ONE: common::IncidenceMatrix<NonSymmetric>
Incidence matrix of codimension one polytopes and maximal polytopes. Rows refer to CODIMENSION_ONE_POLYTOPES, columns to MAXIMAL_POLYTOPES.
These properties are used to define morphisms or rational functions on a Cycle.
Contained in extensionbundled:atint
.- SEPARATED_CODIMENSION_ONE_POLYTOPES: common::IncidenceMatrix<NonSymmetric>
An incidence matrix describing which codimension one polytope in the complex is generated by which vertices. Each row corresponds to a codimension one polytope (More precisely, the i-th element represents the same codim 1 polytope as the i-th element of CODIMENSION_ONE_POLYTOPES). The indices in a row refer to rows of SEPARATED_VERTICES.
- SEPARATED_CONVERSION_VECTOR: common::Vector<Int>
A vector with an entry for each row in SEPARATED_VERTICES. More precisely, the i-th entry gives the row index of the ray in VERTICES that is equal to the i-th row of SEPARATED_VERTICES.
- SEPARATED_MAXIMAL_POLYTOPES: common::IncidenceMatrix<NonSymmetric>
An incidence matrix describing which maximal polytope in the complex us generated by which rays. Each row corresponds to a maximal polytope (More precisely, the i-th element represents the same maximal polytope as the i-th element of MAXIMAL_POLYTOPES). The indices in a row refer to rows of SEPARATED_VERTICES, i.e. the maximal polytope described by the i-th element is generated by the vertices corresponding to these row indices.
- SEPARATED_VERTICES: common::Matrix<Rational, NonSymmetric>
This is a matrix of vertices of the complex. More precisely, each ray r from VERTICES occurs as a row in this matrix... - once, if r_0 = 1 - k times, if r_0 = 0 and k is the number of equivalence classes of maximal cones containing r with respect to the following relation: Two maximal cones m, m' containing r are equivalent, if they are equal or there exists a sequence of maximal cones m = m_1,...m_r = m', such that r is contained in each m_i and each intersection m_i cap m_i+1 contains at least one ray s with s_0 = 1. The reason for this is that, when for example specifying a piecewise affine linear function on a polyhedral complex, the same far ray with x0 = 0 might be assigned two different values, if it is contained in two "non-connected" maximal cones (where connectedness is to be understood as described above). If there is a LOCAL_RESTRICTION the above equivalence relation is changed in such a way that the affine ray s with s_0 = 1 that must be contained in the intersection of two subsequent cones must be a compatible ray
These properties are for input only. They allow redundant information.
- PROJECTIVE_VERTICES: common::Matrix<Rational, NonSymmetric>
A list of (non-redundant) vertices and rays of the complex. Note that VERTICES are supposed to be normalized to a leading 0. If you want to input non-normalized vertices, use this property. VERTICES will be derived from this: Each row of VERTICES is just the corresponding row of PROJECTIVE_VERTICES, with the first coordinate (i.e. column index 1) normalized to 0.
These are general properties related to intersection theory.
Contained in extensionbundled:atint
.- DEGREE: common::Integer
The degree of the tropical variety, i.e. the weight of the intersection product with a uniform tropical linear space of complementary dimension.
These properties are used for doing computations locally around a specified part of a Cycle. ----- These +++ deal with the creation and modification of cycles with nontrivial LOCAL_RESTRICTION.
Contained in extensionbundled:atint
.- LOCAL_RESTRICTION: common::IncidenceMatrix<NonSymmetric>
This contains a list of sets of ray indices (referring to VERTICES). All of these sets should describe polyhedra of the polyhedral complex (though not necessarily maximal ones). A polyhedron is now called compatible with this property, if it contains one of these polyhedra If this list is not empty, all computations will be done only on (or around) compatible cones. The documentation of each property will explain in what way this restriction is enforced. If this list is empty or not defined, there is no restriction. Careful: The implementation assumes that ALL maximal cones are compatible. If in doubt, you can create a complex with a local restriction from a given complex by using one of the "local_..." creation methods This list is assumed to be irredundant, i.e. there are no doubles (though this should not break anything, it is simply less efficient). It is, however, possible that one element is a subset of another.
These properties relate to the weights of a tropical cycle.
- BALANCED_FACES: common::Vector<Bool>
A vector whose entries correspond to the rows of CODIMENSION_ONE_POLYTOPES. The i-th entry is true, if and only if the complex is balanced at that face
Contained in extensionbundled:atint
. - IS_BALANCED: common::Bool
Whether the cycle is balanced. As many functions in a-tint can deal with non-balanced complexes, we include this as a property.
Contained in extensionbundled:atint
. - IS_IRREDUCIBLE: common::Bool
Whether this complex is irreducible.
Contained in extensionbundled:atint
. - LATTICE_BASES: common::IncidenceMatrix<NonSymmetric>
This incidence matrix gives a lattice basis for each maximal polytope. More precisely it gives a lattice basis whose span contains the lattice of the maximal polytope. Row i corresponds to cone i and gives lattice generator indices referring to LATTICE_GENERATORS. If this property is computed via rules, it does indeed give a lattice basis for the cone lattice, but when it is computed during an operation like refinement or divisor it will in general be larger. If this property exists, lattice normals might be computed faster.
Contained in extensionbundled:atint
. - LATTICE_GENERATORS: common::Matrix<Integer, NonSymmetric>
This is an irredundant list of all lattice generators of all maximal polyhedra. If this property exists, lattice normals might be computed faster
Contained in extensionbundled:atint
. - LATTICE_NORMALS: common::Map<Pair<Int, Int>, Vector<Integer>>
The lattice normals of codimension one faces with respect to adjacent maximal cells. It maps a pair of indices (i,j) to the lattice normal of the codimension one face given by row i in CODIMENSION_ONE_POLYTOPES in the maximal cell given by row j in MAXIMAL_POLYTOPES. The lattice normal is a representative of a generator of the quotient of the saturated lattice of the maximal cell by the saturated lattice of the codimension one face. It is chosen such that it "points into the maximal cell" and is only unique modulo the lattice spanned by the codimension one cell.
- LATTICE_NORMAL_FCT_VECTOR: common::Map<Pair<Int, Int>, Vector<Rational>>
For each lattice normal vector, this gives a vector of length (number of rays) + (lineality dim.), such that if a rational function is given by values on the rays and lin space generators, the value of the corresponding normal LATTICE_NORMALS->{i}->{j} can be computed by multiplying the function value vector with the vector LATTICE_NORMAL_FCT_VECTOR->{i}->{j}. This is done in the following way: We use the generating system (and indices refer to SEPARATED_VERTICES) <(r_i-r_0)_i>0, s_j, l_k>, where r_0 is the ray of the maximal cone with the lowest index in SEPARATED_VERTICES, such that it fulfills x0 = 1, r_i are the remaining rays with x0 = 1, ordered according to their index in SEPARATED_VERTICES, s_j are the rays of the cone with x0 = 0 and l_k are the lineality space generators. We will then store the coefficients a_i of (r_i - r_0) at the index of r_i, then - sum(a_i) at the index of r_0 and the remaining coefficients at the appropriate places. In particular, the value of a lattice normal under a rational function can be computed simply by taking the scalar product of RAY_VALUES | LIN_VALUES with this FCT_VECTOR
Contained in extensionbundled:atint
. - LATTICE_NORMAL_SUM: common::Matrix<Rational, NonSymmetric>
Rows of this matrix correspond to CODIMENSION_ONE_POLYTOPES, and each row contains the weighted sum: sum_{cone > codim-1-face}( weight(cone) * LATTICE_NORMALS->{codim-1-face}->{cone})
Contained in extensionbundled:atint
. - LATTICE_NORMAL_SUM_FCT_VECTOR: common::Matrix<Rational, NonSymmetric>
Rows of this matrix correspond to SEPARATED_CODIMENSION_ONE_POLYTOPES and each row contains a function vector for the corresponding row of LATTICE_NORMAL_SUM. This function vector is computed in the same way as described under LATTICE_NORMAL_FCT_VECTOR. Note that for any codim-1-faces at which the complex is not balanced, the corresponding row is a zero row. If a face is balanced can be checked under BALANCED_FACES.
Contained in extensionbundled:atint
. - WEIGHTS: common::Vector<Integer>
These are the integer weights associated to the maximal cells of the complex. Entries correspond to (rows of) MAXIMAL_POLYTOPES.
- WEIGHT_CONE: polytope::Cone<Rational>
The intersection of WEIGHT_SPACE with the positive orthant.
Contained in extensionbundled:atint
. - WEIGHT_SPACE: common::Matrix<Rational, NonSymmetric>
A Z-basis (as rows) for the space of weight distributions on this tropical cycle making it balanced (i.e. this cycle is irreducible, if and only if WEIGHT_SPACE has only one row and the gcd of WEIGHTS is 1.
Contained in extensionbundled:atint
. - WEIGHT_SYSTEM: common::Matrix<Rational, NonSymmetric>
The dual of WEIGHT_SPACE.
Contained in extensionbundled:atint
.
User Methods of Cycle
These methods deal with affine and projective coordinates, conversion between those and properties like dimension that change in projective space.
- affine_chart (chart) → fan::PolyhedralComplex<Rational>
This produces a version of the cycle in the coordinates of a standard tropical chart, i.e. one coordinate is set to 0. It is returned as an ordinary polyhedral complex (which can, for example, be used for visualization).
Parameters
Int chart The coordinate which should be set to 0. Indexed from 0 to AMBIENT_DIM-1 and 0 by default.Returns
fan::PolyhedralComplex<Rational>
These methods provide basic functionality related to polyhedral geometry, but not necessarily to tropical geometry
Contained in extensionbundled:atint
.These methods capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.
- facet_normal (i, j) → Int
Convenience function to ask for FACET_NORMALS_BY_PAIRS->{new Pair<Int,Int>(i,i)}
Parameters
Int i Row index in CODIMENSION_ONE_POLYTOPES.Int j Row index in MAXIMAL_POLYTOPES.Returns
Int Row index in FACET_NORMALS.
These methods are used for doing computations locally around a specified part of a Cycle. ----- These +++ deal with the creation and modification of cycles with nontrivial LOCAL_RESTRICTION.
Contained in extensionbundled:atint
.- delocalize () → Cycle
Returns the cycle without its LOCAL_RESTRICTION (Note that only the defining properties are kept. All derived information is lost).
Returns
Cycle
These methods are for visualization.
- BB_VISUAL ()
Displays a (possibly weighted) polyhedral complex by intersecting it with a bounding box. This bounding box is either defined by the vertices of the complex and the option "BoundingDistance" or explicitly given by "BoundingBox" and by setting "BoundingMode" to "absolute" @options
Contained in extensionbundled:atint
.Options
Int Chart Which affine chart to visualize, i.e. which coordinate to shift to 0. This is 0 by default.String WeightLabels If "hidden", no weight labels are displayed. Not hidden by default.String CoordLabels If "show", coordinate labels are displayed at vertices. Hidden by default.String BoundingMode Can be "relative" (intersects with the bounding box returned by the method boundingBox(BoundingDistance)) or "absolute" (intersects with the given BoundingBox). "relative" by default.Rational BoundingDistance The distance parameter for relative bounding modeMatrix<Rational> BoundingBox The bounding parameter for absolute bounding modeoption list: Visual::Cycle::BoundingDecorations - bounding_box (dist, chart)
Takes a chart and a positive Rational as input and computes the relative bounding box of the cycle, i.e. it takes the coordinate-wise minimum and maximum over the coordinates of the nonfar vertices and adds/subtracts the given number. This is returned as a 2xdim matrix.
Contained in extensionbundled:atint
. - VISUAL () → fan::Visual::PolyhedralFan
Visualizes the tropical cycle. The options are more or less the same as for PolyhedralComplex's VISUAL, except that there is one option for choosing a chart. @options
Options
Int Chart . Which affine chart of the cycle to visualize, i.e. which coordinate should be shifted to 0. This is 0 by default.option list: Visual::Cycle::decorations Returns
fan::Visual::PolyhedralFan
These methods relate to the weights of a tropical cycle.
- lattice_normal (i, j) → Vector<Integer>
Convenience function to ask for LATTICE_NORMALS->{new Pair<Int,Int>(i,i)}
Parameters
Int i Row index in CODIMENSION_ONE_POLYTOPES.Int j Row index in MAXIMAL_POLYTOPES.Returns
Vector<Integer>
Permutations of Cycle
- derived from: Cycle
This is a special instance of a Cycle: It is the tropical locus of a polynomial over the tropical numbers.
Properties of Hypersurface
- COEFFICIENTS: common::Vector
Each row corresponds to one of the monomials in POLYNOMIAL, each column to a variable.
- DOME: polytope::Polytope<Rational>
The dome of a tropical polynomial \(F:\mathbb R^d\to\mathbb R\) (and the corresponding tropical hypersurface) is the set \[D(F)=\left\{(p,s)\in\mathbb R^{d+1}\mid p\in\mathbb R^d, s\in\mathbb R, s \oplus F(p) = s\right\}.\] It is an unbounded convex polyhedron, c.f.
Michael Joswig, Essentials of Tropical Combinatorics, Chapter 1. - MONOMIALS: common::Matrix<Int, NonSymmetric>
Each row corresponds to one of the monomials in POLYNOMIAL, each column to a variable.
User Methods of Hypersurface
- dual_subdivision () → fan::PolyhedralComplex<Rational>
Subdivision of the Newton polytope dual to the tropical hypersurface. The vertices of this PolyhedralComplex are the non-redundant MONOMIALS.
This represents the result of the method lines_in_cubic. It contains: The tropical polynomial representing the surface, the surface itself as a Cycle and lists of lines and families of different types, each starting with LIST_...
The object also has methods, starting with array_... that return the corresponding LIST_... as a perl array. The different (lists of) lines can be visualized nicely with visualize_in_surface.
Contained in extensionbundled:atint
.Properties of LinesInCubic
These count lines or families of lines in the cubic.# @category Defining properties
- N_FAMILIES: common::Int
The total number of families in LIST_FAMILY_FIXED_EDGE, LIST_FAMILY_FIXED_VERTEX, LIST_FAMILY_MOVING_EDGE and LIST_FAMILY_MOVING_VERTEX
- N_ISOLATED: common::Int
The total number of elements in LIST_ISOLATED_EDGE and LIST_ISOLATED_NO_EDGE
The polynomial and the corresponding cubic.
These contain lists of certain (families of) lines.
- LIST_FAMILY_FIXED_VERTEX: Cycle
A list of all families of lines without a bounded edge at a fixed vertex
- LIST_FAMILY_MOVING_EDGE: Cycle
A list of all families of lines with the bounded edge moving in transversal direction
- LIST_FAMILY_MOVING_VERTEX: Cycle
A list of all families of lines without a bounded edge and a moving vertex.
User Methods of LinesInCubic
These contain lists of certain (families of) lines.
- all_families () → Cycle<Addition>
- all_isolated () → Cycle<Addition>
- array_family_fixed_edge () → Cycle<Addition>
- array_family_fixed_vertex () → Cycle<Addition>
- array_family_moving_edge () → Cycle<Addition>
- array_family_moving_vertex () → Cycle<Addition>
- array_isolated_edge () → Cycle<Addition>
- array_isolated_no_edge () → Cycle<Addition>
A morphism is a function between cycles which is locally affine linear and respects the lattices. It is defined by a DOMAIN, which is a cycle, and values on this domain, VERTEX_VALUES and LINEALITY_VALUES, much like RationalFunction. Alternatively, it can be defined as a global affine linear function by giving a matrix and a translation vector.
Contained in extensionbundled:atint
.Properties of Morphism
These properties are used to define morphisms or rational functions on a Cycle.
- IS_GLOBALLY_AFFINE_LINEAR: common::Bool
This is TRUE, if the morphism is defined on the full projective torus by a MATRIX and a TRANSLATE The rules do not actually check for completeness of the DOMAIN. This property will be set to TRUE, if the morphism is only defined by MATRIX and TRANSLATE, otherwise it is false (or you can set it upon creation).
- LINEALITY_VALUES: common::Matrix<Rational, NonSymmetric>
The vector in row i describes the function value (slope) of DOMAIN->LINEALITY_SPACE->row(i)
- MATRIX: common::Matrix<Rational, NonSymmetric>
If the morphism is a global affine linear map x |-> Ax+v, then this contains the matrix A. Note that this must be well-defined on tropical projective coordinates, so the sum of the columns of A must be a multiple of the (1,..,1)-vector. If TRANSLATE is set, but this property is not set, then it is the identity by default.
- TRANSLATE: common::Vector<Rational>
If the morphism is a global affine linear map x |-> Ax+v, then this contains the translation vector v. If MATRIX is set, but this property is not set, then it is the zero vector by default.
- VERTEX_VALUES: common::Matrix<Rational, NonSymmetric>
The vector at row i describes the function value of vertex DOMAIN->SEPARATED_VERTICES->row(i). (In tropical homogenous coordinates, but without leading coordinate). More precisely, if the corresponding vertex is not a far ray, it describes its function value. If it is a directional ray, it describes the slope on that ray.
User Methods of Morphism
These methods deal with affine and projective coordinates, conversion between those and properties like dimension that change in projective space.
- affine_representation (domain_chart, target_chart) → Pair<Matrix<Rational>, Vector<Rational> >
Parameters
Int domain_chart Which coordinate index of the homogenized domain is shifted to zero to identify it with the domain of the affine function. 0 by default.Int target_chart Which coordinate of the homogenized target space is shifted to zero to identify it with the target of the affine function. 0 by default.Returns
Pair<Matrix<Rational>, Vector<Rational> > A matrix and a translate in affine coordinates.
These are general methods that deal with morphisms and their arithmetic.
These methods are for visualization.
An n-marked rational curve, identified by its SETS, i.e. its partitions of {1,...,n} and its COEFFICIENTS, i.e. the lengths of the corresponding edges.
Contained in extensionbundled:atint
.Properties of RationalCurve
These properties capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.
- COEFFS: common::Vector<Rational>
A list of positive rational coefficients. The list should have the same length as SETS and contain only entries > 0. The i-th entry then gives the length of the bounded edge defined by the i-th partition. If you're not sure if all your coefficients are > 0, use INPUT_SETS and INPUT_COEFFS instead. Note that the zero curve (i.e. no bounded edges, only leaves) is represented by one empty set with corresponding lenghth 0.
- GRAPH: graph::Graph<Undirected>
Contains the abstract graph (non-metric) corresponding to the curve. All unbounded leaves are modelled as bounded edges. The vertices at the ends of the "leaves" are always the first N_LEAVES vertices.
- NODES_BY_LEAVES: common::IncidenceMatrix<NonSymmetric>
This incidence matrix gives a list of the vertices of the curve Each row corresponds to a vertex and contains as a set the [[LEAVES] that are attached to that vertex (again, counting from 1!)
- NODES_BY_SETS: common::IncidenceMatrix<NonSymmetric>
This incidence matrix gives a list of the vertices of the curve Each row corresponds to a vertex and contains as a set the row indices of the SETS that correspond to edges attached to that vertex
- NODE_DEGREES: common::Vector<Int>
This gives a list of the vertices of the curve in terms of their valences They appear in the same order as in NODES_BY_LEAVES or NODES_BY_SETS
- SETS: common::IncidenceMatrix<NonSymmetric>
A list of partitions of [n] that define the tree of the curve: For each bounded edge we have the corresponding partition of the n leaves. These should be irredundant. If you want to input a possibly redundant list, use INPUT_SETS and INPUT_COEFFS instead. The number of marked leaves should always be given by N_LEAVES. The sets are subsets of {1,...,n} (NOT {0,..,n-1}!) Note that the zero curve (i.e. no bounded edges, only leaves) is represented by one empty set with corresponding lenghth 0.
These properties are for input only. They allow redundant information.
- INPUT_COEFFS: common::Vector<Rational>
Same as COEFFS, except that entries may be <=0. This should have the same length as INPUT_SETS.
- INPUT_SETS: common::IncidenceMatrix<NonSymmetric>
Same as SETS, except that sets may appear several times.
- INPUT_STRING: common::String
This property can also be used to define a rational curve: A linear combination of partitions is given as a string, using the following syntax: A partition is given as a subset of {1,..,n} and written as a comma-separated list of leaf indices in round brackets, e.g. "(1,2,5)" A linear combination can be created using rational numbers, "+","+" and "-" in the obvious way, e.g. "2*(1,2,5) + 1*(3,4,7) - 2(1,2) (The "*" is optional) Of course, each set should contain at least two elements. If you don't specify N_LEAVES, it is set to be the largest leaf index occuring in the sets. Partitions needn't be irredundant and coefficients can be any rational number. If the resulting element is not in the moduli space, an error is thrown.
User Methods of RationalCurve
These deal with converting the representation of a rational curve between metric vector and matroid fan coordinates.
- metric_vector ()
Returns the (n over 2) metric vector of the rational n-marked curve
These methods are for visualization.
- VISUAL ()
Visualizes a RationalCurve object. This visualization uses the VISUAL method of its GRAPH, so it accepts all the options of Visual::Graph::decorations. In addition it has another option @options
Options
String LengthLabels If "hidden", the edges are not labelled with their lengths. Any other text is ignored. Not set to "hidden" by default.option list: Visual::RationalCurve::decorations
A rational function on a polyhedral complex. It can be described by giving its DOMAIN, a Cycle, and values on this domain - which are encoded in the properties VERTEX_VALUES and LINEALITY_VALUES. Alternatively, it can be defined by a tropical quotient of homogeneous tropical polynomials of the same degree i.e. by giving NUMERATOR and DENOMINATOR. A DOMAIN can be defined additionally (though one should take care that both functions are actually piecewise affine linear on the cells), otherwise it will be computed as the common refinement of the domains of affine linearity of the two polynomials. Note: This has nothing to do with common's RationalFunction (which is univariate). If you want to access that type or use this type from another application, be sure to prepend the appropriate namespace identifier.
Contained in extensionbundled:atint
.Properties of RationalFunction
These properties are used to define morphisms or rational functions on a Cycle.
- DENOMINATOR: common::Polynomial
When representing the function as a quotient of tropical polynomials, this is the denominator. Should be a homogeneous polynomial of the same degree as NUMERATOR.
- IS_GLOBALLY_DEFINED: common::Bool
This is TRUE, if the function is defined on the full projective torus by a NUMERATOR and a DENOMINATOR. The rules do not actually check for completeness of the DOMAIN. This property will be set to true, if the function is created only via NUMERATOR and DENOMINATOR. Otherwise it will be set to FALSE (or you can set it manually upon creation).
- LINEALITY_VALUES: common::Vector<Rational>
The value at index i describes the function value of DOMAIN->LINEALITY_SPACE->row(i)
- NUMERATOR: common::Polynomial
When representing the function as a quotient of tropical polynomials, this is the numerator. Should be a homogeneous polynomial of the same degree as DENOMINATOR.
- POWER: common::Int
This is an internally used property that should not actually be set by the user. When creating a rational function with the ^-operator, this property is set to the exponent. The semantics is that when computing a divisor, this function should be applied so many times The usual application of this is a call to divisor($X, $f^4) or something similar. Warning: This property is not stored if the RationalFunction object is saved. Nor should be assumed to be preserved during any kind of arithmetic or restricting operation.
- VERTEX_VALUES: common::Vector<Rational>
The value at index i describes the function value at DOMAIN->SEPARATED_VERTICES->row(i). More precisely, if the corresponding vertex is not a far ray, it describes its function value. If it is a directional ray, it describes the slope on that ray.
User Methods of RationalFunction
These methods are used to define morphisms or rational functions on a Cycle.
- restrict (C) → RationalFunction<Addition>
These methods are for visualization.
User Functions
These functions deal with abstract rational n-marked curves.
Contained in extensionbundled:atint
.- insert_leaves (curve, nodes)
Takes a RationalCurve and a list of node indices. Then inserts additional leaves (starting from N_LEAVES+1) at these nodes and returns the resulting RationalCurve object
Parameters
RationalCurve curve A RationalCurve objectVector<Int> nodes A list of node indices of the curve - matroid_coordinates_from_curve <Addition> (r) → Vector<Rational>
Takes a rational curve and converts it into the corresponding matroid coordinates in the moduli space of rational curves (including the leading 0 for a ray!)
Type Parameters
Addition Min or Max, i.e. which coordinates to use.Parameters
RationalCurve r A rational n-marked curveReturns
Vector<Rational> The matroid coordinates, tropically homogeneous and with leading coordinate - rational_curve_from_cone (X, n_leaves, coneIndex) → RationalCurve
This takes a weighted complex X that is supposed to be of the form M_0,n x Y for some Y (It assumes that M_0,n occupies the first coordinates) and an index of a maximal cone of that complex. It then computes a rational curve corresponding to an interior point of that cone (ignoring the second component Y)
Parameters
Cycle<Addition> X A weighted complex of the form M_0,n x YInt n_leaves The n in M_0,n. Needed to determine the dimension of the M_0,n componentInt coneIndex The index of the maximal coneReturns
RationalCurve c The curve corresponding to an interior point - rational_curve_from_matroid_coordinates <Addition> (v) → RationalCurve
Takes a vector from Q^((n-1) over 2) that lies in M_0,n (in its matroid coordinates) and computes the corresponding rational curve. In the isomorphism of the metric curve space and the moduli coordinates the last leaf is considered as the special leaf
Type Parameters
Addition Min or Max (i.e. what are the matroid coordinates using)Parameters
Vector<Rational> v A vector in the moduli space (WITH leading coordinate).Returns
RationalCurve - rational_curve_from_metric (v) → RationalCurve
Takes a vector from Q^(n over 2) that describes an n-marked rational abstract curve as a distance vector between its leaves. It then computes the curve corresponding to this vector.
Parameters
Vector<Rational> v A vector of length (n over 2). Its entries are interpreted as the distances d(i,j) ordered lexicographically according to i,j. However, they need not be positive, as long as v is equivalent to a proper metric modulo leaf lengths.Returns
RationalCurve - rational_curve_from_rays <Addition> (rays) → RationalCurve
This takes a matrix of rays of a given cone that is supposed to lie in a moduli space M_0,n and computes the rational curve corresponding to an interior point. More precisely, if there are k vertices in homogeneous coordinates, it computes 1/k * (sum of these vertices), then it adds each directional ray. It then returns the curve corresponding to this point
Type Parameters
Addition Min or Max, where the coordinates live.Parameters
Matrix<Rational> rays The rays of the cone, in tropical homogeneous coordinates.Returns
RationalCurve c The curve corresponding to an interior point - rational_curve_immersion <Addition> (delta, type) → Cycle<Addition>
This function creates an embedding of a rational tropical curve using a given abstract curve and degree
Type Parameters
Addition Min or MaxParameters
Matrix<Rational> delta The degree of the curve in tropical projectve coordinates without leading coordinate. The number of rows should correspond to the number of leaves of type and the number of columns is the dimension of the space in which the curve should be realizedRationalCurve type An abstract rational curveReturns
Cycle<Addition> The corresponding immersed complex. The position of the curve is determined by the first node, which is always placed at the origin - rational_curve_list_from_matroid_coordinates <Addition> (m) → RationalCurve
Takes a matrix whose rows are elements in the moduli space M_0,n in matroid coordinates. Returns a list, where the i-th element is the curve corr. to the i-th row in the matrix
Type Parameters
Addition Mir or Max (i.e. what are the matroid coordinates usingParameters
Matrix<Rational> m A list of vectors in the moduli space (with leading coordinate).Returns
RationalCurve : An array of RationalCurves - rational_curve_list_from_metric (m) → RationalCurve
Takes a matrix whose rows are metrics of rational n-marked curves. Returns a list, where the i-th element is the curve corr. to the i-th row in the matrix
- sum_curves (An, v) → RationalCurve
This function takes a vector of coefficients a_i and a list of RationalCurves c_i and computes sum(a_i * c_i). In particular, it also checks, whether the result lies in M_0,n. If not, it returns undef
Parameters
RationalCurve An arbitrary list of RationalCurve objectsVector<Rational> v A list of coefficients. Superfluous coefficients are ignored, missing ones replaced by +1(!)Returns
RationalCurve The linear combination of the curves defined by the coefficients or undef, if the result is not in M_0,n. The history of the operation is kept in INPUT_SETS and INPUT_COEFFS - testFourPointCondition (v) → Int
Takes a metric vector in Q^{(n over 2)} and checks whether it fulfills the four-point condition, i.e. whether it lies in M_0,n. More precisely it only needs to be equivalent to such a vector
Parameters
Vector<Rational> v The vector to be checkedReturns
Int A quadruple (array) of indices, where the four-point condition is violated or an empty list, if the vector is indeed in M_0,n
These functions deal with affine and projective coordinates, conversion between those and properties like dimension that change in projective space.
- morphism_from_affine <Addition> (A, v, domain_chart, target_chart) → Morphism
Takes a representation of a morphism on affine coordinates and converts it to projective ones.
Contained in extensionbundled:atint
.Type Parameters
Addition Min or MaxParameters
Matrix<Rational> A . The matrix of the morphism x |-> Ax + v in affine coordinates.Vector<Rational> v . The translate of the morphism x |-> Ax + v in affine coordinates.Int domain_chart Which coordinate index of the homogenized domain is shifted to zero to identify it with the domain of the affine function. 0 by default.Int target_chart Which coordinate of the homogenized target space is shifted to zero to identify it with the target of the affine function. 0 by default.Returns
Morphism - rational_fct_from_affine_denominator (p, chart) → RationalFunction
This takes a tropical polynomial p defined on tropical affine coordinates and turns it into the rational function (1/p) on tropical homogeneous coordinates
Contained in extensionbundled:atint
.Parameters
Polynomial<TropicalNumber<Addition> > p A polynomial on affine coordinates.Int chart The index of the homogenizing coordinate. 0 by default.Returns
RationalFunction A rational function, which on the given chart is described by (1/p). - rational_fct_from_affine_denominator (p, chart) → RationalFunction
Same as rational_fct_from_affine_denominator(Polynomial), except that it takes a string which it converts to a tropical polynomial using toTropicalPolynomial.
Contained in extensionbundled:atint
.Parameters
String p A string that will be converted to a tropical polynomialInt chart The index of the homogenizing coordinate. 0 by default.Returns
RationalFunction - rational_fct_from_affine_numerator (p, chart) → RationalFunction
This takes a tropical polynomial defined on tropical affine coordinates and turns it into a rational function on tropical homogeneous coordinates
Contained in extensionbundled:atint
.Parameters
Polynomial<TropicalNumber<Addition> > p A polynomial on affine coordinates.Int chart The index of the homogenizing coordinate. 0 by default.Returns
RationalFunction A rational function, which on the given chart is described by p. - rational_fct_from_affine_numerator (p, chart) → RationalFunction
Same as rational_fct_from_affine_numerator(Polynomial), except that it takes a string which it converts to a tropical polynomial using toTropicalPolynomial.
Contained in extensionbundled:atint
.Parameters
String p A string that will be converted to a tropical polynomialInt chart The index of the homogenizing coordinate. 0 by default.Returns
RationalFunction - tdehomog (A, chart, has_leading_coordinate) → Matrix<Rational>
This is the inverse operation of thomog. It assumes a list of rays and vertices is given in tropical projective coordinates and returns a conversion into affine coordinates.
Parameters
Matrix<Rational> A The matrix. Can also be given as an anonymous array.Int chart Optional. Indicates which coordinate should be shifted to 0. If there is a leading coordinate, the first column of the matrix will remain untouched and the subsequent ones are numbered from 0. The default value for this is 0.Bool has_leading_coordinate Whether the matrix has a leading 1/0 to indicate whether a row is a vertex or a ray. In that case, this coordinate is not touched. This is true by default.Returns
Matrix<Rational> - thomog (A, chart, has_leading_coordinate) → Matrix<Rational>
Converts tropical affine to tropical projective coordinates. It takes a matrix of row vectors in Rn-1 and identifies the latter with Rn mod (1,..,1) by assuming a certain coordinate has been set to 0. I.e. it will return the matrix with a 0 column inserted at the position indicated by chart
Parameters
Matrix<Rational> A The matrix. Can also be given as an anonymous array [[..],[..],..]Int chart Optional. Indicates, which coordinate of Rn mod (1,..,1) should be set to 0 to identify it with Rn-1. Note that if there is a leading coordinate, the first column is supposed to contain the 1/0-coordinate indicating whether a row is a vertex or a ray and the remaining coordinates are then labelled 0,..,n-1. This option is 0 by default.Bool has_leading_coordinate Whether the matrix has a leading 1/0 to indicate whether a row is a vertex or a ray. In that case, this coordinate is not touched. This is true by default.Returns
Matrix<Rational>
These functions provide basic functionality related to polyhedral geometry, but not necessarily to tropical geometry
Contained in extensionbundled:atint
.- affine_transform (C, M, T) → Cycle<Addition>
Computes the affine transform of a cycle under an affine linear map. This function assumes that the map is a lattice isomorphism on the cycle, i.e. no push-forward computations are performed, in particular the weights remain unchanged
Parameters
Cycle<Addition> C a tropical cycleMatrix<Rational> M The transformation matrix. Should be given in tropical projective coordinates and be homogeneous, i.e. the sum over all rows should be the same.Vector<Rational> T The translate. Optional and zero vector by default. Should be given in tropical projective coordinates (but without leading coordinate for vertices or rays). If you only want to shift a cycle, use shift_cycle.Returns
Cycle<Addition> The transform M*C + T - affine_transform (C, M) → Cycle<Addition>
Computes the affine transform of a cycle under an affine linear map. This function assumes that the map is a lattice isomorphism on the cycle, i.e. no push-forward computations are performed, in particular the weights remain unchanged
Parameters
Cycle<Addition> C a tropical cycleMorphism<Addition> M Returns
Cycle<Addition> The transform M(C) - cartesian_product (cycles) → Cycle
Computes the cartesian product of a set of cycles. If any of them has weights, so will the product (all non-weighted cycles will be treated as if they had constant weight 1)
Parameters
Cycle cycles a list of CyclesReturns
Cycle The cartesian product. Note that the representation is noncanonical, as it identifies the product of two projective tori of dimensions d and e with a projective torus of dimension d+e by dehomogenizing and then later rehomogenizing after the first coordinate. - check_cycle_equality (X, Y, check_weights) → Bool
This takes two pure-dimensional polyhedral complexes and checks if they are equal i.e. if they have the same lineality space, the same rays (modulo lineality space) and the same cones. Optionally, it can also check if the weights are equal
Parameters
Cycle<Addition> X A weighted complexCycle<Addition> Y A weighted complexBool check_weights Whether the algorithm should check for equality of weights. This parameter is optional and true by defaultReturns
Bool Whether the cycles are equal - coarsen (complex, testFan) → Cycle<Addition>
Takes a tropical variety on which a coarsest polyhedral structure exists and computes this structure.
Parameters
Cycle<Addition> complex A tropical variety which has a unique coarsest polyhedral structreBool testFan (Optional, FALSE by default). Whether the algorithm should perform some consistency checks on the result. If true, it will check the following: - That equivalence classes of cones have convex support - That all equivalence classes have the same lineality space If any condition is violated, the algorithm throws an exception Note that it does not check whether equivalence classes form a fan This can be done via fan::check_fan afterwards, but it is potentially slow.Returns
Cycle<Addition> The corresponding coarse complex. Throws an exception if testFan = True and consistency checks fail. - contains_point (A, point) → Bool
Takes a weighted complex and a point and computed whether that point lies in the complex
Parameters
Cycle A weighted complexVector<Rational> point An arbitrary vector in the same ambient dimension as complex. Given in tropical projective coordinates with leading coordinate.Returns
Bool Whether the point lies in the support of complex - fan_decomposition (C) → Cycle<Addition>
This computes the local fans at all (nonfar) vertices of a tropical cycle
- insert_rays (F, R) → Cycle<Addition>
Takes a cycle and a list of rays/vertices in tropical projective coordinates with leading coordinate and triangulates the fan such that it contains these rays
Parameters
Cycle<Addition> F A cycle (not necessarily weighted).Matrix<Rational> R A list of normalized vertices or rays Note that the function will NOT subdivide the lineality space, i.e. rays that are equal to an existing ray modulo lineality space will be ignored.Returns
Cycle<Addition> A triangulation of F that contains all the original rays of F plus the ones in R - intersect_container (cycle, container, forceLatticeComputation) → Cycle
Takes two Cycles and computes the intersection of both. The function relies on the fact that the second cycle contains the first cycle to compute the refinement correctly The function copies WEIGHTS, LATTICE_BASES and LATTICE_GENERATORS in the obvious manner if they exist.
Parameters
Cycle cycle An arbitrary CycleCycle container A cycle containing the first one (as a set) Doesn't need to have any weights and its tropical addition is irrelevant.Bool forceLatticeComputation Whether the properties LATTICE_BASES and LATTICE_GENERATORS of cycle should be computed before refining. False by default.Returns
Cycle The intersection of both complexes (whose support is equal to the support of cycle). It uses the same tropical addition as cycle. - recession_fan (complex) → Cycle
- set_theoretic_intersection (A, B) → fan::PolyhedralComplex
Computes the set-theoretic intersection of two cycles and returns it as a polyhedral complex. The cycles need not use the same tropical addition
Parameters
Cycle A Cycle B Returns
fan::PolyhedralComplex The set-theoretic intersection of the supports of A and B - shift_cycle (C, T) → Cycle<Addition>
Computes the shift of a tropical cycle by a given vector
Parameters
Cycle<Addition> C a tropical cycleVector<Rational> T The translate. Optional and zero vector by default. Should be given in tropical projective coordinates (but without leading coordinate for vertices or rays).Returns
Cycle<Addition> The shifted cycle - skeleton_complex (C, k, preserveRays) → Cycle<Addition>
Takes a polyhedral complex and computes the k-skeleton. Will return an empty cycle, if k is larger then the dimension of the given complex or smaller than 0.
Parameters
Cycle<Addition> C A polyhedral complex.Int k The dimension of the skeleton that should be computedBool preserveRays When true, the function assumes that all rays of the fan remain in the k-skeleton, so it just copies the VERTICES, instead of computing an irredundant list. By default, this property is false.Returns
Cycle<Addition> The k-skeleton (without any weights, except if k is the dimension of C - triangulate_cycle (F) → Cycle<Addition>
Takes a cycle and computes a triangulation
Parameters
Cycle<Addition> F A cycle (not necessarily weighted)Returns
Cycle<Addition> A simplicial refinement of F
These functions deal with the conversion of tropical objects between Min and Max.
- dual_addition_version (cone, strong_conversion) → Cone
This function takes a tropical cone and returns a tropical cone that uses the opposite tropical addition. By default, the signs of the POINTS are inverted.
Parameters
Cone<Addition,Scalar> cone Bool strong_conversion This is optional and TRUE by default. It indicates, whether the signs of the vertices should be inverted.Returns
Cone - dual_addition_version (cycle, strong_conversion) → Cycle
This function takes a tropical cycle and returns a tropical cycle that uses the opposite tropical addition. By default, the signs of the vertices are inverted.
Parameters
Cycle<Addition> cycle Bool strong_conversion This is optional and TRUE by default. It indicates, whether the signs of the vertices should be inverted.Returns
Cycle - dual_addition_version (number, strong_conversion) → TropicalNumber
This function takes a tropical number and returns a tropical number that uses the opposite tropical addition. By default, the sign is inverted.
Parameters
TropicalNumber<Addition,Scalar> number Bool strong_conversion This is optional and TRUE by default. It indicates, whether the sign of the number should be inverted.Returns
TropicalNumber - dual_addition_version (vector, strong_conversion) → Vector<TropicalNumber>
This function takes a vector of tropical numbers and returns a vector that uses the opposite tropical addition. By default, the signs of the entries are inverted.
Parameters
Vector<TropicalNumber<Addition,Scalar> > vector Bool strong_conversion This is optional and TRUE by default. It indicates, whether the signs of the entries should be inverted.Returns
Vector<TropicalNumber> - dual_addition_version (matrix, strong_conversion) → Matrix<TropicalNumber>
This function takes a matrix of tropical numbers and returns a matrix that uses the opposite tropical addition. By default, the signs of the entries are inverted.
Parameters
Matrix<TropicalNumber<Addition,Scalar> > matrix Bool strong_conversion This is optional and TRUE by default. It indicates, whether the signs of the entries should be inverted.Returns
Matrix<TropicalNumber> - dual_addition_version (ring) → Ring<TropicalNumber>
This function takes a ring over the tropical numbers and returns a ring that uses the opposite tropical addition. Variable names are preserved
- dual_addition_version (polynomial, strong_conversion) → Polynomial<TropicalNumber>
This function takes a tropical polynomial and returns a tropical polynomial that uses the opposite tropical addition. By default, the signs of the coefficients are inverted.
Parameters
Polynomial<TropicalNumber<Addition,Scalar> > polynomial Bool strong_conversion This is optional and TRUE by default. It indicates, whether the signs of the coefficients should be inverted.Returns
Polynomial<TropicalNumber>
These functions create specific morphisms and functions.
Contained in extensionbundled:atint
.- projection_map <Addition> (n, s) → Morphism<Addition>
This creates a linear projection from the projective torus of dimension n to a given set of coordinates.
Type Parameters
Addition Min or MaxParameters
Int n The dimension of the projective torus which is the domain of the projection.Set<Int> s The set of coordinaes to which the map should project. Should be a subset of (0,..,n)Returns
Morphism<Addition> The projection map. - projection_map (n, m) → Morphism
These functions are special +++ for creating special tropical cycles.
Contained in extensionbundled:atint
.- affine_linear_space <Addition> (lineality, translate, weight) → Cycle<Addition>
This creates a true affine linear space.
Type Parameters
Addition Min or MaxParameters
Matrix<Rational> lineality (Row) generators of the lineality space, in tropical homogeneous coordinates, but without the leading zeroVector<Rational> translate Optional. The vertex of the space. By default this is the originInteger weight Optional. The weight of the space. By default, this is 1.Returns
Cycle<Addition> - cross_variety <Addition> (n, k, h, weight) → Cycle<Addition>
This creates the k-skeleton of the tropical variety dual to the cross polytope
Type Parameters
Addition Min or MaxParameters
Int n The (projective) ambient dimensionInt k The (projective) dimension of the variety.Rational h Optional, 1 by default. It is a nonnegative number, describing the height of the one interior lattice point of the cross polytope.Integer weight Optional, 1 by default. The (global) weight of the varietyReturns
Cycle<Addition> The k-skeleton of the tropical hypersurface dual to the cross polytope. It is a smooth (for weight 1), irreducible (for h > 0) variety, which is invariant under reflection. - empty_cycle <Addition> (ambient_dim) → Cycle
Creates the empty cycle in a given ambient dimension (i.e. it will set the property PROJECTIVE_AMBIENT_DIM.
- halfspace_subdivision <Addition> (a, g, The) → Cycle
Creates a subdivision of the tropical projective torus along an affine hyperplane into two halfspaces. This hyperplane is defined by an equation gx = a
Type Parameters
Addition Max or MinParameters
Rational a The constant coefficient of the equationVector<Rational> g The linear coefficients of the equation Note that the equation must be homogeneous in the sense that (1,..1) is in its kernel, i.e. all entries of g add up to 0.Integer The (constant) weight this cycle should haveReturns
Cycle The halfspace subdivision - orthant_subdivision <Addition> (point, chart, weight)
Creates the orthant subdivision around a given point on a given chart, i.e. the corresponding affine chart of this cycle consists of all 2^n fulldimensional orthants
Type Parameters
Addition Min or MaxParameters
Vector<Rational> point The vertex of the subdivision. Should be given in tropical homogeneous coordinates with leading coordinate.Int chart On which chart the cones should be orthants, 0 by default.Integer weight The constant weight of the cycle, 1 by default. - point_collection <Addition> (points, weights) → Cycle
Creates a cycle consisting of a collection of points with given weights
Type Parameters
Addition Max or MinParameters
Matrix<Rational> points The points, in tropical homogeneous coordinates (though not with leading ones for vertices).Vector<Integer> weights The list of weights for the pointsReturns
Cycle The point collection. - projective_torus <Addition> (n, w) → Cycle
Creates the tropical projective torus of a given dimension. In less fancy words, the cycle is the complete complex of given (tropical projective) dimension n, i.e. Rn
- uniform_linear_space <Addition> (n, k, weight) → Cycle
Creates the linear space of the uniform matroid of rank k+1 on n+1 variables.
These functions test cycles for degeneracy, i.e. whether a cycle is the empty cycle
Contained in extensionbundled:atint
.- is_empty ()
This tests wheter a cycle is the empty cycle.
These functions deal with the computation of divisors
Contained in extensionbundled:atint
.- divisor (C, F) → Cycle
This function computes the divisor of one or more rational functions on a tropical cycle.
Parameters
Cycle C A tropical cycleRationalFunction F An arbitrary list of rational functions (r_1,...r_n). The DOMAIN of r_i should contain the support of r_{i-1} * ... * r_1 * C. Note that using the ^-operator on these rational functions is allowed and will result in applying the corresponding function several times.Returns
Cycle The divisor r_n * ... * r_1 * C - divisor_nr (C, F) → Cycle
Parameters
Cycle C A tropical cycleRationalFunction F An arbitrary list of rational functions (r_1,...r_n). The DOMAIN of each function should be equal (in terms of VERTICES and MAXIMAL_POLYTOPES) to the cycle. Note that using the ^-operator on these rational functions is allowed and will result in applying the corresponding function several times.Returns
Cycle The divisor r_n * ... * r_1 * C - piecewise_divisor (F, cones, coefficients) → Cycle<Addition>
Computes a divisor of a linear sum of certain piecewise polynomials on a simplicial fan.
Parameters
Cycle<Addition> F A simplicial fan without lineality space in non-homog. coordinatesIncidenceMatrix cones A list of cones of F (not maximal, but all of the same dimension). Each cone t corresponds to a piecewise polynomial psi_t, defined by subsequently applying the rational functions that are 1 one exactly one ray of t and 0 elsewhere. Note that cones should refer to indices in SEPARATED_VERTICES, which may have a different orderVector<Integer> coefficients A list of coefficients a_t corresponding to the cones.Returns
Cycle<Addition> The divisor sum_t a_t psi_t * F
These functions are wrappers for gfan functions.
- gfan_tropicalbruteforce (I) → fan::SymmetricFan
Calls gfan_tropicalbruteforce for a homogeneous ideal. If the ideal contains a monomial, gfan will return an empty object and the xslt parsing fails. We do not catch this for you.
- gfan_tropicalhypersurface (p) → Cycle<Max>
Calls gfan_tropicalhypersurface for a single polynomial. If the polynomial is a monomial, gfan will return an empty object and the xslt parsing fails. We do not catch this for you.
- gfan_tropicalintersection (I) → Cycle<Max>
Calls gfan_tropicalintersection for a homogeneous ideal.
- gfan_tropicalvariety_of_prime (I) → Cycle<Max>
Calls gfan_tropicalstartingcone | gfan_tropicaltraverse for a homogeneous prime ideal. If the ideal contains a monomial, gfan will return an empty object and the xslt parsing fails. We do not catch this for you.
These functions deal with the creation and study of tropical Hurwitz cycles.
Contained in extensionbundled:atint
.- hurwitz_cycle <Addition> (k, degree, points) → Cycle<Addition>
This function computes the Hurwitz cycle H_k(x), x = (x_1,...,x_n)
Type Parameters
Addition Min or Max, where the coordinates live.Parameters
Int k The dimension of the Hurwitz cycle, i.e. the number of moving verticesVector<Int> degree The degree x. Should add up to 0Vector<Rational> points Optional. Should have length n-3-k. Gives the images of the fixed vertices (besides 0). If not given all fixed vertices are mapped to 0 and the function computes the recession fan of H_k(x)Options
Bool Verbose If true, the function outputs some progress information. True by default.Returns
Cycle<Addition> H_k(x), in homogeneous coordinates - hurwitz_marked_cycle <Addition> (k, degree, pullback_points) → Cycle<Addition>
Computes the marked k-dimensional tropical Hurwitz cycle H_k(degree)
Type Parameters
Addition Min or MaxParameters
Int k The dimension of the Hurwitz cycleVector<Int> degree The degree of the covering. The sum over all entries should be 0 and if n := degree.dim, then 0 <= k <= n-3Vector<Rational> pullback_points The points p_i that should be pulled back to determine the Hurwitz cycle (in addition to 0). Should have length n-3-k. If it is not given, all p_i are by default equal to 0 (same for missing points)Returns
Cycle<Addition> The marked Hurwitz cycle H~_k(degree) - hurwitz_pair <Addition> (k, degree, points) → List
This function computes hurwitz_subdivision and hurwitz_cycle at the same time, returning the result in an array
Type Parameters
Addition Min or Max, where the coordinates live.Parameters
Int k The dimension of the Hurwitz cycle, i.e. the number of moving verticesVector<Int> degree The degree x. Should add up to 0Vector<Rational> points Optional. Should have length n-3-k. Gives the images of the fixed vertices (besides 0). If not given all fixed vertices are mapped to 0 and the function computes the subdivision of M_0,n containing the recession fan of H_k(x)Options
Bool Verbose If true, the function outputs some progress information. True by default.Returns
List ( Cycle subdivision of M_0,n, Cycle Hurwitz cycle ) - hurwitz_pair_local <Addition> (k, degree, local_curve)
Does the same as hurwitz_pair, except that no points are given and the user can give a RationalCurve object representing a ray. If given, the computation will be performed locally around the ray.
Type Parameters
Addition Min or Max, where the coordinates live.Parameters
Int k Vector<Int> degree RationalCurve local_curve Options
Bool Verbose If true, the function outputs some progress information. True by default. - hurwitz_subdivision <Addition> (k, degree, points) → Cycle
This function computes a subdivision of M_0,n containing the Hurwitz cycle H_k(x), x = (x_1,...,x_n) as a subfan. If k = n-4, this subdivision is the unique coarsest subdivision fulfilling this property
Type Parameters
Addition Min or Max, where the coordinates live.Parameters
Int k The dimension of the Hurwitz cycle, i.e. the number of moving verticesVector<Int> degree The degree x. Should add up to 0Vector<Rational> points Optional. Should have length n-3-k. Gives the images of the fixed vertices (besides the first one, which always goes to 0) as elements of R. If not given, all fixed vertices are mapped to 0 and the function computes the subdivision of M_0,n containing the recession fan of H_k(x)Options
Bool Verbose If true, the function outputs some progress information. True by default.Returns
Cycle A subdivision of M_0,n
These are general functions related to intersection theory.
Contained in extensionbundled:atint
.- intersect (X, Y) → Cycle
Computes the intersection product of two tropical cycles in the projective torus Use intersect_check_transversality to check for transversal intersections
- intersect_check_transversality (X, Y, ensure_transversality) → List
Computes the intersection product of two tropical cycles in R^n and tests whether the intersection is transversal (in the sense that the cycles intersect set-theoretically in the right dimension).
Parameters
Cycle X A tropical cycleCycle Y A tropical cycle, living in the same space as XBool ensure_transversality Whether non-transversal intersections should not be computed. Optional and false by default. If true, returns the zero cycle if it detects a non-transversal intersectionReturns
List ( Cycle intersection product, Bool is_transversal). Intersection product is a zero cycle if ensure_transversality is true and the intersection is not transversal. is_transversal is false if the codimensions of the varieties add up to more than the ambient dimension. - intersect_in_smooth_surface (surface, A, B) → Cycle<Addition>
Computes the intersection product of two cycles in a smooth surface
Parameters
Cycle<Addition> surface A smooth surfaceCycle<Addition> A any cycle in the surfaceCycle<Addition> B any cycle in the surfaceReturns
Cycle<Addition> The intersection product of A and B in the surface - point_functions <Addition> (A) → RationalFunction
Constructs a list of rational functions that cut out a single point in the projective torus
Type Parameters
Addition Min or Max. Determines the type of the rational functions.Parameters
Vector<Rational> A point in the projective torus, given in tropical homogeneous coordinates, but without leading coordinate.Returns
RationalFunction . A perl array of rational functions of the form (v_i*x_0 + x_i)/(x_0), i = 1,..,n - pullback (m, r) → RationalFunction
This computes the pullback of a rational function via a morphism Due to the implementation of composition of maps, the DOMAIN of the rational function need not be contained in the image of the morphism The pullback will be defined in the preimage of the domain.
Parameters
Morphism m A morphism.RationalFunction r A rational function.Returns
RationalFunction The pullback m*r.
These functions deal with finding rational functions to given divisors.
Contained in extensionbundled:atint
.- cutting_functions (F, weight_aim) → Matrix<Rational>
Takes a weighted complex and a list of desired weights on its codimension one faces and computes all possible rational functions on (this subdivision of ) the complex
Parameters
Cycle<Addition> F A tropical variety, assumed to be simplicial.Vector<Integer> weight_aim A list of weights, whose length should be equal to the number of CODIMENSION_ONE_POLYTOPES. Gives the desired weight on each codimension one faceReturns
Matrix<Rational> The space of rational functions defined on this particular subdivision. Each row is a generator. The columns correspond to values on SEPARATED_VERTICES and LINEALITY_SPACE, except the last one, which is either 0 (then this function cuts out zero and can be added to any solution) or non-zero (then normalizing this entry to -1 gives a function cutting out the desired weights on the codimension one skeleton Note that the function does not test if these generators actually define piecewise linear functions, as it assumes the cycle is simplicial - simplicial_diagonal_system (fan) → Matrix<Rational>
This function computes the inhomogeneous version of simplicial_piecewise_system in the sense that it computes the result of the above mentioned function (i.e. which coefficients for the piecewise polynomials yield the zero divisor) and adds another column at the end where only the entries corresponding to the diagonal cones are 1, the rest is zero. This can be seen as asking for a solution to the system that cuts out the diagonal (all solutions whose last entry is 1)
- simplicial_piecewise_system (F) → Matrix<Rational>
This function takes a d-dimensional simplicial fan F and computes the linear system defined in the following way: For each d-dimensional cone t in the diagonal subdivision of FxF, let psi_t be the piecewise polynomial defined by subsequently applying the rational functions that are 1 one exactly one ray of t and 0 elsewhere. Now for which coefficients a_t is sum_t a_t psi_t * (FxF) = 0?
Parameters
Cycle<Addition> F A simplicial fan without lineality spaceReturns
Matrix<Rational> The above mentioned linear system. The rows are equations, the columns correspond to d-dimensional cones of FxF in the order given by skeleton_complex(simplicial_with_diagonal(F), d, 1) - simplicial_with_diagonal (F) → Cycle<Addition>
This function takes a simplicial fan F (without lineality space) and computes the coarsest subdivision of F x F containing all diagonal rays (r,r)
Parameters
Cycle<Addition> F A simplicial fan without lineality space.Returns
Cycle<Addition> The product complex FxF subdivided such that it contains all diagonal rays
These functions deal with lattices (meaning free abelian, finitely generated groups).
Contained in extensionbundled:atint
.- lattice_index (m) → Integer
This computes the index of a lattice in its saturation.
Parameters
Matrix<Integer> m A list of (row) generators of the lattice.Returns
Integer The index of the lattice in its saturation. - randomInteger (max_arg, n) → Array<Integer>
Returns n random integers in the range 0.. (max_arg-1),inclusive Note that this algorithm is not optimal for real randomness: If you change the range parameter and then change it back, you will usually get the exact same sequence as the first time
Parameters
Int max_arg The upper bound for the random integersInt n The number of integers to be createdReturns
Array<Integer>
These functions deal with the computation and representation of (families of) lines in surfaces.
Contained in extensionbundled:atint
.- lines_in_cubic (p) → LinesInCubic<Addition>
This takes either: - A homogeneous polynomial of degree 3 in 4 variables or - A polynomial of degree 3 in 3 variables and computes the corresponding cubic and finds all tropical lines and families thereof in the cubic. The result is returned as a LinesInCubic object. Note that the function has some heuristics for recognizing families, but might still return a single family as split up into two.
Parameters
Polynomial<TropicalNumber<Addition> > p A homogeneous tropical polynomial of degree 3 in four variables.Returns
LinesInCubic<Addition>
These functions are used for doing computations locally around a specified part of a Cycle. ----- These +++ deal with the creation and modification of cycles with nontrivial LOCAL_RESTRICTION.
Contained in extensionbundled:atint
.- local_codim_one (complex, face) → Cycle<Addition>
This takes a weighted complex and an index of one of its codimension one faces (The index is in CODIMENSION_ONE_POLYTOPES) and computes the complex locally restricted to that face
Parameters
Cycle<Addition> complex An arbitrary weighted complexInt face An index of a face in CODIMENSION_ONE_POLYTOPESReturns
Cycle<Addition> The complex locally restricted to the given face - local_point (complex, v) → Cycle<Addition>
This takes a weighted complex and an arbitrary vertex in homogeneous coordinates (including the leading coordinate) that is supposed to lie in the support of the complex. It then refines the complex such that the vertex is a cell in the polyhedral structure and returns the complex localized at this vertex
Parameters
Cycle<Addition> complex An arbitrary weighted complexVector<Rational> v A vertex in homogeneous coordinates and with leading coordinate. It should lie in the support of the complex (otherwise an error is thrown)Returns
Cycle<Addition> The complex localized at the vertex - local_restrict (complex, cones) → Cycle<Addition>
This takes a tropical variety and an IncidenceMatrix describing a set of cones (not necessarily maximal ones) of this variety. It will then create a variety that contains all compatible maximal cones and is locally restricted to the given cone set.
Parameters
Cycle<Addition> complex An arbitrary weighted complexIncidenceMatrix cones A set of cones, indices refer to VERTICESReturns
Cycle<Addition> The same complex, locally restricted to the given cones - local_vertex (complex, ray) → Cycle<Addition>
This takes a weighted complex and an index of one of its vertices (the index is to be understood in VERTICES) It then localizes the variety at this vertex. The index should never correspond to a far vertex in a complex, since this would not be a cone
Parameters
Cycle<Addition> complex An arbitrary weighted complexInt ray The index of a ray/vertex in RAYSReturns
Cycle<Addition> The complex locally restricted to the given vertex
These functions deal with matroids and matroidal fans.
Contained in extensionbundled:atint
.- is_smooth (a) → List
Takes a weighted fan and returns if it is smooth (i.e. isomorphic to a Bergman fan B(M)/L for some matroid M) or not. The algorithm works for fans of dimension 0,1,2 and codimension 0,1! For other dimensions the algorithm could give an answer but it is not guaranteed.
Parameters
Cycle<Addition> a tropical fan FReturns
List ( Int s, Matroid M, Morphism<Addition> A ). If s=1 then F is smooth, the corresponding matroid fan is Z-isomorphic to the matroid fan associated to M. The Z-isomorphism is given by A, i.e. B(M)/L = affine_transform(F,A) If s=0, F is not smooth. If s=2 the algorithm is not able to determine if F is smooth or not. - matroid_fan <Addition> (m) → Cycle
Uses an algorithm by Felipe Rincón to compute the matroidal fan of a given matroid. If you have a matrix at hand that represents this matroid, it is recommended to call this function with that matrix as an argument - it is significantly faster.
Type Parameters
Addition Min or Max - determines the coordinates.Parameters
matroid::Matroid m A matroidReturns
Cycle The matroidal fan or Bergman fan of the matroid. - matroid_fan <Addition> (m) → Cycle
Uses an algorithm by Felipe Rincón to compute the bergman fan of the column matroid of the given matrix. Calling the function in this manner is significantly faster than calling it on the matroid.
Type Parameters
Addition Min or Max - determines the coordinates.Parameters
Matrix<Rational> m A matrix, whose column matroid is considered.Returns
Cycle The matroidal fan or Bergman fan of the matroid. - matroid_fan_from_flats <Addition> (A) → Cycle<Addition>
Computes the fan of a matroid in its chains-of-flats subdivision. Note that this is potentially very slow for large matroids.
Type Parameters
Addition Min or max, determines the matroid fan coordinates.Parameters
matroid::Matroid A matroid. Should be loopfree.Returns
Cycle<Addition> - matroid_from_fan (A) → matroid::Matroid
Takes the bergman fan of a matroid and reconstructs the corresponding matroid The fan has to be given in its actual matroid coordinates, not as an isomorphic transform. The actual subdivision is not relevant.
These functions deal with moduli spaces of abstract or parametrized rational curves.
Contained in extensionbundled:atint
.- count_mn_cones (n, k) → Integer
- count_mn_rays (n) → Integer
- evaluation_map <Addition> (n, Delta, i) → Morphism<Addition>
This creates the i-th evaluation function on M_0,n^(lab)(R^r,Delta) (which is actually realized as M_0,(n+|Delta|) x R^r and can be created via space_of_stable_maps).
Type Parameters
Addition Min or MaxParameters
Int n The number of marked (contracted) pointsMatrix<Rational> Delta The directions of the unbounded edges (given as row vectors in tropical projective coordinates without leading coordinate, i.e. have r+1 columns)Int i The index of the marked point that should be evaluated. Should lie in between 1 and n Note that the i-th marked point is realized as the |Delta|+i-th leaf in M_0,(n+|Delta|) and that the R^r - coordinate is interpreted as the position of the n-th leaf. In particular, ev_n is just the projection to the R^r-coordinatesReturns
Morphism<Addition> ev_i. Its domain is the ambient space of the moduli space as created by space_of_stable_maps. The target space is the tropical projective torus of dimension r - evaluation_map <Addition> (n, r, d, i) → Morphism<Addition>
This creates the i-th evaluation function on M_0,n^(lab)(R^r,d) (which is actually realized as M_0,(n+d(r+1)) x R^r) This is the same as calling the function evaluation_map(Int,Int,Matrix<Rational>,Int) with the standard d-fold degree as matrix (i.e. each (inverted) unit vector of R^(r+1) occuring d times).
Type Parameters
Addition Min or MaxParameters
Int n The number of marked (contracted) pointsInt r The dimension of the target spaceInt d The degree of the embedding. The direction matrix will be the standard d-fold directions, i.e. each unit vector (inverted for Max), occuring d times.Int i The index of the marked point that should be evaluated. i should lie in between 1 and nReturns
Morphism<Addition> ev_i. Its domain is the ambient space of the moduli space as created by space_of_stable_maps. The target space is the tropical projective torus of dimension r - forgetful_map <Addition> (n, S) → Morphism
This computes the forgetful map from the moduli space M_0,n to M_0,(n-|S|)
Type Parameters
Addition Min or MaxParameters
Int n The number of leaves in the moduli space M_0,nSet<Int> S The set of leaves to be forgotten. Should be a subset of (1,..,n)Returns
Morphism The forgetful map. It will identify the remaining leaves i_1,..,i_(n-|S|) with the leaves of M_0,(n-|S|) in canonical order. The domain of the morphism is the ambient space of the morphism in matroid coordinates, as created by m0n. - local_m0n <Addition> (R ...) → Cycle<Addition>
Computes the moduli space M_0,n locally around a given list of combinatorial types. More precisely: It computes the weighted complex consisting of all maximal cones containing any of the given combinatorial types and localizes at these types This should only be used for curves of small codimension. What the function actually does, is that it combinatorially computes the cartesian products of M_0,v's, where v runs over the possible valences of vertices in the curves For max(v) <= 8 this should terminate in a reasonable time (depending on the number of curves) The coordinates are the same that would be produced by the function m0n
Type Parameters
Addition Min or Max, determines the coordinatesParameters
RationalCurve R ... A list of rational curves (preferrably in the same M_0,n)Returns
Cycle<Addition> The local complex - m0n <Addition> (n) → Cycle
Creates the moduli space of abstract rational n-marked curves. Its coordinates are given as the coordinates of the bergman fan of the matroid of the complete graph on n-1 nodes (but not computed as such) The isomorphism to the space of curve metrics is obtained by choosing the last leaf as special leaf
- psi_class <Addition> (n, i) → Cycle
- psi_product <Addition> (n, exponents) → Cycle
Computes a product of psi classes psi_1^k_1 * ... * psi_n^k_n on the moduli space of rational n-marked tropical curves M_0,n
Type Parameters
Addition Min or MaxParameters
Int n The number of leaves in M_0,nVector<Int> exponents The exponents of the psi classes k_1,..,k_n. If the vector does not have length n or if some entries are negative, an error is thrownReturns
Cycle The corresponding psi class divisor - space_of_stable_maps <Addition> (n, d, r) → Cycle
Creates the moduli space of stable maps of rational n-marked curves into a projective torus. It is given as the cartesian product of M_{0,n+d} and R^r, where n is the number of contracted leaves, d the number of non-contracted leaves and r is the dimension of the target torus. The R^r - coordinate is interpreted as the image of the last (n-th) contracted leaf. Due to the implementation of cartesian_product, the projective coordinates are non-canonical: Both M_{0,n+d} and R^r are dehomogenized after the first coordinate, then the product is taken and homogenized after the first coordinate again. Note that functions in a-tint will usually treat this space in such a way that the first d leaves are the non-contracted ones and the remaining n leaves are the contracted ones.
These are general functions that deal with morphisms and their arithmetic.
Contained in extensionbundled:atint
.- add_morphisms (f, g) → Morphism
Special purpose functions.
- lifted_pluecker (V) → Vector<TropicalNumber<Addition> >
Compute the tropical Pluecker vector from a matrix representing points in the tropical torus. This can be used to lift regular subdivisions of a product of simplices to a matroid decomposition of hypersimplices.
- points_in_pseudovertices (points, pseudovertices) → Set<Int>
This function takes a Matrix of tropical vectors in projective coordinates (e.g. the POINTS of a Cone) and a Matrix of Scalar vectors in extended tropical projective coordinates (e.g. the PSEUDOVERTICES of a tropical Cone). It returns the set of row indices of the second matrix such that the corresponding row starts with a 1 and the remaining vector occurs in the first matrix.
These functions produce a tropical hypersurface from other objects.
- hyperplane <Addition> (coeffs) → Hypersurface<Addition>
Create a tropical hyperplane as object of type Hypersurface.
Type Parameters
Addition Parameters
Vector<TropicalNumber<Addition> > coeffs coefficients of the tropical linear form (can also be specified as anonymous array).Returns
Hypersurface<Addition> - points2hypersurface (points) → Hypersurface
Constructs a tropical hypersurface defined by the linear hypersurfaces associated to the points. If the points are min-tropical points then the output is a max-tropical hypersurface, and conversely.
These functions produce an object of type Cone from other objects.
- cyclic <Addition> (d, n) → Cone<Addition>
Produces a tropical cyclic d-polytope with n vertices. Cf.
Josephine Yu & Florian Block, arXiv: math.MG/0503279.Type Parameters
Addition Min or Max.Parameters
Int d the dimensionInt n the number of generatorsReturns
Cone<Addition> - hypersimplex <Addition> (d, k) → Cone<Addition>
Produce the tropical hypersimplex Δ(k,d). Cf.
M. Joswig math/0312068v3, Ex. 2.10.The value of k defaults to 1, yielding a tropical standard simplex.
Type Parameters
Addition Max or MinParameters
Int d the dimensionInt k the number of +/-1 entriesReturns
Cone<Addition> - matroid_polytope <Addition, Scalar> (m, v) → Cone<Addition,Scalar>
Produce the tropical matroid polytope from a matroid m. Each vertex corresponds to a basis of the matroid, the non-bases coordinates get value 0, the bases coordinates get value v, default is -orientation.
Type Parameters
Addition Min or MaxScalar coordinate typeParameters
matroid::Matroid m Scalar v value for the basesReturns
Cone<Addition,Scalar> - minkowski_sum (lambda, P, mu, Q) → Cone<Addition,Scalar>
Produces the tropical polytope (lambda \( \otimes \) P) \( \oplus \) (mu \( \otimes \) Q), where \( \otimes \) and \( \oplus \) are tropical scalar multiplication and tropical addition, respectively.
These functions deal with covectors of subdivision of tropical point configurations.
- coarse_covectors (points, generators) → Matrix<int>
This computes the coarse covector of a list of points relative to a list of generators.
Parameters
Matrix<TropicalNumber<Addition,Scalar> > points Matrix<TropicalNumber<Addition,Scalar> > generators Returns
Matrix<int> . Rows correspond to points, columns to coordinates. Each entry encodes, how many generators contain a given point in a certain coordinate. - coarse_covectors_of_scalar_vertices (points, maximal_cells, generators) → Matrix<int>
Computes the coarse covectors of a list of scalar points, as described in covectors_of_scalar_vertices
Parameters
Matrix<Scalar> points IncidenceMatrix maximal_cells Matrix<TropicalNumber<Addition,Scalar> > generators Returns
Matrix<int> . Rows correspond to points, columns to coordinates. Each entry encodes, how many generators contain a given point in a certain coordinate. - covectors (points, generators) → Array<IncidenceMatrix>
This computes the (fine) covector of a list of points relative to a list of generators.
Parameters
Matrix<TropicalNumber<Addition,Scalar> > points Matrix<TropicalNumber<Addition,Scalar> > generators Returns
Array<IncidenceMatrix> . Each IncidenceMatrix corresponds to a point. Rows of a matrix correspond to coordinates and columns to generators. Each row indicates, which generators contain the point in the sector corresponding to the coordinate. - covectors_of_scalar_vertices (points, maximal_cells, generators) → Array<IncidenceMatrix>
This computes the (fine) covector of a list of points relative to a list of generators. The points are scalar points and they are supposed to be normalized in the following sense: - All bounded vertices have a leading 1 - All unbounded vertices have a leading 0 and all nonzero entries are either +1 or -1. (but not both) Furthermore, the points make up a polyhedral complex - in particular, every maximal cell has a bounded vertex. For the bounded vertices, covectors are computed as usual. For unbounded vertices, the nonzero entries are replaced by tropical zero, the complementary entries are copied from a bounded vertex sharing a cell and then the covector is computed.
Parameters
Matrix<Scalar> points IncidenceMatrix maximal_cells Matrix<TropicalNumber<Addition,Scalar> > generators Returns
Array<IncidenceMatrix> . Each IncidenceMatrix corresponds to a point. Rows of a matrix correspond to coordinates and columns to generators. Each row indicates, which generators contain the point in the sector corresponding to the coordinate.
These functions deal with general tropical arithmetic.
- nearest_point (C, x) → Vector<TropicalNumber<Addition,Scalar> >
Compute the projection of a point x in tropical projective space onto a tropical cone C. Cf.
Develin & Sturmfels math.MG/0308254v2, Proposition 9. - norm (v) → Scalar
The tropical norm of a vector v in the tropical torus is the difference between the maximal and minimal coordinate in any coordinate representation of the vector.
- principal_solution (A, b) → Vector<TropicalNumber<Addition, Scalar> >
Compute the solution of the tropical equation A * x = b. If there is no solution, the return value is 'near' a solution. Cf. Butkovič 'Max-linear systems: theory and algorithms' (MR2681232), Theorem 3.1.1
- tdet (matrix) → TropicalNumber<Addition,Scalar>
The tropical determinant of a matrix.
These functions are for visualization.
Contained in extensionbundled:atint
.- visualize_in_surface ()
This visualizes a surface in R^3 and an arbitrary list of (possibly non-pure) Cycle objects. A common bounding box is computed for all objects and a random color is chosen for each object (except the surface)
These functions deal with the weight space of a cycle, i.e. the space of weights that make it balanced and related properties.
Contained in extensionbundled:atint
.- decomposition_polytope (A) → polytope::Polytope
Computes the possible positive decompositions into irreducible subvarieties of the same weight positivity signature (i.e. the weight on a cone has to have the same sign as in the cycle) To be precise, it computes the irreducible varieties as rays of the weight cone (where the corresponding orthant is taken such that the weight vector of X lies in that orthant). It then computes the polytope of all positive linear combinations of those irreducible varieties that produce the original weight vector.
- weight_cone (X, negative)
Takes a polyhedral complex and computes a weight cone, i.e. intersects the WEIGHT_SPACE with a chosen orthant (by default the positive orthant)
Parameters
Cycle X A polyhedral complexSet<int> negative A subset of the coordinates {0,..,N-1}, where N is the number of maximal cells of X. Determines the orthant to intersect the weight space with: All coordinates in the set are negative, the others positive If the set is not given, it is empty by default (i.e. we take the positive orthant)
These functions relate to the weights of a tropical cycle.
- is_balanced (C)
This computes whether a given cycle is balanced.
Parameters
Cycle C The cycle for which to check balancing.$ @return Bool Whether the cycle is balanced.
Common Option Lists
- Visual::CovectorLattice::decorationsimports from: Visual::Lattice::decorations, Visual::Graph::decorations, Visual::Wire::decorations, Visual::PointSet::decorations
UNDOCUMENTED
Options
String Covectors if set to "hidden", the covectors are not displayed. - Visual::Cycle::BoundingDecorationsimports from: Visual::Polygons::decorations, Visual::PointSet::decorations
UNDOCUMENTED
Options
Flexible<Int> Chart The visualization is always a visualization in affine coordinates. This chooses the chartString WeightLabels if set to "hidden", the labels indicating the weights are hiddenString CoordLabels If set to "show", the labels indicating the vertex coordinates are displayed, otherwise they are not. Note that this is expensive and significantly increases computation time.Flexible<Rational> BoundingDistance The distance of the border of the bounding box from the smallest box containing the affine points of the complex. This is only relevant, if BoundingMode is "relative"Matrix<Rational> BoundingBox A fixed bounding box, determined by two row vectors that specify two of its vertices (on "on top" and one "at the bottom"). Is only relevant, if BoundingMode is "absolute"String BoundingMode If set to "relative", the function determines the smallest possible box containing all affine points of the complex and then enlarges the box by BoundingDistance to all sides. If set to "absolute", BoundingBox must be specified and the complex will be intersected with that box. "cube" does a similar thing as relative, except that the resulting bounding box is always a cube. By default this is set to "cube".Array<String> ConeLabels A list of strings to be displayed as labels for the maximal cones. If this is empty, the weight labels (if present and not suppressed by WeightLabels=>"hidden") are displayed - Visual::Cycle::decorationsimports from: Visual::Polygons::decorations, Visual::PointSet::decorations
UNDOCUMENTED
Options
Flexible<Int> Chart The visualization is always a visualization in affine coordinates. This chooses the chart - Visual::Cycle::FunctionDecorationsimports from: Visual::Cycle::BoundingDecorations, Visual::Polygons::decorations, Visual::PointSet::decorations
Visualization options for RationalFunction/Morphism->BB_VISUAL
Options
String FunctionLabels , if set to "show", the function labels indicatingt the affine linear representation of each function on each cone are computed - Visual::RationalCurve::decorationsimports from: Visual::Graph::decorations, Visual::Wire::decorations, Visual::PointSet::decorations
UNDOCUMENTED
Options
String LengthLabels if set to "hidden", the labels indicating the lengths are hidden - Visual::TropicalCone::decorationsimports from: Visual::Polygons::decorations, Visual::PointSet::decorations
UNDOCUMENTED
Options
Flexible<Int> Chart The visualization is always in affine coordinates. This chooses the chart.Flexible<Color> GeneratorColor color of the finite POINTS of a coneFlexible<Color> PseudovertexColor color of the finite PSEUDOVERTICES of a cone