Polygons

In RationalPolygons.jl, we represent a polygon $P \subseteq \mathbb{R}^2$ by two pieces of data: An integral matrix $V \in \mathbb{Z}^{2\times N}$, called the vertex matrix and an integer $k \in \mathbb{Z}$, called the rationality. The associated polygon has as vertices the columns of $V$ divided by $k$. To represent $V$, we use a static matrix, which are faster than Julia's internal matrices for many common operations. However, this implies that the type RationalPolygon will depend on the number of vertices N as a type parameter, which affects Julia's dispatch mechanism at runtime. As long as the number of distinct values of N occurring during a computation remains relatively small, this should not cause performance issues.

There are two ways in which our encoding of rational polygons is not unique: First, scaling $V$ and $k$ by the same factor does not change the polygon, e.g. $(V,k)$ describes the same polygon as $(2V,2k)$. Even though they are mathematically the same polygon, RationalPolygon.jl views them as different objects: once as a $k$-rational polygon and once as a $2k$-rational polygon. The second way in which this encoding is not unique is that we can change the order of the columns. While we require them to be sorted counterclockwise, we may use any vertex as the first column. This problem is addressed in the section on normal forms.

RationalPolygons.RationalPolygonType
RationalPolygon{T<:Integer,N,M}

The type of rational polygons in two-dimensional space. T is the type of integers to be used. N is the number of vertices of the polygon and M equals 2*N. It has the following fields:

  • rationality :: T: The rationality of the polygon, e.g. 1 for lattice polygons, 2 for half-integral polygons etc.
  • vertex_matrix :: SMatrix{2,N,T,M}: An integral 2xN matrix. The vertices of the polygon are understood to be the columns of this matrix divided by rationality.
  • number_of_vertices :: Int: The number of vertices of the polygon. This is redundant information, since the number of vertices is already available as the type parameter N. However, getting the number of vertices of a polygon through the type paremeter means lots of work for Julia's dispatch algorithm. Therefore, we found it to improve performance to put it as a variable into the struct as well.
  • is_unimodular_normal_form :: Bool: A flag to remember that the polygon is already in unimodular normal form. This is used internally to avoid redundant computations of the normal form.
  • is_affine_normal_form :: Bool: A flag variable to remember that the polygon is already in affine normal form. This is used internally to avoid redundant computations of the normal form.
source

Constructors

Besides the type constructor methods, we provide the functions convex_hull and intersect_halfplanes to create a polygon from an unstructured collection of points or affine halfplanes.

RationalPolygons.RationalPolygonMethod
RationalPolygon(vertex_matrix :: SMatrix{2,N,T,M}, rationality :: T) where {N, M, T <: Integer}
RationalPolygon(scaled_vertices :: Vector{LatticePoint{T}}, rationality :: T)
RationalPolygon(vertices :: Vector{RationalPoint{T}})
RationalPolygon(vertices :: Vector{RationalPoint{T}}, rationality :: T)

Construct a rational polygon from a given set of vertices and possibly prescribed rationality. If vertex_matrix or scaled_vertices is provided, these are understood to be the integral vertices of the scaled polygon rationality(P) * P.

When using this constructor, no consistency checks are done on the input. In particular, the user must be sure that the given points are truly vertices of the polygon and that they are ordered counterclockwise. If this is not known ahead of contraction, convex_hull should be used instead of this constructor.

All constructors accept the optional keyword arguments is_unimodular_normal_form and is_affine_normal_form, which are set to false by default. If they are set to true, this means the user is certain that the given polygon is already in the respective normal form. This information will be used to prevent additional computations of the normal form and thus speed up equivalence checking.

Example:

The standard lattice triangle.

julia> P = RationalPolygon(LatticePoint{Int}[(0,0),(1,0),(0,1)], 1)
Rational polygon of rationality 1 with 3 vertices.   

Example:

A half-integral polygon.

