Constructions

The standard way to define a polyhedron is by either giving a $V$-representation or an $H$-representation. But polyhedra may also be constructed through other means: by name, via operations on other polyhedra, or from other objects in OSCAR.

$H$- and $V$-representations

Intersecting halfspaces: $H$-representation

polyhedronMethod
polyhedron([::Union{Type{T}, Field},] A::AnyVecOrMat, b) where T<:scalar_types

The (convex) polyhedron defined by

\[P(A,b) = \{ x | Ax ≤ b \}.\]

see Def. 3.35 and Section 4.1. of [JT13]

The first argument either specifies the Type of its coefficients or their parent Field.

Examples

The following lines define the square $[0,1]^2 \subset \mathbb{R}^2$:

julia> A = [1 0; 0 1; -1 0 ; 0 -1];

julia> b = [1, 1, 0, 0];

julia> polyhedron(A,b)
Polyhedron in ambient dimension 2
source
polyhedronFunction
polyhedron(::Union{Type{T}, Field}, I::Union{Nothing, AbstractCollection[AffineHalfspace]}, E::Union{Nothing, AbstractCollection[AffineHyperplane]} = nothing) where T<:scalar_types

The (convex) polyhedron obtained intersecting the halfspaces I (inequalities) and the hyperplanes E (equations). The first argument either specifies the Type of its coefficients or their parent Field.

Examples

The following lines define the square $[0,1]^2 \subset \mathbb{R}^2$:

julia> A = [1 0; 0 1; -1 0 ; 0 -1];

julia> b = [1, 1, 0, 0];

julia> polyhedron((A,b))
Polyhedron in ambient dimension 2

As an example for a polyhedron constructed from both inequalities and equations, we construct the polytope $[0,1]\times\{0\}\subset\mathbb{R}^2$

julia> P = polyhedron(([-1 0; 1 0], [0,1]), ([0 1], [0]))
Polyhedron in ambient dimension 2

julia> is_feasible(P)
true

julia> dim(P)
1

julia> vertices(P)
2-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 0]
 [0, 0]
source

The complete $H$-representation can be retrieved using facets and affine_hull:

julia> P = polyhedron(([-1 0; 1 0], [0,1]), ([0 1], [0]))
Polyhedron in ambient dimension 2

julia> facets(P)
2-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^2 described by:
-x_1 <= 0
x_1 <= 1


julia> affine_hull(P)
1-element SubObjectIterator{AffineHyperplane{QQFieldElem}} over the Hyperplanes of R^2 described by:
x_2 = 0


julia> Q0 = polyhedron(facets(P))
Polyhedron in ambient dimension 2

julia> P == Q0
false

julia> Q1 = polyhedron(facets(P), affine_hull(P))
Polyhedron in ambient dimension 2

julia> P == Q1
true

Computing convex hulls: $V$-representation

convex_hullMethod
convex_hull([::Union{Type{T}, Field} = QQFieldElem,] V [, R [, L]]; non_redundant::Bool = false)

Construct the convex hull of the vertices V, rays R, and lineality L. If R or L are omitted, then they are assumed to be zero.

Arguments

  • The first argument either specifies the Type of its coefficients or their

parent Field.

  • V::AbstractCollection[PointVector]: Points whose convex hull is to be computed.
  • R::AbstractCollection[RayVector]: Rays completing the set of points.
  • L::AbstractCollection[RayVector]: Generators of the Lineality space.

If an argument is given as a matrix, its content has to be encoded row-wise.

R can be given as an empty matrix or as nothing if the user wants to compute the convex hull only from V and L.

If it is known that V and R only contain extremal points and that the description of the lineality space is complete, set non_redundant = true to avoid unnecessary redundancy checks.

See Def. 2.11 and Def. 3.1 of [JT13].

Examples

The following lines define the square $[0,1]^2 \subset \mathbb{R}^2$:

julia> Square = convex_hull([0 0; 0 1; 1 0; 1 1])
Polyhedron in ambient dimension 2

To construct the positive orthant, rays have to be passed:

julia> V = [0 0];

julia> R = [1 0; 0 1];

julia> PO = convex_hull(V, R)
Polyhedron in ambient dimension 2

The closed-upper half plane can be constructed by passing rays and a lineality space:

julia> V = [0 0];

julia> R = [0 1];

julia> L = [1 0];

julia> UH = convex_hull(V, R, L)
Polyhedron in ambient dimension 2

To obtain the x-axis in $\mathbb{R}^2$:

julia> V = [0 0];

julia> R = nothing;

julia> L = [1 0];

julia> XA = convex_hull(V, R, L)
Polyhedron in ambient dimension 2
source

This is a standard triangle, defined via a (redundant) $V$-representation and its unique minimal $H$-representation:

