Rational Polygons
The RationalPolygon
type
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}$ 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.
The standard lattice triangle can be created as follows:
julia> using RationalPolygons, StaticArrays
julia> P = RationalPolygon(SMatrix{2,3}([0 1 0 ; 0 0 1]), 1)
Rational polygon of rationality 1 with 3 vertices.
When creating a RationalPolygon
from a constructor, the user has to be certain that the columns of $V$ truly are vertices of a polygons, i.e. no column is contained in the convex hull of the other columns and the columns are sorted by angle (both clockwise and counterclockwise is allowed). If this is not known, use convex_hull
to create the polygon instead.
There are two ways in which this encoding 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 viewed as a $k$-rational polygon and once viewed as a $2k$-rational polygon. The second way in which this encoding is not unique is that there is no canonical "first vertex" of a polygon, i.e. we can shift the columns of $V$ around and still describe the same polygon. Moreover, we choose to order them clockwise or counterclockwise. This problem is adressed in the section on normal forms.
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 in the plane of affine halfplanes.
RationalPolygons.RationalPolygon
— TypeRationalPolygon{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 byrationality
.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 parameterN
. 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 variable to remember that the polygon is already in unimodular normal form.is_affine_normal_form :: Bool
: A flag variable to remember that the polygon is already in affine normal form.
RationalPolygons.convex_hull
— Functionconvex_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
.
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.
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.
RationalPolygons.intersect_halfplanes
— Functionintersect_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.
RationalPolygons.empty_polygon
— Functionempty_polygon(rationality :: T) where {T <: Integer}
Return the empty polygon of given rationality.
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.
Properties
Basic Properties
RationalPolygons.number_of_vertices
— Functionnumber_of_vertices(P :: RationalPolygon)
Return the number of vertices of a polygon.
RationalPolygons.rationality
— Methodrationality(P :: RationalPolygon)
Return the rationality of P
. Note that this does not need to be 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 2
, even though all vertices are integral.
RationalPolygons.vertex_matrix
— Functionvertex_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
RationalPolygons.is_unimodular_normal_form
— Functionis_unimodular_normal_form(P :: RationalPolygon)
Return whether it is already known that P
is in unimodular normal form.
RationalPolygons.is_affine_normal_form
— Functionis_affine_normal_form(P :: RationalPolygon)
Return whether it is already known that P
is in affine normal form.
RationalPolygons.scaled_vertex
— Functionscaled_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.
RationalPolygons.vertex
— Functionvertex(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.
RationalPolygons.vertices
— Functionvertices(P :: RationalPolygon)
Return the vertices of P
.
RationalPolygons.affine_halfplane
— Methodaffine_halfplane(P :: RationalPolygon, i :: Int)
Return the i
-th describing halfplane of P
.
RationalPolygons.affine_halfplanes
— Methodaffine_halfplanes(P :: RationalPolygon)
Return the describing halfplanes of P
.
Base.in
— MethodBase.in(x :: Point{T}, P :: RationalPolygon{T}) where {T <: Integer}
Check whether a point x
is contained in the polygon P
.
RationalPolygons.contains_in_interior
— Methodcontains_in_interior(x :: Point{T}, P :: RationalPolygon{T}) where {T <: Integer}
Check whether a point x
is contained in the interior of P
.
RationalPolygons.dim
— Functiondim(P :: RationalPolygon)
Return the dimension of P
. For empty polygons, this returns -1. Otherwise, it returns 0, 1 or 2.
RationalPolygons.normalized_area
— Functionnormalized_area(P :: RationalPolygon)
Return the normalized area of a k
-rational polygon, i.e. 2k^2
times its euclidian area. The result is always an integer, counting the number of standard k
-rational triangles contained in P
.
RationalPolygons.euclidian_area
— Functioneuclidian_area(P :: RationalPolygon)
Return the euclidian area of a rational polygon.
RationalPolygons.is_maximal
— Functionis_maximal(P :: RationalPolygon)
Check whether a rational polygon is maximal among all polygons sharing the same rationality and number of interior lattice points.
RationalPolygons.dual
— Functiondual(P :: RationalPolygon{T}) where {T <: Integer}
Return the dual of a polygon P
. Note that P
must contain the origin in its interior.
LDP polygons
An LDP polygon is a lattice polygon with primitive vertices containing the origin in its interior. LDP polygons correspond to toric log del Pezzo surfaces. In this package, we call more generally a k
-rational polygon LDP, if it contains the origin in its interior and its k
-fold dilation has primitive vertices. Here, we list some properties of polygons that are primarily meant to be used for LDP polygons, where they correspond to some invariant of the associated toric del Pezzo surface.
RationalPolygons.contains_origin_in_interior
— Functioncontains_origin_in_interior(P :: RationalPolygon)
Check whether P
contains the origin in its interior.
RationalPolygons.is_primitive
— Functionis_primitive(p :: Point)
Checks whether a point is primitive, i.e. is integral and its coordinates are coprime.
is_primitive(P :: RationalPolygon)
Check whether all vertices of the scaled lattice polygon rationality(P) * P
are primitive.
RationalPolygons.is_ldp
— Functionis_ldp(P :: RationalPolygon)
Check whether P
is a ldp polygon, i.e. is primitive and contains the origin in its interior.
RationalPolygons.multiplicity
— Functionmultiplicity(p :: Point)
The unique rational number x
such that x*p
is primitive and integral.
Example
julia> multiplicity(RationalPoint(4//3,2//3))
3//2
multiplicity(P :: RationalPolygon{T}, i :: Int) where {T <: Integer}
The index of the sublattice spanned by the i
-th and i+1
-th scaled vertex of P
(i.e. the determinant of those two vertices). For ldp polygons, this equals the order of the local class group associated with the toric fixed point associated to the i
-th and i+1
-th ray.
multiplicity(P :: RationalPolygon{T,N}) where {N,T <: Integer}
The order of the sublattice spanned by the scaled vertices of P
. For ldp polygons, this equals the order of the torsion part of the divisors class group of the associated toric surface.
RationalPolygons.picard_index
— Functionpicard_index(P :: RationalPolygon{T,N}) where {N,T <: Integer}
The product of all local multiplicities of P
divided by the global multiplicity. For ldp polygons, this equals the index of the Picard group inside the divisor class group of the associated toric surface, see [2].
RationalPolygons.gorenstein_index
— Functiongorenstein_index(P :: RationalPolygon{T}, i :: Int)
The multiplicity of P
divided by gcd(w[2] - v[2], v[1] - w[1])
, where v
and w
are the i
-th and i+1
-th scaled vertices of P
respectively. For ldp polygons, this equals the local gorenstein at the toric fixed point associated to the i
-th and i+1
-th ray of P
, see e.g. Lemma 3.9 of [3].
gorenstein_index(P :: RationalPolygon{T}) where {T <: Integer}
The least common multiple of the local gorenstein indices of P
. For ldp polygons, this equals the gorenstein index of the associated toric surface.
RationalPolygons.log_canonicity
— Functionlog_canonicity(P :: RationalPolygon)
Given a k
-rational polygon P
, return the maximal rational number 0 < ϵ ≤ 1 such that ε*P contains only one k
-rational point in its interior (the origin). For an ldp polygon, this equals the maximal rational number 0 < ε ≤ 1 such that the associated toric surface is ε-log canonical.
RationalPolygons.toric_prime_divisor_self_intersection
— Functiontoric_prime_divisor_self_intersection(P :: RationalPolygon, i :: Int)
Writing u
, v
and w
for the i-1
-th, i
-th and i+1
-th scaled vertex of P
respectively, return det(w,u) // (det(u,v) * det(v,w))
. For ldp polygons, this equals the self intersection number of the i
-th toric prime divisor, see e.g. Summary 3.2 of [4].
RationalPolygons.toric_prime_divisor_adjacent_intersection
— Functiontoric_prime_divisor_adjacent_intersection(P :: RationalPolygon, i :: Int)
Writing v
and w
for the i
-th and i+1
-th scaled vertex of P
, return 1 // det(v,w)
. For ldp polygons, this equals the intersection number between the i
-th and i+1
-th toric prime divisors.
RationalPolygons.degree
— Functiondegree(P :: RationalPolygon)
For ldp polygons, return the self intersection number of an anticanonocal divisor of the associated toric surface.
Example
The projective plane has degree 9.
julia> P = convex_hull(LatticePoint{Int}[(1,0),(0,1),(-1,-1)])
Rational polygon of rationality 1 with 3 vertices.
julia> degree(P)
9//1
Ehrhart Theory
Consider a $k$-rational polygon $P$. The main result of Ehrhart Theory is that the the number of lattice points in integral multiples of $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 euclidian 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. We then 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 polyomial 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_points
— Functiongeneric_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.
RationalPolygons.boundary_k_rational_points
— Functionboundary_k_rational_points(P :: RationalPolygon{T,N}, k :: T) where {N,T <: Integer}
Return all k
-rational points on the boundary of P
.
RationalPolygons.number_of_boundary_k_rational_points
— Functionnumber_of_boundary_k_rational_points(P :: RationalPolygon{T}, k :: T) where {T <: Integer}
Return the number of k
-rational points on the boundary of P
.
RationalPolygons.boundary_lattice_points
— Functionboundary_lattice_points(P :: RationalPolygon{T}) where {T <: Integer}
Return all lattice points on the boundary of P
.
RationalPolygons.number_of_boundary_lattice_points
— Functionnumber_of_boundary_lattice_points(P :: RationalPolygon)
Return the number of lattice points on the boundary of P
.
RationalPolygons.interior_k_rational_points
— Functioninterior_k_rational_points(P :: RationalPolygon{T}, k :: T) where {T <: Integer}
Return all k
-rational points in the interior of P
.
RationalPolygons.number_of_interior_k_rational_points
— Functionnumber_of_interior_k_rational_points(k :: T, P :: RationalPolygon{T}) where {T <: Integer}
Return the number of k
-rational points in the interior of P
.
RationalPolygons.interior_lattice_points
— Functioninterior_lattice_points(P :: RationalPolygon{T}) where {T <: Integer}
Return all lattice points in the interior of P
.
RationalPolygons.number_of_interior_lattice_points
— Functionnumber_of_interior_lattice_points(P :: RationalPolygon{T}) where {T <: Integer}
Return the number of lattice points in the interior of P
.
RationalPolygons.k_rational_points
— Functionk_rational_points(k :: T, P :: RationalPolygon{T}) where {T <: Integer}
Return all k
-rational points in P
.
RationalPolygons.number_of_k_rational_points
— Functionnumber_of_k_rational_points(P :: RationalPolygon{T}, k :: T) where {T <: Integer}
Return the number of k
-rational points in P
.
RationalPolygons.lattice_points
— Functionlattice_points(P :: RationalPolygon{T}) where {T <: Integer}
Return all lattice points in P
.
RationalPolygons.number_of_lattice_points
— Functionnumber_of_lattice_points(P :: RationalPolygon{T}) where {T <: Integer}
Return the number of lattice points in P
.
RationalPolygons.k_rational_hull
— Functionk_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.
RationalPolygons.interior_k_rational_hull
— Functioninterior_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.
RationalPolygons.integer_hull
— Functioninteger_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.
RationalPolygons.interior_integer_hull
— Functioninterior_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.
Ehrhart quasipolynomial
RationalPolygons.is_periodic
— Functionis_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)
.
RationalPolygons.period
— Functionperiod(v :: Vector)
Return the smallest positive integer k
such that v
is k
-periodic.
RationalPolygons.ehrhart_quasipolynomial_with_periods
— Functionehrhart_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.
RationalPolygons.ehrhart_quasipolynomial
— Functionehrhart_quasipolynomial(P :: RationalPolygon)
Return a k x 3-matrix of normalized coefficients of the Ehrhart quasipolynomial of a k
-rational polygon P
.
RationalPolygons.ehrhart_quasipolynomial_periods
— Functionehrhart_quasipolynomial_periods(P :: RationalPolygon)
Return the periods of the Ehrhart quasipolynomial of a k
-rational polygon P
.
RationalPolygons.ehrhart_quasipolynomial_period
— Functionehrhart_quasipolynomial_period(P :: RationalPolygon)
Return the period of the Ehrhart quasipolynomial of a k
-rational polygon P
.
RationalPolygons.is_quasiintegral
— Functionis_quasiintegral(P :: RationalPolygon)
Check whether a rational polygon is quasiintegral, i.e. has an Ehrhart polynomial.
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, we refer to [1].
Normal forms
RationalPolygons.unimodular_normal_form
— Functionunimodular_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.
RationalPolygons.are_unimodular_equivalent
— Functionare_unimodular_equivalent(P :: RationalPolygon, Q :: RationalPolygon)
Checks whether two rational polygons are equivalent by a unimodular transformation.
RationalPolygons.affine_normal_form
— Functionaffine_normal_form(P :: RationalPolygon{T,N}) where {N,T <: Integer}
Return a 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.
RationalPolygons.are_affine_equivalent
— Functionare_affine_equivalent(P :: RationalPolygon, Q :: RationalPolygon)
Checks whether two rational polygons are equivalent by an affine unimodular transformation.
Automorphism groups
RationalPolygons.PolygonAutomorphismGroup
— TypePolygonAutomorphismGroup
A struct holding information about the automorphism group of a rational polygon. It has two fields is_cyclic :: Bool
and n :: Int
.
RationalPolygons.CyclicGroup
— FunctionCyclicGroup(n :: Int)
Return the cyclic group of order n
, as an PolygonAutomorphismGroup
.
RationalPolygons.DihedralGroup
— FunctionDihedralGroup(n :: Int)
Return the dihedral group of order 2n
, as an PolygonAutomorphismGroup
.
RationalPolygons.is_cyclic
— Functionis_cyclic(G :: PolygonAutomorphismGroup)
Check whether G
is a cyclic group.
RationalPolygons.order
— Functionorder(G :: PolygonAutomorphismGroup)
Return the order of the group G
.
RationalPolygons.unimodular_automorphism_group
— Functionunimodular_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
RationalPolygons.affine_automorphism_group
— Functionaffine_automorphism_group(P :: RationalPolygon{T,N}) where {N,T <: Integer}
Return the automorphism group of P
with respect to affine unimodular transformations.
Lattice width
In [5], Bohnert describes the concept of lattice width data, which captures information about the slicing lengths of a polygon with respect to a given direction vectors. RationalPolygons.jl implements this concept, following his Definition 2.11.
RationalPolygons.width
— Functionwidth(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
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
RationalPolygons.all_direction_vectors_with_width_less_than
— Functionall_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]
RationalPolygons.width_direction_vectors
— Functionwidth_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]
RationalPolygons.adjust_to_width_direction
— Functionadjust_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 [5].
RationalPolygons.number_of_interior_integral_lines
— Functionnumber_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
RationalPolygons.minimal_number_of_interior_integral_lines
— Functionminimal_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
RationalPolygons.is_realizable_in_interval
— Functionis_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
RationalPolygons.LatticeWidthData
— TypeLatticeWidthData{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 [5]. 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}}
.
RationalPolygons.lattice_width_data
— Functionlattice_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 [5]. 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))
RationalPolygons.number_of_interior_integral_vertical_lines
— Functionnumber_of_interior_integral_vertical_lines(P :: RationalPolygon, w :: Point{T}) where {T <: Integer}
Return the number of interior integral vertical lines of P
with respect to a given lattice width direction vector, see Definition 2.11 of [5].
RationalPolygons.position_of_longest_vertical_slice_length
— Functionposition_of_longest_vertical_slice_length(P :: RationalPolygon, w :: Point{T}) where {T <: Integer}
Return the position of the longest vertical slicing length of P
with respect to a given lattice width direction vector, see Definition 2.11 of [5].
RationalPolygons.lattice_width_datas
— Functionlattice_width_datas(P :: RationalPolygon{T}) where {T <: Integer}
Return the lattice with datas for all lattice width direction vectors of P
.
RationalPolygons.numbers_of_interior_integral_vertical_lines
— Functionnumbers_of_interior_integral_vertical_lines(P :: RationalPolygon{T}) where {T <: Integer}
Return the number of interior integral vertical lines for all width direction vectors of P
.
RationalPolygons.positions_of_longest_vertical_slice_length
— Functionpositions_of_longest_vertical_slice_length(P :: RationalPolygon{T}) where {T <: Integer}
Return the positions of the longest vertical slicing lengths for all width direction vectors of P
.
IO
RationalPolygons.jl provides two ways to save polygons to a file: The first is text-based, where polygons can be written and read to 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-efficients, as they contain lots of redundant control characters. Hence we provide another way to store polygons in binary form, which uses the HDF5 format and is more suitable for large datasets. For an example session, see write_polygon_dataset
.
RationalPolygons.parse_rational_polygons
— Functionparse_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]] ....
RationalPolygons.write_rational_polygons
— Functionwrite_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]] ....
RationalPolygons.create_polygon_dataset
— Functioncreate_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 isInt
.chunk_size :: Int
: Chunk size for HDF5. The default is1000
.deflate :: Int
: Deflate parameter for HDF5. The default is3
.
RationalPolygons.write_polygon_dataset
— Functionwrite_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.
RationalPolygons.read_polygon_dataset
— Functionread_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.