julia> P = RationalPolygon(RationalPoint{Int}[(1,0),(0,1//2),(-1,0),(0,-1//2)])
Rational polygon of rationality 2 with 4 vertices.
source
RationalPolygons.convex_hullFunction
convex_hull(points :: Vector{LatticePoint{T}}, k :: T = one(T)) where {T <: Integer}

Return the k-rational polygon given by the convex hull of p // k, where p ∈ points.

source
convex_hull(points :: Vector{RationalPoint{T}}) where {T <: Integer}

Return the convex hull of a given set of rational points. The rationality will be inferred from the input.

Example:

julia> convex_hull(RationalPoint{Int}[(1,0),(0,1//2),(0,-1//3)])
Rational polygon of rationality 6 with 3 vertices.
source
convex_hull(points :: Vector{RationalPoint{T}}, k :: T) where {T <: Integer}

Return the convex hull of a given set of rational points, viewed as a k-rational polygon.

Example:

julia> convex_hull(RationalPoint{Int}[(1,0),(0,1//2),(0,-1//3)], 12)
Rational polygon of rationality 12 with 3 vertices.
source
RationalPolygons.intersect_halfplanesFunction
intersect_halfplanes(halfplanes :: Vector{AffineHalfplane{T}}; rationality :: Union{Missing,T} = missing) where {T <: Integer}

Return the polygon described by the intersection of the given affine halfplanes. Throws an error if the intersection is unbounded. If rationality is not provided, it will be inferred from the input.

source
RationalPolygons.empty_polygonFunction
empty_polygon(rationality :: T) where {T <: Integer}

Return the empty polygon of given rationality.

source
empty_polygon(::Type{T}) where {T <: Integer}

Return the empty polygon of integer type T. The rationality is understood to be one, i.e. it is a lattice polygon.

source

Basic Properties

We provide basic properties and checks for rational polygons. Note that indices corresponding to vertices are always considered cyclic, i.e. the $N+1$-th vertex of a polygon with $N$ vertices cycles back to its first vertex.

RationalPolygons.rationalityMethod
rationality(P :: RationalPolygon)

Return the rationality of P. Note that this does not need to be the denominator of P, which is the smallest positive integer k such that k*P is a lattice polygon: The standard lattice triangle may also be viewed as a half-integral polygon, in which case the rationality would be two, but the denominator is one.

source
Base.denominatorMethod
Base.denominator(P :: RationalPolygon{T,N}) where {N,T <: Integer}

Return the smallest positive integer k such that k * P is a lattice polygon.

source
RationalPolygons.vertex_matrixFunction
vertex_matrix(P :: RationalPolygon)

The vertex matrix of P is the 2xN integral matrix containing the vertices of rationality(P) * P as its columns.

Example:

julia> P = convex_hull(RationalPoint{Int}[(1,0),(0,1//2),(0,-1//3)])
Rational polygon of rationality 6 with 3 vertices.

julia> vertex_matrix(P)
2×3 StaticArraysCore.SMatrix{2, 3, Int64, 6} with indices SOneTo(2)×SOneTo(3):
  0  6  0
 -2  0  3
source
RationalPolygons.scaled_vertexFunction
scaled_vertex(P :: RationalPolygon, i :: Int)

Return the i-th integral vertex of the polygon rationality(P) * P. The index i is regarded as a cyclic index, e.g. the N+1-th vertex is equal to the first vertex.

source
RationalPolygons.vertexFunction
vertex(P :: RationalPolygon, i :: Int)

Return the i-th vertex of the polygon P. The index i is regarded as a cyclic index, e.g. the N+1-th vertex is equal to the first vertex.

source
Base.inMethod
Base.in(x :: Point{T}, P :: RationalPolygon{T}) where {T <: Integer}

Check whether a point x is contained in the polygon P.

source
RationalPolygons.dimFunction
dim(P :: RationalPolygon)

Return the dimension of P. For empty polygons, this returns -1. Otherwise, it returns 0, 1 or 2.

source
RationalPolygons.normalized_areaFunction
normalized_area(P :: RationalPolygon)

Return the normalized area of a k-rational polygon, i.e. 2k^2 times its euclidean area. The result is always an integer, counting the number of standard k-rational triangles contained in P.

source
RationalPolygons.is_maximalFunction
is_maximal(P :: RationalPolygon)

Check whether a rational polygon is maximal among all polygons sharing the same rationality and number of interior lattice points.

source
RationalPolygons.dualFunction
dual(P :: RationalPolygon{T}) where {T <: Integer}

Return the dual of a polygon P. Throws an error if P does not contain the origin in its interior.

source

Ehrhart Theory

The main result of Ehrhart Theory is that the number of lattice points in integral multiples of a $k$-rational polygon $P$ is a quasipolynomial, called its Ehrhart quasipolynomial:

\[\mathrm{ehr_P}(t) = |tP \cap \mathbb{Z}^2| = At^2 + a(t)t+b(t), \qquad t \in \mathbb{Z}.\]

Here, $A$ is the euclidean area of $P$ and $a, b\colon \mathbb{Z} \to \mathbb{Q}$ are $k$-periodic functions. These can be computed by

\[a(t) = -(2t+k)\cdot A + \frac{\mathrm{ehr}_P(t+k)-\mathrm{ehr}_P(t)}{k},\]

\[b(t) = (t^2+tk)\cdot A + \frac{(t+k)\mathrm{ehr}_P(t)-t\mathrm{ehr}_P(t+k)}{k}.\]

Setting $\tilde{A} := 2k^2A,\ \tilde{a} := 2k^2 a$ and $\tilde{b} := 2k^2b$, we get integer valued functions $\tilde a$ and $\tilde b$ which we call the normalized Ehrhart coefficients. Hence we can encode the Ehrhart quasipolynomial by the $3\times k$-integral matrix of its normalized Ehrhart coefficients:

\[\begin{bmatrix} \tilde{A} & \tilde{a}(1) & \tilde{b}(1) \\ \tilde{A} & \tilde{a}(2) & \tilde{b}(2) \\ \vdots & \vdots & \vdots \\ \tilde{A} & \tilde{a}(k) & \tilde{b}(k) \end{bmatrix} \in \mathbb{Z}^{3\times k}.\]

If $P$ is integral, we have $k=1$ and its Ehrhart quasipolynomial is a regular polynomial of degree 2. In general, the periods of $a$ and $b$ are divisors of $k$. If they are strictly smaller than $k$, we speak of quasiperiod collapse. A rational polygon is called quasiintegral if the periods of $a$ and $b$ are both 1, hence it has an Ehrhart polynomial.

RationalPolygons.jl comes with many methods for counting the (interior, boundary) lattice points of a rational polygon as well as computing its Ehrhart quasipolynomial and its periods.

Counting lattice points

RationalPolygons.generic_lattice_pointsFunction
generic_lattice_points(P :: RationalPolygon{T,N}, k :: T)

Return the lattice points in the k-fold of the rational polygon P. If the keyword argument interior = true is passed, only the lattice points in the relative interior are computed. If only_count = true is passed, only the total number of lattice points is returned.

source
RationalPolygons.k_rational_hullFunction
k_rational_hull(P :: RationalPolygon{T}, k :: T) where {T <: Integer}

Return the convex hull of all k-rational points contained in P. If primitive = true is passed, only the primitive k-rational points are taken.

source
RationalPolygons.interior_k_rational_hullFunction
interior_k_rational_hull(P :: RationalPolygon{T}, k :: T) where {T <: Integer}

Return the convex hull of all interior k-rational points of P. If primitive = true is passed, only the primitive k-rational points are taken.

source
RationalPolygons.integer_hullFunction
integer_hull(P :: RationalPolygon{T}) where {T <: Integer}

Return the convex hull of all lattice points in the interior of P. If primitive = true is passed, only the primitive interior lattice points are taken.

source
RationalPolygons.interior_integer_hullFunction
interior_integer_hull(P :: RationalPolygon{T}) where {T <: Integer}

Return the convex hull of all interior lattice points of P. If primitive = true is passed, only the primitive lattice points are taken.

source

Ehrhart quasipolynomial

RationalPolygons.is_periodicFunction
is_periodic(v :: Vector, k :: Int)

Check if a vector v is k-periodic, i.e. v[i] == v[mod(i+k,1:n)] for all i, where n = length(v).

source
RationalPolygons.ehrhart_quasipolynomial_with_periodsFunction
ehrhart_quasipolynomial_with_periods(P :: RationalPolygon{T}) where {T <: Integer}

Return a k x 3-matrix of normalized coefficients of the Ehrhart quasipolynomial of a k-rational polygon P, together with a vector of it's three periods.

source

Normal forms and automorphism groups

Two $k$-rational polygons are called (affine) unimodular equivalent if they can be transformed into each other by an (affine) unimodular transformation. The purpose of a normal form is to provide a unique representative for every equivalence class, i.e. two polygons should be (affine) unimodular equivalent to each other if and only if their (affine) unimodular normal forms coincide.

For details about the normal form used in RationalPolygons.jl, see [2].

RationalPolygons.unimodular_normal_formFunction
unimodular_normal_form(P :: RationalPolygon{T,N}) where {N,T <: Integer}

Return a unimodular normal form of a rational polygon. Two rational polygons have the same unimodular normal form if and only if the can be transformed into each other by applying a unimodular transformation.

source
RationalPolygons.affine_normal_formFunction
affine_normal_form(P :: RationalPolygon{T,N}) where {N,T <: Integer}

Return an affine normal form of a rational polygon. Two rational polygons have the same affine normal form if and only if the can be transformed into each other by applying an affine unimodular transformation.

source
RationalPolygons.unimodular_automorphism_groupFunction
unimodular_automorphism_group(P :: RationalPolygon)

Return the automorphism group of P with respect to unimodular transformations.

Example:

julia> P = convex_hull(LatticePoint{Int}[(1,0),(0,1),(-1,1),(-1,0),(0,-1),(1,-1)])
Rational polygon of rationality 1 with 6 vertices.

julia> unimodular_automorphism_group(P)
D6
source

Lattice width

In [3], Bohnert describes the concept of lattice width data, which captures information about the slicing lengths of a polygon with respect to given direction vectors. RationalPolygons.jl implements this concept, following his Definition 2.11.

RationalPolygons.widthFunction
width(P :: RationalPolygon{T}, w :: Point{T}) where {T <: Integer}

Return the lattice width of P in direction w.

Example:

julia> P = convex_hull(LatticePoint{Int}[(1,1),(1,-2),(-4,2),(-2,2)],2)
Rational polygon of rationality 2 with 4 vertices.

julia> width(P, Point(0,1))
2//1

julia> width(P, Point(1,0))
5//2
source
width(P :: RationalPolygon)

Return the lattice width of P.

Example:

julia> P = convex_hull(LatticePoint{Int}[(1,1),(1,-2),(-4,2),(-2,2)],2)
Rational polygon of rationality 2 with 4 vertices.

julia> width(P)
2//1
source
RationalPolygons.all_direction_vectors_with_width_less_thanFunction
all_direction_vectors_with_width_less_than(P :: RationalPolygon{T}, c :: Rational{T}) where {T <: Integer}

Return all direction vectors in which the width of P is less than or equal to a given constant.

Example:

julia> P = convex_hull(LatticePoint{Int}[(1,1),(1,-2),(-4,2),(-2,2)],2)
Rational polygon of rationality 2 with 4 vertices.

julia> all_direction_vectors_with_width_less_than(P, 3//1)
4-element Vector{StaticArraysCore.SVector{2, Int64}}:
 [1, 0]
 [0, 1]
 [1, 1]
 [1, 2]
source
RationalPolygons.width_direction_vectorsFunction
width_direction_vectors(P :: RationalPolygon)

Return the lattice width direction vectors of P, i.e. those directions that realize the lattice width of P.

Example:

julia> P = convex_hull(LatticePoint{Int}[(1,1),(1,-2),(-4,2),(-2,2)],2)
Rational polygon of rationality 2 with 4 vertices.

julia> width_direction_vectors(P)
2-element Vector{StaticArraysCore.SVector{2, Int64}}:
 [0, 1]
 [1, 1]
source
RationalPolygons.adjust_to_width_directionFunction
adjust_to_width_direction(P :: RationalPolygon{T}, w :: Point{T}) where {T <: Integer}

Apply an affine unimodular transformation to P that transforms the given width direction vector to (1,0), see Lemma 2.10 of [3].

source
RationalPolygons.number_of_interior_integral_linesFunction
number_of_interior_integral_lines(P :: RationalPolygon{T}, w :: Point{T}) where {T <: Integer}

Return the number of interior integral lines of P with respect to a given direction vector w.

Example:

julia> P = convex_hull(RationalPoint{Int}[(1//3,-1),(4//3,2),(2//3,2),(-4//3,-1)])
Rational polygon of rationality 3 with 4 vertices.

julia> number_of_interior_integral_lines(P, LatticePoint(1,0))
3

julia> number_of_interior_integral_lines(P, LatticePoint(0,1))
2
source
RationalPolygons.minimal_number_of_interior_integral_linesFunction
minimal_number_of_interior_integral_lines(P :: RationalPolygon{T}) where {T <: Integer}

Return the minimal number of interior integral lines of P.

Example:

julia> P = convex_hull(RationalPoint{Int}[(1//3,-1),(4//3,2),(2//3,2),(-4//3,-1)])
Rational polygon of rationality 3 with 4 vertices.

julia> minimal_number_of_interior_integral_lines(P)
2
source
RationalPolygons.is_realizable_in_intervalFunction
is_realizable_in_interval(P :: RationalPolygon{T}, h :: T) where {T <: Integer}

Check whether P is realizable in $\mathbb{R} \times [0,h]$. This is true if and only if minimal_number_of_interior_integral_lines is less than or equal to $h-1$.

Example:

julia> P = convex_hull(RationalPoint{Int}[(1//3,-1),(4//3,2),(2//3,2),(-4//3,-1)])
Rational polygon of rationality 3 with 4 vertices.

julia> is_realizable_in_interval(P,2)
false

julia> is_realizable_in_interval(P,3)
true
source
RationalPolygons.LatticeWidthDataType
LatticeWidthData{T <: Integer}

A struct capturing information about the lattice width of a rational polygon with respect to a some width direction vector, see Definition 2.11 of [3]. It has two fields:

  • interval_of_nonzero_vertical_slice_length :: Tuple{Rational{T},Rational{T}},
  • interval_of_longest_vertical_slice_length :: Tuple{Rational{T},Rational{T}}.
source
RationalPolygons.lattice_width_dataFunction
lattice_width_data(P :: RationalPolygon{T}, w :: Point{T}) where {T <: Integer}

Compute the lattice width data of a rational polygon with respect to a given width direction vector, see Definition 2.11 of [3]. This function returns a value of type LatticeWidthData.

Example:

julia> P = convex_hull(LatticePoint{Int}[(1,1),(1,-2),(-4,2),(-2,2)],2)
Rational polygon of rationality 2 with 4 vertices.

julia> ws = width_direction_vectors(P)
2-element Vector{StaticArraysCore.SVector{2, Int64}}:
 [0, 1]
 [1, 1]

julia> lattice_width_data(P,ws[1])
LatticeWidthData{Int64}((0//1, 2//1), (3//2, 3//2))

julia> lattice_width_data(P,ws[2])
LatticeWidthData{Int64}((0//1, 2//1), (1//2, 1//2))
source

IO

RationalPolygons.jl provides two ways to save and retrieve polygons from a file. The first is text-based. Polygons can be read from files containing one polygon per line like this:

[[2, 0], [1, 3], [-1, 0], [-3, -4]]
[[1, 0], [2, 6], [-4, -9]]
[[1, 0], [3, 5], [0, 1], [-5, -8]]
....

This text-based format has the advantage of being universally understandable and easy to use. However, storing polygons as ASCII strings is not very space-efficient, as they contain lots of redundant control characters. Hence we provide another way to store polygons in binary and compressed form, which uses the HDF5 file format and is more suitable for large datasets. For an example session, see write_polygon_dataset.

RationalPolygons.parse_rational_polygonsFunction
parse_rational_polygons(k :: T, files :: AbstractVector{String}) where {T <: Integer}
parse_rational_polygons(k :: T, file :: String) where {T <: Integer}

Parse a list of files containing the vertices of a k-rational polygon. The files must have one polygon per line and the vertices must be given as a list of lists of integers, i.e. as in the following example:

[[2, 0], [1, 3], [-1, 0], [-3, -4]]
[[1, 0], [2, 6], [-4, -9]]
[[1, 0], [3, 5], [0, 1], [-5, -8]]
....
source
RationalPolygons.write_rational_polygonsFunction
write_rational_polygons(Ps :: Vector{<:RationalPolygon{T}}, filepath :: String) where {T <: Integer}

Write a list of polygons to a text file. Each polygon will be written as one line containing its vertices, as in the following example:

[[2, 0], [1, 3], [-1, 0], [-3, -4]]
[[1, 0], [2, 6], [-4, -9]]
[[1, 0], [3, 5], [0, 1], [-5, -8]]
....
source
RationalPolygons.create_polygon_datasetFunction
create_polygon_dataset(f :: Union{HDF5.File, HDF5.Group}, path :: String, N :: Int)

Create an HDF5 dataset named path for storing rational polygons with N vertices. The dataset will have an HDF5 compound datatype with 2*N integers, which are the entries of the vertex matrices stored in a column major layout. This function takes three keyword arguments:

  • T :: Type{<:Integer}: The integer type to be used, e.g. Int64, Int32, ... The default is Int.
  • chunk_size :: Int: Chunk size for HDF5. The default is 1000.
  • deflate :: Int: Deflate parameter for HDF5. The default is 3.
source
RationalPolygons.write_polygon_datasetFunction
write_polygon_dataset(f :: Union{HDF5.File, HDF5.Group}, path :: String, Ps :: Vector{RationalPolygon{T,N,M}}) where {N,M,T <: Integer}

Write the polygons Ps to an HDF5 dataset named path. Creates the dataset if it does not exist already. If it does exist, the data will be appended to it.

Example:

Write the reflexive lattice triangles to an HDF5 file:

julia> using RationalPolygons, StaticArrays, HDF5

julia> Vs = SMatrix{2,3}[[1 0 -2 ; 0 1 -1], [1 1 -3 ; 0 2 -4], [1 0 -1 ; 0 1 -1],[1 0 -3 ; 0 1 -2], [1 1 -2 ; 0 3 -3]];

julia> Ps = [RationalPolygon(V,1) for V ∈ Vs];

julia> f = h5open("/tmp/reflexive_triangles.h5", "cw")
🗂️ HDF5.File: (read-write) /tmp/reflexive_triangles.h5

julia> write_polygon_dataset(f, "my_triangles", Ps)

julia> close(f)

The resulting HDF5 file can be inspected using a command line tool such as h5dump:

$ h5dump /tmp/reflexive_triangles.h5 
HDF5 "/tmp/reflexive_triangles.h5" {
GROUP "/" {
   DATASET "my_triangles" {
      DATATYPE  H5T_COMPOUND {
         H5T_STD_I64LE "1";
         H5T_STD_I64LE "2";
         H5T_STD_I64LE "3";
         H5T_STD_I64LE "4";
         H5T_STD_I64LE "5";
         H5T_STD_I64LE "6";
      }
      DATASPACE  SIMPLE { ( 5 ) / ( H5S_UNLIMITED ) }
      DATA {
      (0): {
            1,
            0,
            0,
            1,
            -2,
            -1
         },
      (1): {
            1,
            0,
            1,
            2,
            -3,
            -4
         },
      (2): {
            1,
            0,
            0,
            1,
            -1,
            -1
         },
      (3): {
            1,
            0,
            0,
            1,
            -3,
            -2
         },
      (4): {
            1,
            0,
            1,
            3,
            -2,
            -3
         }
      }
   }
}
}

We can read them back into RationalPolygons.jl at any time using read_polygon_dataset:

julia> using RationalPolygons, HDF5

julia> f = h5open("/tmp/reflexive_triangles.h5", "r")
🗂️ HDF5.File: (read-only) /tmp/reflexive_triangles.h5
└─ 🔢 my_triangles

julia> Ps = read_polygon_dataset(1, f, "my_triangles")
5-element Vector{RationalPolygon{Int64, 3, 6}}:
 Rational polygon of rationality 1 with 3 vertices.
 Rational polygon of rationality 1 with 3 vertices.
 Rational polygon of rationality 1 with 3 vertices.
 Rational polygon of rationality 1 with 3 vertices.
 Rational polygon of rationality 1 with 3 vertices.
source
RationalPolygons.read_polygon_datasetFunction
read_polygon_dataset(k :: T, f :: Union{HDF5.File, HDF5.Group}, path :: String, I...) where {T <: Integer}

Read from an HDF5 dataset containing k-rational polygons.

source