julia> T = convex_hull([ 0 0 ; 1 0 ; 0 1; 0 1//2 ])
Polyhedron in ambient dimension 2

julia> halfspace_matrix_pair(facets(T))
(A = [-1 0; 0 -1; 1 1], b = QQFieldElem[0, 0, 1])

The complete $V$-representation can be retrieved using minimal_faces, rays_modulo_lineality and lineality_space:

julia> P = convex_hull([0 0], [1 0], [0 1])
Polyhedron in ambient dimension 2

julia> Q0 = convex_hull(vertices(P))
Polyhedron in ambient dimension 2

julia> P == Q0
false

julia> mfP = minimal_faces(P)
(base_points = PointVector{QQFieldElem}[[0, 0]], lineality_basis = RayVector{QQFieldElem}[[0, 1]])

julia> rmlP = rays_modulo_lineality(P)
(rays_modulo_lineality = RayVector{QQFieldElem}[[1, 0]], lineality_basis = RayVector{QQFieldElem}[[0, 1]])

julia> Q1 = convex_hull(mfP.base_points, rmlP.rays_modulo_lineality)
Polyhedron in ambient dimension 2

julia> P == Q1
false

julia> Q0 == Q1
false

julia> Q2 = convex_hull(mfP.base_points, rmlP.rays_modulo_lineality, lineality_space(P))
Polyhedron in ambient dimension 2

julia> P == Q2
true

Regular polytopes

A polytope is regular, in the strict sense, if it admits a flag-transtive group of (linear) automorphisms. There are three infinite families of regular polytopes which exist in each dimension: the (regular) simplices, cubes and cross polytopes. In addition there are two exceptional regular 3-polytopes (dodecahedron and icosahedron) plus three exceptional regular 4-polytopes (24-cell, 120-cell and 600-cell).

The regular 3-polytopes are also known as the Platonic solids. Here we also list the Archimedean, Catalan and Johnson solids, which form various generalizations of the Platonic solids. However, here we implement "disjoint families", i.e., the proper Archimedean solids exclude the Platonic solids; similarly, the proper Johnson solids exclude the Archmidean solids.

simplexFunction
simplex([::Union{Type{T}, Field} = QQFieldElem,] d::Int [,n])

Construct the simplex which is the convex hull of the standard basis vectors along with the origin in $\mathbb{R}^d$, scaled by $n$. The first argument either specifies the Type of its coefficients or their parent Field.

Examples

Here we take a look at the facets of the 7-simplex and a scaled 7-simplex:

julia> s = simplex(7)
Polytope in ambient dimension 7

julia> facets(s)
8-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^7 described by:
-x_1 <= 0
-x_2 <= 0
-x_3 <= 0
-x_4 <= 0
-x_5 <= 0
-x_6 <= 0
-x_7 <= 0
x_1 + x_2 + x_3 + x_4 + x_5 + x_6 + x_7 <= 1

julia> t = simplex(7, 5)
Polytope in ambient dimension 7

julia> facets(t)
8-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^7 described by:
-x_1 <= 0
-x_2 <= 0
-x_3 <= 0
-x_4 <= 0
-x_5 <= 0
-x_6 <= 0
-x_7 <= 0
x_1 + x_2 + x_3 + x_4 + x_5 + x_6 + x_7 <= 5
source
cross_polytopeFunction
cross_polytope([::Union{Type{T}, Field} = QQFieldElem,] d::Int [,n])

Construct a $d$-dimensional cross polytope around origin with vertices located at $\pm e_i$ for each unit vector $e_i$ of $R^d$, scaled by $n$. The first argument either specifies the Type of its coefficients or their parent Field.

Examples

Here we print the facets of a non-scaled and a scaled 3-dimensional cross polytope:

julia> C = cross_polytope(3)
Polytope in ambient dimension 3

julia> facets(C)
8-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^3 described by:
x_1 + x_2 + x_3 <= 1
-x_1 + x_2 + x_3 <= 1
x_1 - x_2 + x_3 <= 1
-x_1 - x_2 + x_3 <= 1
x_1 + x_2 - x_3 <= 1
-x_1 + x_2 - x_3 <= 1
x_1 - x_2 - x_3 <= 1
-x_1 - x_2 - x_3 <= 1

julia> D = cross_polytope(3, 2)
Polytope in ambient dimension 3

julia> facets(D)
8-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^3 described by:
x_1 + x_2 + x_3 <= 2
-x_1 + x_2 + x_3 <= 2
x_1 - x_2 + x_3 <= 2
-x_1 - x_2 + x_3 <= 2
x_1 + x_2 - x_3 <= 2
-x_1 + x_2 - x_3 <= 2
x_1 - x_2 - x_3 <= 2
-x_1 - x_2 - x_3 <= 2
source
cubeFunction
cube([::Union{Type{T}, Field} = QQFieldElem,] d::Int , [l::Rational = -1, u::Rational = 1])

Construct the $[l,u]$-cube in dimension $d$. The first argument either specifies the Type of its coefficients or their parent Field.

Examples

In this example the 5-dimensional unit cube is constructed to ask for one of its properties:

julia> C = cube(5,0,1);

julia> normalized_volume(C)
120
source
tetrahedronFunction
tetrahedron()

Construct the regular tetrahedron, one of the Platonic solids.

source
dodecahedronFunction
dodecahedron()

Construct the regular dodecahedron, one out of two Platonic solids.

source
icosahedronFunction
icosahedron()

Construct the regular icosahedron, one out of two exceptional Platonic solids.

source
platonic_solidFunction
platonic_solid(s)

Construct a Platonic solid with the name given by String s from the list below.

Arguments

  • s::String: The name of the desired Archimedean solid. Possible values:
    • "tetrahedron" : Tetrahedron. Regular polytope with four triangular facets.
    • "cube" : Cube. Regular polytope with six square facets.
    • "octahedron" : Octahedron. Regular polytope with eight triangular facets.
    • "dodecahedron" : Dodecahedron. Regular polytope with 12 pentagonal facets.
    • "icosahedron" : Icosahedron. Regular polytope with 20 triangular facets.

Examples

julia> T = platonic_solid("icosahedron")
Polytope in ambient dimension 3 with EmbeddedAbsSimpleNumFieldElem type coefficients

julia> n_facets(T)
20
source
archimedean_solidFunction
archimedean_solid(s)

Construct an Archimedean solid with the name given by String s from the list below. The polytopes are realized with floating point numbers and thus not exact; Vertex-facet-incidences are correct in all cases.

Arguments

  • s::String: The name of the desired Archimedean solid. Possible values:
    • "truncated_tetrahedron" : Truncated tetrahedron. Regular polytope with four triangular and four hexagonal facets.
    • "cuboctahedron" : Cuboctahedron. Regular polytope with eight triangular and six square facets.
    • "truncated_cube" : Truncated cube. Regular polytope with eight triangular and six octagonal facets.
    • "truncated_octahedron" : Truncated Octahedron. Regular polytope with six square and eight hexagonal facets.
    • "rhombicuboctahedron" : Rhombicuboctahedron. Regular polytope with eight triangular and 18 square facets.
    • "truncated_cuboctahedron" : Truncated Cuboctahedron. Regular polytope with 12 square, eight hexagonal and six octagonal facets.
    • "snub_cube" : Snub Cube. Regular polytope with 32 triangular and six square facets. The vertices are realized as floating point numbers. This is a chiral polytope.
    • "icosidodecahedron" : Icosidodecahedon. Regular polytope with 20 triangular and 12 pentagonal facets.
    • "truncated_dodecahedron" : Truncated Dodecahedron. Regular polytope with 20 triangular and 12 decagonal facets.
    • "truncated_icosahedron" : Truncated Icosahedron. Regular polytope with 12 pentagonal and 20 hexagonal facets.
    • "rhombicosidodecahedron" : Rhombicosidodecahedron. Regular polytope with 20 triangular, 30 square and 12 pentagonal facets.
    • "truncated_icosidodecahedron" : Truncated Icosidodecahedron. Regular polytope with 30 square, 20 hexagonal and 12 decagonal facets.
    • "snub_dodecahedron" : Snub Dodecahedron. Regular polytope with 80 triangular and 12 pentagonal facets. The vertices are realized as floating point numbers. This is a chiral polytope.

Examples

julia> T = archimedean_solid("cuboctahedron")
Polytope in ambient dimension 3

julia> sum([n_vertices(F) for F in faces(T, 2)] .== 3)
8

julia> sum([n_vertices(F) for F in faces(T, 2)] .== 4)
6

julia> n_facets(T)
14
source
johnson_solidFunction
johnson_solid(i::Int)

Construct the i-th proper Johnson solid.

A Johnson solid is a 3-polytope whose facets are regular polygons, of various gonalities. It is proper if it is not an Archimedean solid. Up to scaling there are exactly 92 proper Johnson solids.

source
catalan_solidFunction
catalan_solid(s::String)

Construct a Catalan solid with the name s from the list below. The polytopes are realized with floating point coordinates and thus are not exact. However, vertex-facet-incidences are correct in all cases.

Arguments

  • s::String: The name of the desired Archimedean solid. Possible values:
    • "triakis_tetrahedron" : Triakis Tetrahedron. Dual polytope to the Truncated Tetrahedron, made of 12 isosceles triangular facets.
    • "triakis_octahedron" : Triakis Octahedron. Dual polytope to the Truncated Cube, made of 24 isosceles triangular facets.
    • "rhombic_dodecahedron" : Rhombic dodecahedron. Dual polytope to the cuboctahedron, made of 12 rhombic facets.
    • "tetrakis_hexahedron" : Tetrakis hexahedron. Dual polytope to the truncated octahedron, made of 24 isosceles triangluar facets.
    • "disdyakis_dodecahedron" : Disdyakis dodecahedron. Dual polytope to the truncated cuboctahedron, made of 48 scalene triangular facets.
    • "pentagonal_icositetrahedron" : Pentagonal Icositetrahedron. Dual polytope to the snub cube, made of 24 irregular pentagonal facets. The vertices are realized as floating point numbers.
    • "pentagonal_hexecontahedron" : Pentagonal Hexecontahedron. Dual polytope to the snub dodecahedron, made of 60 irregular pentagonal facets. The vertices are realized as floating point numbers.
    • "rhombic_triacontahedron" : Rhombic triacontahedron. Dual polytope to the icosidodecahedron, made of 30 rhombic facets.
    • "triakis_icosahedron" : Triakis icosahedron. Dual polytope to the icosidodecahedron, made of 30 rhombic facets.
    • "deltoidal_icositetrahedron" : Deltoidal Icositetrahedron. Dual polytope to the rhombicubaoctahedron, made of 24 kite facets.
    • "pentakis_dodecahedron" : Pentakis dodecahedron. Dual polytope to the truncated icosahedron, made of 60 isosceles triangular facets.
    • "deltoidal_hexecontahedron" : Deltoidal hexecontahedron. Dual polytope to the rhombicosidodecahedron, made of 60 kite facets.
    • "disdyakis_triacontahedron" : Disdyakis triacontahedron. Dual polytope to the truncated icosidodecahedron, made of 120 scalene triangular facets.

Examples

julia> T = catalan_solid("triakis_tetrahedron");

julia> count(F -> n_vertices(F) == 3, faces(T, 2))
12

julia> n_facets(T)
12
source
regular_24_cellFunction
regular_24_cell()

Construct the regular 24-cell, one out of three exceptional regular 4-polytopes.

source
regular_120_cellFunction
regular_120_cell()

Construct the regular 120-cell, one out of three exceptional regular 4-polytopes.

source
regular_600_cellFunction
regular_600_cell()

Construct the regular 600-cell, one out of three exceptional regular 4-polytopes.

source

Other polytope constructions

SIM_body_polytopeFunction
SIM_body_polytope(alpha::AbstractVector)

Produce an $n$-dimensional SIM-body as generalized permutahedron in $(n+1)$-space. SIM-bodies are defined in [GK14], but the input needs to be descending instead of ascending, as used in [JKS22], i.e. alpha has parameters $(a_1,\dots,a_n)$ such that $a_1 \geq \dots \geq a_n \geq 0$.

Example

To produce a $2$-dimensional SIM-body, use for example the following code. Note that the polytope lives in $3$-space, so we project it down to $2$-space by eliminating the last coordinate.

julia> s = SIM_body_polytope([3,1])
Polyhedron in ambient dimension 3

julia> p = convex_hull(map(x->x[1:dim(s)],vertices(s)))
Polyhedron in ambient dimension 2

julia> vertices(p) 
5-element SubObjectIterator{PointVector{QQFieldElem}}:
 [0, 0]
 [3, 0]
 [3, 1]
 [0, 3]
 [1, 3]
source
associahedronFunction
associahedron(d::Int)

Produce a $d$-dimensional associahedron (or Stasheff polytope). We use the facet description given in section 9.2. of [Zie95].

Note that in polymake, this function has an optional Boolean parameter group, to also construct the symmetry group of the polytope. For details, see [CSZ15].

Example

Produce the $2$-dimensional associahedron is a polygon in $\mathbb{R}⁴$ having $5$ vertices and $5$ facets.

julia> A =  associahedron(2)
Polyhedron in ambient dimension 4

julia> vertices(A)
5-element SubObjectIterator{PointVector{QQFieldElem}}:
 [9, 4, 1, 10]
 [10, 1, 4, 9]
 [1, 10, 1, 9]
 [1, 4, 9, 6]
 [4, 1, 10, 6]

julia> facets(A)
5-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^4 described by:
-x_1 <= -1
-2*x_1 - 2*x_2 <= -10
-x_2 <= -1
-2*x_2 - 2*x_3 <= -10
-x_3 <= -1
source
billera_lee_polytopeFunction
billera_lee_polytope(h::AbstractVector)

Construct a simplicial polytope whose h-vector is $h$. The corresponding g-vector must be an M-sequence. The ambient dimension equals the length of $h$, and the polytope lives in codimension one.

Examples

julia> BL = billera_lee_polytope([1,3,3,1])
Polyhedron in ambient dimension 4

julia> f_vector(BL)
3-element Vector{ZZRingElem}:
 6
 12
 8
source
binary_markov_graph_polytopeFunction
binary_markov_graph_polytope(observation::AbstractVector)

Defines a very simple graph for a polytope propagation related to a Hidden Markov Model. The length of observation is the number of possible oberservations. Its elements are of types Bool or Int. The propagated polytope is always a polygon. For a detailed description see [Jos05].

Examples

julia> P = binary_markov_graph_polytope([1,1,1,1])
Polyhedron in ambient dimension 2

julia> vertices(P)
4-element SubObjectIterator{PointVector{QQFieldElem}}:
 [3, 0]
 [1, 1]
 [0, 2]
 [0, 7]
source
birkhoff_polytopeFunction
birkhoff_polytope(n::Integer, even::Bool = false)

Construct the Birkhoff polytope of dimension $n^2$.

This is the polytope of $n \times n$ stochastic matrices (encoded as row vectors of length $n^2$), i.e., the matrices with non-negative real entries whose row and column entries sum up to one. Its vertices are the permutation matrices.

Use even = true to get the vertices only for the even permutation matrices.

Examples

julia> b = birkhoff_polytope(3)
Polytope in ambient dimension 9

julia> vertices(b)
6-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 0, 0, 0, 1, 0, 0, 0, 1]
 [0, 1, 0, 1, 0, 0, 0, 0, 1]
 [0, 0, 1, 1, 0, 0, 0, 1, 0]
 [1, 0, 0, 0, 0, 1, 0, 1, 0]
 [0, 1, 0, 0, 0, 1, 1, 0, 0]
 [0, 0, 1, 0, 1, 0, 1, 0, 0]
source
cyclic_caratheodory_polytopeFunction
cyclic_caratheodory_polytope(d::Int, n::Int)

Produce a $d$-dimensional cyclic polytope with $n$ points. Clearly $n\geq d$ is required. It is a prototypical example of a neighborly polytope whose combinatorics completely known due to Gale's evenness criterion. The coordinates are chosen on the trigonometric moment curve.

Example

julia> C= cyclic_caratheodory_polytope(4,5)
Polytope in ambient dimension 4

julia> vertices(C)
5-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 0, 1, 0]
 [347922205179541//1125899906842624, 8566355544790271//9007199254740992, -7286977268806823//9007199254740992, 5294298886396511//9007199254740992]
 [-7286977268806823//9007199254740992, 5294298886396511//9007199254740992, 1391688820718163//4503599627370496, -33462326346837//35184372088832]
 [-7286977268806825//9007199254740992, -5294298886396509//9007199254740992, 5566755282872661//18014398509481984, 8566355544790271//9007199254740992]
 [1391688820718163//4503599627370496, -33462326346837//35184372088832, -3643488634403413//4503599627370496, -5294298886396507//9007199254740992]
source
cyclic_polytopeFunction
cyclic_polytope(d::Int, n::Int)

Construct the cyclic polytope that is the convex hull of $n$ points on the moment curve in dimension $d$.

Examples

julia> cp = cyclic_polytope(3, 20)
Polytope in ambient dimension 3

julia> n_vertices(cp)
20
source
del_pezzo_polytopeFunction
del_pezzo_polytope(d::Int)

Produce the d-dimensional del Pezzo polytope, which is the convex hull of the cross polytope together with the all-ones and minus all-ones vector.

Examples

julia> DP = del_pezzo_polytope(4)
Polytope in ambient dimension 4

julia> f_vector(DP)
4-element Vector{ZZRingElem}:
 10
 40
 60
 30
source
dwarfed_cubeFunction
dwarfed_cube(d::Int)

Produce the $d$-dimensional dwarfed cube as defined in [ABS97].

Example

The $3$-dimensional dwarfed cube is illustrated in [Jos03].

julia> c = dwarfed_cube(3)
Polytope in ambient dimension 3

julia> vertices(c)
10-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1//2, 0, 1]
 [1//2, 1, 0]
 [1, 0, 1//2]
 [1, 1//2, 0]
 [1, 0, 0]
 [0, 0, 0]
 [0, 1, 0]
 [0, 1//2, 1]
 [0, 0, 1]
 [0, 1, 1//2]
source
dwarfed_product_polygonsFunction
dwarfed_product_polygons(d::Int, s::Int)

Produce a $d$-dimensional dwarfed product of polygons of size $s$ as defined in [ABS97]. It must be $d\geq4$ and even as well as $s\geq 3$.

Example

julia> p = dwarfed_product_polygons(4,3)
Polytope in ambient dimension 4

julia> vertices(p)
11-element SubObjectIterator{PointVector{QQFieldElem}}:
 [5, 3, 0, 0]
 [5, 0, 0, 0]
 [2, 0, 3, 9]
 [0, 0, 5, 3]
 [0, 0, 3, 9]
 [2, 6, 3, 9]
 [0, 0, 5, 0]
 [0, 0, 0, 0]
 [3, 9, 2, 6]
 [3, 9, 2, 0]
 [3, 9, 0, 0]
source
explicit_zonotopeFunction
explicit_zonotope(zones::Matrix; rows_are_points::Bool=true)

Produce the points of a zonotope as the iterated Minkowski sum of all intervals $[-x,x]$, where $x$ ranges over the rows of the input matrix zones. If rows_are_points is true (default), the rows of the input matrix represent affine points, otherwise they represent linear vectors.

Examples

julia> Z = explicit_zonotope([1 1; 1 -1], rows_are_points=false)
Polyhedron in ambient dimension 2

julia> vertices(Z)
4-element SubObjectIterator{PointVector{QQFieldElem}}:
 [2, 0]
 [0, -2]
 [0, 2]
 [-2, 0]
source
fano_simplexFunction
fano_simplex(d::Int)

Construct a lattice simplex such that the origin is the unique interior lattice point. The normal toric variety associated with its face fan is smooth.

Keywords

  • d::Int: the dimension.

Examples

julia> S = fano_simplex(3)
Polytope in ambient dimension 3

julia> X = normal_toric_variety(face_fan(S))
Normal toric variety

julia> is_smooth(X)
true
source
fractional_cut_polytopeFunction
fractional_cut_polytope(G::Graph{Undirected})

Construct the fractional cut polytope of the graph $G$.

Examples

julia> G = complete_graph(4);

julia> fractional_cut_polytope(G)
Polytope in ambient dimension 6
source
fractional_knapsack_polytopeFunction
fractional_knapsack_polytope(b::AbstractVector{<:Base.Number})

Produce a knapsack polytope defined by one linear inequality (and non-negativity constraints).

Example

julia> f = fractional_knapsack_polytope([10,-2,-3,-5])
Polytope in ambient dimension 3

julia> print_constraints(f)
2*x_1 + 3*x_2 + 5*x_3 <= 10
-x_1 <= 0
-x_2 <= 0
-x_3 <= 0
source
fractional_matching_polytopeFunction
fractional_matching_polytope(G::Graph{Undirected})

Construct the fractional matching polytope of the graph $G$.

Examples

julia> G = complete_graph(4);

julia> fractional_matching_polytope(G)
Polytope in ambient dimension 6
source
gelfand_tsetlin_polytopeFunction
gelfand_tsetlin_polytope(lambda::AbstractVector)

Construct the Gelfand-Tsetlin polytope indexed by a weakly decreasing vector lambda.

Examples

julia> P = gelfand_tsetlin_polytope([5,3,2])
Polyhedron in ambient dimension 6

julia> is_fulldimensional(P)
false

julia> p = project_full(P)
Polyhedron in ambient dimension 3

julia> is_fulldimensional(p)
true

julia> volume(p)
3
source
gelfand_tsetlin_polytope(lambda::AbstractVector, sigma::PermGroupElem)

Construct the generalized Gelfand-Tsetlin polytope indexed by a weakly decreasing vector lambda and a permutation sigma.

julia> P = gelfand_tsetlin_polytope([5,3,2], @perm (1,3,2))
Polyhedron in ambient dimension 6
source
goldfarb_cubeFunction
goldfarb_cube(d::Int, e::Number, g::Number)

Produce a d-dimensional Goldfarb cube. The first parameter of deformation e must be $<\frac{1}{2}$, the second parameter of deformation d must be $\geq \frac{\texttt{e}}{4}$.

The Goldfarb cube is a combinatorial cube and yields a bad example for the Simplex Algorithm using the Shadow Vertex Pivoting Strategy. Here we use the description as a deformed product due to [AZ99]. For $g=0$ we obtain a Klee-Minty cube, in particular for $e=g=0$ we obtain the standard cube.

Example

The following produces a $3$-dimensional Klee-Minty cube for $e=\frac{1}{3}$.

julia> c = goldfarb_cube(3,1//3,0)
Polytope in ambient dimension 3

julia> vertices(c)
8-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 1//3, 8//9]
 [1, 2//3, 7//9]
 [1, 2//3, 2//9]
 [1, 1//3, 1//9]
 [0, 0, 0]
 [0, 1, 1//3]
 [0, 1, 2//3]
 [0, 0, 1]
source
goldfarb_sit_cubeFunction
goldfarb_sit_cube(d::Int, eps::Number, delta::Number)

Produces a d-dimensional variation of the Klee-Minty cube, which is scaled in direction $x_{d-i}$ by eps*delta^i. The first parameter of deformation eps must be $<\frac{1}{2}$, the second parameter of deformation delta must be $\geq \frac{1}{2}$. This cube is a combinatorial cube and yields a bad example for the Simplex Algorithm using the Steepest Edge Pivoting Strategy. Here we use a scaled description of the construction of Goldfarb and Sit, see [GS79].

Examples

julia> c = goldfarb_sit_cube(3,1//3,1//2)
Polytope in ambient dimension 3

julia> vertices(c) 
8-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1//36, 1//18, 8//9]
 [1//36, 1//9, 7//9]
 [1//36, 1//9, 2//9]
 [1//36, 1//18, 1//9]
 [0, 0, 0]
 [0, 1//6, 1//3]
 [0, 1//6, 2//3]
 [0, 0, 1]
source
hypersimplexFunction
hypersimplex(k::Int, d::Int; no_vertices::Bool=false, no_facets::Bool=false, no_vif::Bool=false)

Produce the hypersimplex $\Delta(k,d)$, that is the the convex hull of all $0/1$-vector in $\mathbb{R}^d$ with exactly $k$ ones. Note that the output is never full-dimensional.

Optional Arguments

  • no_vertices::Bool: If set equal to true, vertices of the underlying polymake object are not computed.
  • no_facets::Bool: If set equal to true, facets of the underlying polymake object are not computed.
  • no_vif::Bool: If set equal to true, vertices in facets of the underlying polymake object are not computed.

Example

julia> H = hypersimplex(3,4)
Polytope in ambient dimension 4

julia> G = hypersimplex(3,4,no_facets=true)
Polytope in ambient dimension 4

julia> facets(G)
4-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^4 described by:
x_4 <= 1
x_3 <= 1
-x_1 - x_3 - x_4 <= -2
x_1 <= 1

julia> facets(H)
4-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^4 described by:
x_4 <= 1
x_3 <= 1
-x_1 - x_3 - x_4 <= -2
x_1 <= 1
source
hypertruncated_cubeFunction
hypertruncated_cube(d::Int, k::Number, lambda::Number)

Produce a $d$-dimensional hypertruncated cube with symmetric linear objective function $(1,1,…,1)$.

Arguments

  • k: cutoff parameter
  • lambda: scaling of extra vertex

Example

julia> H = hypertruncated_cube(3,2,3)
Polytope in ambient dimension 3

julia> print_constraints(H)
-x_1 <= 0
-x_2 <= 0
-x_3 <= 0
x_1 <= 1
x_2 <= 1
x_3 <= 1
5*x_1 - 2*x_2 - 2*x_3 <= 3
-2*x_1 + 5*x_2 - 2*x_3 <= 3
-2*x_1 - 2*x_2 + 5*x_3 <= 3
source
k_cyclic_polytopeFunction
k_cyclic_polytope(n::Int, s::Vector)

Produce a (rounded) $2*k$-dimensional $k$-cyclic polytope with n points, where $k$ is the length of the input vector s. Special cases are the bicyclic ($k=2$) and tricyclic ($k=3$) polytopes. Only possible in even dimensions.

The parameters $\texttt{s}_i$ can be integers, floating-points or rational numbers. The $i$-th vertex then is: $(\cos(\texttt{s}_1 * 2\pi i/\texttt{n}), \sin(\texttt{s}_1 * 2\pi i/\texttt{n}), ... , \cos(\texttt{s}_k * 2\pi i/\texttt{n}), \sin(\texttt{s}_k * 2\pi i/\texttt{n}))$.

Warning: Some of the $k-$cyclic polytopes are not simplicial. Since the components are rounded, this function might output a polytope which is not a $k-$cyclic polytope! More information see [Sch95].

Example

To produce a (not exactly) regular pentagon, type this:

julia> p = k_cyclic_polytope(5,[1])
Polytope in ambient dimension 2

julia> dim(p) 
2

julia> n_vertices(p)
5
source
klee_minty_cubeFunction
klee_minty_cube(d::Int, e::Number)

Produces a $d-$dimensional Klee-Minty-cube if $\texttt{e} < 1/2$. Uses the goldfarb_cube method with the argument $\texttt{g} = 0$.

#Example

julia> k = klee_minty_cube(3,1//8)
Polytope in ambient dimension 3

julia> print_constraints(k)
-x_1 <= 0
x_1 <= 1
1//8*x_1 - x_2 <= 0
1//8*x_1 + x_2 <= 1
1//8*x_2 - x_3 <= 0
1//8*x_2 + x_3 <= 1
source
lecture_hall_simplexFunction
lecture_hall_simplex(d::Int)

Produce the $d$-dimensional lecture hall simplex for the sequence $(s_i)=i$ for $1\geq i \geq d$ as defined in [SS12].

Note that in polymake, this function has an optional Boolean parameter group, to also construct the symmetry group of the simplex.

Example

The $3$-dimensional lecture hall simplex:

julia> S = lecture_hall_simplex(3) 
Polytope in ambient dimension 3

julia> vertices(S)
4-element SubObjectIterator{PointVector{QQFieldElem}}:
 [0, 0, 0]
 [0, 0, 3]
 [0, 2, 3]
 [1, 2, 3]
source
max_GC_rank_polytopeFunction
max_GC_rank_polytope(d::Int)

Produce a d-dimensional polytope of maximal Gomory-Chvatal rank $\Omega(d/\log(d))$, integrally infeasible. With symmetric linear objective function $(1,1..,1)$. Construction due to Pokutta and Schulz, see [PS11].

Example

julia> c = max_GC_rank_polytope(3)
Polytope in ambient dimension 3

julia> vertices(c)
6-element SubObjectIterator{PointVector{QQFieldElem}}:
 [0, 1//2, 1//2]
 [1//2, 0, 1//2]
 [1//2, 1//2, 0]
 [1//2, 1, 1//2]
 [1//2, 1//2, 1]
 [1, 1//2, 1//2]
source
n_gonFunction
n_gon(n::Int; r::RationalUnion=1, alpha_0::RationalUnion=0)

Produce a regular n-gon. All vertices lie on a circle of radius r (defaults to $1$) and initial angle divided by pi alpha_0 (defaults to $0$).

Examples

To store the regular pentagon in the variable p, do this:

julia> p = n_gon(3)
Polytope in ambient dimension 2

julia> n_gon(3,r=3)
Polytope in ambient dimension 2
source
newton_polytopeFunction
newton_polytope(poly::Polynomial)

Compute the Newton polytope of the multivariate polynomial poly.

Examples

julia> S, (x, y) = polynomial_ring(ZZ, ["x", "y"])
(Multivariate polynomial ring in 2 variables over ZZ, ZZMPolyRingElem[x, y])

julia> f = x^3*y + 3x*y^2 + 1
x^3*y + 3*x*y^2 + 1

julia> NP = newton_polytope(f)
Polyhedron in ambient dimension 2

julia> vertices(NP)
3-element SubObjectIterator{PointVector{QQFieldElem}}:
 [3, 1]
 [1, 2]
 [0, 0]
source
orbit_polytopeFunction
orbit_polytope(V::AbstractCollection[PointVector], G::PermGroup)

Construct the convex hull of the orbit of one or several points (given row-wise in V) under the action of G.

Examples

This will construct the $3$-dimensional permutahedron:

julia> V = [1 2 3];

julia> G = symmetric_group(3);

julia> P = orbit_polytope(V, G)
Polyhedron in ambient dimension 3

julia> vertices(P)
6-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 2, 3]
 [1, 3, 2]
 [2, 1, 3]
 [2, 3, 1]
 [3, 1, 2]
 [3, 2, 1]
source
permutahedronFunction
permutahedron(d::Int)

Produce a d-dimensional permutahedron. The vertices correspond to the elements of the symmetric group of degree d$+1$.

#Example

julia> p = permutahedron(2)
Polytope in ambient dimension 3

julia> vertices(p)
6-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 2, 3]
 [1, 3, 2]
 [2, 1, 3]
 [2, 3, 1]
 [3, 1, 2]
 [3, 2, 1]
source
pile_polytopeFunction
pile_polytope(sizes::Vector{Int})

Produce a $($d$+1)$-dimensional polytope from a pile of cubes. Start with a d-dimensional pile of cubes. Take a generic convex function to lift this polytopal complex to the boundary of a $($d$+1)$–polytope. The argument sizes is a vector $(s_1,…,s_d)$ where $s_i$ specifies the number of boxes in the $i$-th dimension.

source
pitman_stanley_polytopeFunction
pitman_stanley_polytope(y::AbstractVector)

Produce a Pitman-Stanley polytope of dimension $n-1$, where y is a Vector of $n$ positive parameters. Does not check if the parameters are actually positive; negative values are legal but that do not yield a Pitman-Stanley polytope. Zeros just reduce the dimension; negative numbers may produce unbounded polyhedra.

Example:

Pitman-Stanley polytopes are combinatorial cubes:

julia> p = pitman_stanley_polytope([1,2,3])
Polyhedron in ambient dimension 3

julia> f_vector(p) 
2-element Vector{ZZRingElem}:
 4
 4
source
perles_nonrational_8_polytopeFunction
perles_nonrational_8_polytope()

Create an $8$-dimensional polytope without rational realizations due to Perles. See [Gru03].

Example

julia> perles_nonrational_8_polytope()
Polytope in ambient dimension 8 with EmbeddedAbsSimpleNumFieldElem type coefficients
source
pseudo_del_pezzo_polytopeFunction
pseudo_del_pezzo_polytope(d::Int)

Produce a d-dimensional del-Pezzo polytope, which is the convex hull of the cross polytope together with the all-ones vector. All coordinates are plus or minus one.

Example

julia> DP = pseudo_del_pezzo_polytope(4)
Polytope in ambient dimension 4

julia> f_vector(DP)
4-element Vector{ZZRingElem}:
 9
 32
 46
 23
source
rand01_polytopeFunction
rand01_polytope(d::Int, n::Int; seed=nothing)

Produce a d-dimensional $0/1$-polytope with n random vertices. Uniform distribution.

Optional Argument

-seed::Int: Seed for random number generation

Example

julia> s = rand01_polytope(2, 4; seed=3)
Polytope in ambient dimension 2

julia> vertices(s)
4-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 1]
 [1, 0]
 [0, 0]
 [0, 1]
source
rand_box_polytopeFunction
rand_box_polytope(d::Int, n::Int, b::Int; seed::Int=nothing)

Computes the convex hull of n points sampled uniformly at random from the integer points in the cube $\[0,\texttt{b}\]^{\textttt{d}}$.

Optional Argument

-seed: Seed for random number generation.

Example

julia> r = rand_box_polytope(3, 10, 3, seed=1)
Polyhedron in ambient dimension 3

julia> vertices(r) 
8-element SubObjectIterator{PointVector{QQFieldElem}}:
 [3, 2, 3]
 [0, 3, 0]
 [3, 3, 0]
 [3, 0, 1]
 [1, 1, 0]
 [2, 0, 3]
 [0, 3, 3]
 [0, 1, 2]
source
rand_cyclic_polytopeFunction
rand_cyclic_polytope(d::Int, n::Int; seed::Int=nothing)

Computes a random instance of a cyclic polytope of dimension d on n vertices by randomly generating a Gale diagram whose cocircuits have alternating signs.

Optional Argument

-seed: Seed for random number generation

Examples

julia> r = rand_cyclic_polytope(3, 5)
Polytope in ambient dimension 3

julia> f_vector(r)
3-element Vector{ZZRingElem}:
 5
 9
 6
source
rand_metricFunction
rand_metric(n::Int; seed=nothing)

Produce a rational n-point metric with random distances. The values are uniformily distributed in $\[1,2\]$.

Examples

julia> rand_metric(3, seed=132)
[                               0   260222460282405//140737488355328   371474612593257//281474976710656]
[260222460282405//140737488355328                                  0   388326899436839//281474976710656]
[371474612593257//281474976710656   388326899436839//281474976710656                                  0]
source
rand_metric_intFunction
rand_metric_int(n::Int, digits::Int; seed=nothing)

Produce a n-point metric with random integral distances. The values are uniformily distributed in $\[1,2\]$. The distances are integers and lie in $[10^digits, 10^(digits+1)[$.

source
rand_normal_polytopeFunction
rand_normal_polytope(d::Int, n::Int; seed=nothing, precision=nothing)

Produce a rational d-dimensional polytope from n random points approximately normally distributed in the unit ball.

Optional Arguments

-seed: controls the outcome of the random number generator; fixing a seed number guarantees the same outcome -precision: number of bits for MPFR sphere approximation

Example

julia> rnp = rand_normal_polytope(2,4; seed=42, precision=4)
Polytope in ambient dimension 2

julia> is_simplicial(rnp)
true

julia> sort(map(x->dot(x,x), vertices(rnp)))
4-element Vector{QQFieldElem}:
 1417//4096
 481//1024
 225//256
 101//32
source
rand_spherical_polytopeFunction
rand_spherical_polytope([rng::AbstractRNG,] d::Int, n::Int;
distribution=:uniform, precision=nothing, seed=nothing)

Construct the convex hull of $n$ points on the unit sphere in $\mathbb{R}^d$. Almost surely this is a simplicial polytope.

Keywords

  • distribution::Symbol: One of the following two options:
    • :uniform (default): Use intermediate floating point numbers for an almost uniform distribution on the sphere. The points will not be exactly on the sphere.
    • :exact: Create exact rational points on the unit sphere, this works at the expense of both uniformity and log-height of the points.
  • precision::Int64: Precision in bits during floating point approximation for uniform distribution.
  • seed::Int64: Seed for random number generation. Cannot be used together with the AbstractRNG argument.

Examples

julia> rsph = rand_spherical_polytope(3, 20)
Polytope in ambient dimension 3

julia> is_simplicial(rsph)
true

julia> rsph = rand_spherical_polytope(3, 4; precision=5, seed=132)
Polytope in ambient dimension 3

julia> map(x->dot(x,x), vertices(rsph))
4-element Vector{QQFieldElem}:
 4306545//4194304
 15849//16384
 4165//4096
 8281//8192

julia> rsph = rand_spherical_polytope(3, 4; distribution=:exact)
Polytope in ambient dimension 3

julia> map(x->dot(x,x), vertices(rsph))
4-element Vector{QQFieldElem}:
 1
 1
 1
 1
source
rand_subpolytopeFunction
rand_subpolytope(P::Polyhedron, n::Int; seed=nothing)

Construct a subpolytope of $P$ as the convex hull of $n$ vertices, chosen uniformly at random. The polyhedron $P$ must be bounded, and the number $n$ must not exceed the number of vertices.

Keywords

  • seed::Int64: Seed for random number generation.

Examples

julia> n_vertices(rand_subpolytope(cube(3), 5))
5
source
rss_associahedronFunction

rss_associahedron(n::Int)

Produce a polytope of constrained expansions in ambient dimension n according to [RSS03].

Examples:

To produce a $3$-dimensional associahedron in $5$-space, do:

julia> a= rss_associahedron(5)
Polyhedron in ambient dimension 5

julia> vertices(a)
14-element SubObjectIterator{PointVector{QQFieldElem}}:
 [0, 1, 12, 13, 16]
 [0, 7, 8, 11, 16]
 [0, 1, 4, 9, 16]
 [0, 3, 4, 9, 16]
 [0, 5, 6, 9, 16]
 [0, 5, 8, 9, 16]
 [0, 1, 8, 9, 16]
 [0, 7, 10, 11, 16]
 [0, 7, 12, 13, 16]
 [0, 7, 12, 15, 16]
 [0, 1, 12, 15, 16]
 [0, 7, 8, 15, 16]
 [0, 1, 4, 15, 16]
 [0, 3, 4, 15, 16]

julia> facets(a) 
9-element SubObjectIterator{AffineHalfspace{QQFieldElem}} over the Halfspaces of R^5 described by:
x_1 - x_2 <= -1
x_1 - x_3 <= -4
x_1 - x_4 <= -9
x_2 - x_3 <= -1
x_2 - x_4 <= -4
x_2 - x_5 <= -9
x_3 - x_4 <= -1
x_3 - x_5 <= -4
x_4 - x_5 <= -1
source
signed_permutahedronFunction
signed_permutahedron(d::Int)

Produce the d-dimensional signed permutahedron. I.e. for all possible permutations of the vector $(1,\dots,d)$, all possible sign patterns define vertices of this polytope. Contrary to the classical permutahedron, the signed permutahedron is full-dimensional.

Examples:

To produce the $2$-dimensional signed permutahedron, do:

julia> P = signed_permutahedron(2)
Polytope in ambient dimension 2

julia> vertices(P)
8-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 2]
 [-1, 2]
 [1, -2]
 [-1, -2]
 [2, 1]
 [-2, 1]
 [2, -1]
 [-2, -1]
source
stable_set_polytopeFunction
stable_set_polytope(G::Graph{Undirected})

Produces the stable set polytope from an undirected graph G$=(V,E)$. The stable set Polytope has the following inequalities: $x_i + x_j \leq 1 \forall \{i,j\} \in E$, $x_i \geq 0 \forall i \in V$ and $x_i \leq 1 \forall i \in V \text{ with } \mathrm{deg}(i)=0$

Example:

The following produces first the standard cube in $3$ dimensions, and then a bipyramid over the convex hull of the unit vectors.

julia> G = Graph{Undirected}(3)
Undirected graph with 3 nodes and no edges

julia> S = stable_set_polytope(G)
Polytope in ambient dimension 3

julia> vertices(S)
8-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 0, 0]
 [1, 1, 0]
 [1, 1, 1]
 [1, 0, 1]
 [0, 0, 1]
 [0, 0, 0]
 [0, 1, 0]
 [0, 1, 1]

julia> add_edge!(G, 1, 2);

julia> add_edge!(G, 1, 3);

julia> add_edge!(G, 2, 3);

julia> S = stable_set_polytope(G)
Polytope in ambient dimension 3

julia> vertices(S)
5-element SubObjectIterator{PointVector{QQFieldElem}}:
 [0, 0, 1]
 [0, 1, 0]
 [1, 0, 0]
 [1//2, 1//2, 1//2]
 [0, 0, 0]
source
transportation_polytopeFunction
transportation_polytope(r::AbstractVector, c::AbstractVector)

Produce the transportation polytope from two vectors r of length $m$ and c of length $n$, i.e. all positive $m\times n$ Matrizes with row sums equal to r and column sums equal to c.

Example:

We can see that the set of $3\times 3$ magic squares with magic constant $15$ is a $4$-dimensional polytope.

julia> r = c = [15,15,15]
3-element Vector{Int64}:
 15
 15
 15

julia> t = transportation_polytope(r,c) 
Polytope in ambient dimension 9

julia> dim(t) 
4

julia> is_bounded(t) 
true
source
zonotopeFunction
zonotope(M::Matrix{<:Number}; centered::Bool=true)

Create a zonotope from a matrix whose rows are input points.

Optional Arguments

  • centered::Bool: This is true if the output should be centered; the default is true.

Examples

The following produces a parallelogram with the origin as its vertex barycenter:

julia> Z = zonotope([1 0; 1 1])
Polyhedron in ambient dimension 2

julia> vertices(Z)
4-element SubObjectIterator{PointVector{QQFieldElem}}:
 [-1, -1//2]
 [0, -1//2]
 [0, 1//2]
 [1, 1//2]

The following produces a parallelogram with the origin being a vertex (not centered case):

julia> Z = zonotope([1 0; 1 1], centered = false)
Polyhedron in ambient dimension 2

julia> vertices(Z)
4-element SubObjectIterator{PointVector{QQFieldElem}}:
 [0, 0]
 [1, 0]
 [1, 1]
 [2, 1]
source
zonotope_vertices_fukuda_matrixFunction
zonotope_vertices_fukuda(M::Matrix)

Create the vertices of a zonotope from a matrix whose rows are input points or vectors.

Examples

The following creates the vertices of a parallelogram with the origin as its vertex barycenter.

julia> zonotope_vertices_fukuda_matrix([1 1 0; 1 1 1])
pm::Matrix<pm::Rational>
-1 -1 -1/2
0 0 -1/2
0 0 1/2
1 1 1/2
source

Operations on polyhedra

Polyhedra can be produced through operations on other polyhedra. For example, they can be added using Minkowski addition or scaled; each of which results in a new polyhedron.

+Method
+(P::Polyhedron, Q::Polyhedron)

Return the Minkowski sum $P + Q = \{ x+y\ |\ x∈P, y∈Q\}$ of P and Q (see also minkowski_sum).

Examples

The Minkowski sum of a square and the 2-dimensional cross-polytope is an octagon:

julia> P = cube(2);

julia> Q = cross_polytope(2);

julia> M = minkowski_sum(P, Q)
Polyhedron in ambient dimension 2

julia> n_vertices(M)
8
source
*Method
*(k::Union{Number, FieldElem}, Q::Polyhedron)

Return the scaled polyhedron $kQ = \{ kx\ |\ x∈Q\}$.

Note that k*Q = Q*k.

Examples

Scaling an $n$-dimensional bounded polyhedron by the factor $k$ results in the volume being scaled by $k^n$. This example confirms the statement for the 6-dimensional cube and $k = 2$.

julia> C = cube(6);

julia> SC = 2*C
Polyhedron in ambient dimension 6

julia> volume(SC)//volume(C)
64
source
*Method
*(P::Polyhedron, Q::Polyhedron)

Return the Cartesian product of P and Q (see also product).

Examples

The Cartesian product of a triangle and a line segment is a triangular prism.

julia> T=simplex(2)
Polytope in ambient dimension 2

julia> S=cube(1)
Polytope in ambient dimension 1

julia> length(vertices(T*S))
6
source
bipyramidFunction
bipyramid(P::Polyhedron, z::Union{Number, FieldElem} = 1, z_prime::Union{Number, FieldElem} = -z)

Make a bipyramid over a pointed polyhedron P.

The bipyramid is the convex hull of the input polyhedron P and two apexes (v, z), (v, z_prime) on both sides of the affine span of P. For bounded polyhedra, the projections of the apexes v to the affine span of P is the vertex barycenter of P.

Examples

julia> c = cube(2)
Polytope in ambient dimension 2

julia> vertices(bipyramid(c,2))
6-element SubObjectIterator{PointVector{QQFieldElem}}:
 [-1, -1, 0]
 [1, -1, 0]
 [-1, 1, 0]
 [1, 1, 0]
 [0, 0, 2]
 [0, 0, -2]
source
intersectMethod
intersect(P::Polyhedron...)

Return the intersection $\bigcap\limits_{p \in P} p$.

Examples

The positive orthant of the plane is the intersection of the two halfspaces with $x≥0$ and $y≥0$ respectively.

julia> UH1 = convex_hull([0 0],[1 0],[0 1]);

julia> UH2 = convex_hull([0 0],[0 1],[1 0]);

julia> PO = intersect(UH1, UH2)
Polyhedron in ambient dimension 2

julia> rays(PO)
2-element SubObjectIterator{RayVector{QQFieldElem}}:
 [1, 0]
 [0, 1]
source
pyramidFunction
pyramid(P::Polyhedron, z::Union{Number, FieldElem} = 1)

Make a pyramid over the given polyhedron P.

The pyramid is the convex hull of the input polyhedron P and a point v outside the affine span of P. For bounded polyhedra, the projection of v to the affine span of P coincides with the vertex barycenter of P. The scalar z is the distance between the vertex barycenter and v.

Examples

julia> c = cube(2)
Polytope in ambient dimension 2

julia> vertices(pyramid(c,5))
5-element SubObjectIterator{PointVector{QQFieldElem}}:
 [-1, -1, 0]
 [1, -1, 0]
 [-1, 1, 0]
 [1, 1, 0]
 [0, 0, 5]
source
vertex_figureFunction
vertex_figure(P::Polyhedron, n::Int; cutoff=1//2)

Construct the vertex figure of the vertex n of a bounded polytope. The vertex figure is dual to a facet of the dual polytope.

Optional Arguments

  • cutoff::Number: controls the exact location of the cutting hyperplane. It should lie in the open Interval $(0,1)$. Value $0$ would let the hyperplane go through the chosen vertex, thus degenerating the vertex figure to a single point. Value $1$ would let the hyperplane touch the nearest neighbor vertex of a polyhedron. Default value is $\frac{1}{2}$.

Example

To produce a triangular vertex figure of a $3$-dimensional cube in the positive orthant, do:

julia> T = vertex_figure(cube(3), 8) 
Polyhedron in ambient dimension 3

julia> vertices(T)
3-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 1, 0]
 [1, 0, 1]
 [0, 1, 1]

julia> T = vertex_figure(cube(3), 8, cutoff = 1/4)
Polyhedron in ambient dimension 3

julia> vertices(T)
3-element SubObjectIterator{PointVector{QQFieldElem}}:
 [1, 1, 1//2]
 [1, 1//2, 1]
 [1//2, 1, 1]
source

The convex hull of two polytopes can be computed via convex_hull.

convex_hullMethod
convex_hull(P::Polyhedron, Q::Polyhedron)

Return the convex_hull of P and Q.

Examples

The convex hull of the following two line segments in $R^3$ is a tetrahedron.

julia> L₁ = convex_hull([-1 0 0; 1 0 0])
Polyhedron in ambient dimension 3

julia> L₂ = convex_hull([0 -1 0; 0 1 0])
Polyhedron in ambient dimension 3

julia> T=convex_hull(L₁,L₂);

julia> f_vector(T)
2-element Vector{ZZRingElem}:
 4
 4
source