Theta.jl Documentation

Theta.jl Documentation

A Julia package for computing the Riemann theta function and its derivatives.

Overview

The Riemann theta function is the holomorphic function

\[\theta(z,\tau) = \sum_{n\in \mathbb{Z}^g} \exp\left( \pi i n^t \tau n + 2\pi i n^t z \right)\,.\]

where $z\in\mathbb{C}^g$ and $\tau$ belongs to the Siegel upper-half space $\mathbb{H}_g$, which consists of all complex symmetric $g\times g$ matrices with positive definite imaginary part.

A characteristic is a vector of length $2g$ whose entries are 0 or 1, which we write as $\begin{bmatrix}\varepsilon \\\delta\end{bmatrix}$, where $\varepsilon,\delta \in \{0,1\}^g$. The Riemann theta function with characteristic $\begin{bmatrix}\varepsilon \\\delta\end{bmatrix}$ is

\[\theta\begin{bmatrix}\varepsilon \\\delta\end{bmatrix}(z,\tau)\,\,=\,\,\,\sum_{n\in\mathbb{Z}^g}\exp\left[\pi i \left(n+\frac{\varepsilon}{2}\right)^t \tau \left(n+\frac{\varepsilon}{2}\right)+2\pi i\left(n+\frac{\varepsilon}{2}\right)^t\left(z+\frac{\delta}{2}\right)\right]\,.\]

This package computes the Riemann theta function in Julia, with derivatives and characteristics. We optimize the implementation for working with a fixed Riemann matrix $\tau$ and computing the theta function for multiple choices of $z$, as well as multiple derivatives and characteristics. For more details, refer to our upcoming preprint.

Installation

Download Julia 1.1 and above. Start Julia and run

julia> import Pkg
julia> Pkg.add("Theta")

Examples

We start with a matrix M in the Siegel upper-half space.

M = [0.794612+1.9986im 0.815524+1.95836im 0.190195+1.21249im 0.647434+1.66208im 0.820857+1.68942im; 
0.0948191+1.95836im 0.808422+2.66492im 0.857778+1.14274im 0.754323+1.72747im 0.74972+1.95821im; 
0.177874+1.21249im 0.420423+1.14274im 0.445617+1.44248im 0.732018+0.966489im 0.564779+1.57559im; 
0.440969+1.66208im 0.562332+1.72747im 0.292166+0.966489im 0.433763+1.91571im 0.805161+1.46982im; 
0.471487+1.68942im 0.0946854+1.95821im 0.837648+1.57559im 0.311332+1.46982im 0.521253+2.29221im];      

We construct a RiemannMatrix using M.

R = RiemannMatrix(M);

We can then compute the theta function on inputs z and M as follows.

z = [0.30657351+0.34017115im; 0.71945631+0.87045964im; 0.19963849+0.71709398im; 0.64390182+0.97413482im; 0.02747232+0.59071266im];
theta(z, R)

We can also compute first derivatives of theta functions by specifying the direction using the optional argument derivs. The following code computes the partial derivative of the theta function with respect to the first coordinate of z.

theta(z, R, derivs=[[1,0,0,0,0]])

We specify higher order derivatives by adding more elements into the input to derivs, where each element specifies the direction of the derivative. For instance, to compute the partial derivative of the theta function with respect to the first, second and fifth coordinates of z, we run

theta(z, R, derivs=[[1,0,0,0,0], [0,1,0,0,0], [0,0,0,0,1]])

We can compute theta functions with characteristics using the optional argument char.

theta(z, R, char=[[0,1,0,1,1],[0,1,1,0,0]])

We can also compute derivatives of theta functions with characteristics.

theta(z, R, derivs=[[1,0,0,0,0]], char=[[0,1,0,1,1],[0,1,1,0,0]])

Functions

Computing theta functions

Theta.thetaFunction.
theta(z, R, char=[], derivs=[], derivs_t=[])

Compute the Riemann theta function of the matrix in R at z, with optional characteristics or derivatives.

The input derivatives must be unit vectors, and the input z must be of the form $z = a + τb, b∈[0,1)^g, a∈ \mathbb{R}^g$.

Arguments

  • z::Array{<:Number}: Array of size g, each entry is a complex number.
  • R::RiemannMatrix: RiemannMatrix type containing matrix.
  • char::Array{}=[]: 2-element array containing 2 arrays of size g, where entries are either 0 or 1.
  • derivs::Array{}=[]: N-element array specifying the derivative to the N-th order with respect to z, where each element is an array representing a unit vector in the direction of the derivative.
  • derivs_t::Array{}=[]: M-element array specifying the derivative to the M-th order with respect to τ, where each element is an array of size 2 whose entries are the index in τ.

Examples

julia> theta([0,0], RiemannMatrix([1+im -1; -1 1+im]), char=[[1,0],[1,1]], derivs=[[1,0]], derivs_t=[[2,1]])
source

Riemann matrix

RiemannMatrix

Type containing information about the Riemann Matrix, for use in computing the Riemann theta function.

source
Theta.random_siegelFunction.
random_siegel(g)

Sample a random matrix in the Siegel upper half space with genus g.

Arguments

  • g::Integer: genus.

Examples

julia> random_siegel(5)
source
siegel_transform(τ)

Compute Siegel transformation on input Riemann matrix.

Arguments

  • τ::Array{<:Number}: 2-dimensional array of complex numbers

Examples

julia> siegel_transform([1+im -1; -1 1+im])
source
symplectic_transform(Γ, τ)

Compute symplectic group action of Γ on τ.

Arguments

  • Γ::Array{<:Number}: 2-dimensional array of numbers
  • τ::Array{<:Number}: 2-dimensional array of complex numbers

Examples

julia> τ = [1+im -1; -1 1+im];
julia> Γ = siegel_transform(τ)[1];
julia> symplectic_transform(Γ, τ)
source

Theta characteristics

Theta.theta_charFunction.
theta_char(g)

Compute all theta characteristics of genus g.

Arguments

  • g::Integer: genus

Examples

julia> theta_char(4)
source
Theta.even_theta_charFunction.
even_theta_char(g)

Compute all even theta characteristics of genus g.

Arguments

  • g::Integer: genus

Examples

julia> even_theta_char(4)
source
Theta.odd_theta_charFunction.
odd_theta_char(g)

Compute all odd theta characteristics of genus g.

Arguments

  • g::Integer: genus

Examples

julia> odd_theta_char(4)
source
Theta.check_azygeticFunction.
check_azygetic(chars)

Check if a list of characteristics is azygetic.

Arguments

  • chars::Array{}: array where each entry is an array consisting of two arrays of the same length with 0 or 1 entries.

Examples

```julia julia> check_azygetic([[[1,0,1,0], [1,0,1,0]], [[0,0,0,1], [1,0,0,0]], [[0,0,1,1], [1,0,1,1]]])

source

Schottky computations

Schottky problem in genus 4

schottky_genus_4(τ)

Compute the value of the genus 4 Schottky-Igusa polynomial, given as input a 4x4 matrix in the Siegel upper half space.

Arguments

  • τ::Array{<:Number}: 2-dimensional array of complex numbers

Examples

julia> τ = random_siegel(4);
julia> schottky_genus_4(τ)
source
random_nonschottky_genus_4(tol=0.1, trials=100)

Find a random 4x4 matrix in the Siegel upper half space which is not in the Schottky locus, up to the input tolerance and number of trials.

Arguments

  • tol::Real=0.1: tolerance for deciding whether absolute value of Schottky-Igusa polynomial is nonzero.
  • trials::Integer=100: number of trials

Examples

julia> random_nonschottky_genus_4(tol=0.1, trials=10)
source

Schottky problem in genus 5

Theta.accola_charsFunction.
accola_chars()

Compute all characteristics in the 8 Accola special theta relations, for a fixed fundamental system.

Examples

julia> accola_chars()
source
Theta.accolaFunction.
accola(τ, chars=accola_chars())

Compute the 8 Accola special theta relations, without taking the product over the signs of the square roots. Return the largest absolute value of the special theta relation.

Arguments

  • τ::Array{<:Number}: 2-dimensional array of complex numbers
  • chars=accola_chars(): characteristics in the 8 Accola special theta relations

Examples

julia> τ = random_siegel(5);
julia> accola(τ)
source
accola(R, chars=accola_chars())

Compute the 8 Accola special theta relations, without taking the product over the signs of the square roots. Return the largest absolute value of the special theta relation.

Arguments

  • R::RiemannMatrix: RiemannMatrix type containing matrix.
  • chars=accola_chars(): characteristics in the 8 Accola special theta relations

Examples

julia> R = RiemannMatrix(random_siegel(5));
julia> accola(R)
source
random_nonaccola(tol=0.1, trials=100)

Find a random genus 5 matrix in the Siegel upper half space which is not in the Accola locus, such that the largest Accola relation has absolute value at least tol, using input number of trials.

Arguments

  • tol::Real=0.1: tolerance for deciding whether absolute value of largest Accola relation is nonzero.
  • trials::Integer=100: number of trials

Examples

julia> random_nonaccola(tol=0.1, trials=10)
source
Theta.fgsm_charsFunction.
fgsm_chars()

Compute the characteristics used in the FGSM relations in genus 5.

Examples

julia> fgsm_chars()
source
Theta.fgsmFunction.
fgsm(τ, chars=fgsm_chars())

Compute the FGSM relations in genus 5. Return the largest absolute value of the relation.

Arguments

  • τ::Array{<:Number}: 2-dimensional array of complex numbers
  • chars=fgsm_chars(): characteristics in the FGSM relations

Examples

julia> τ = random_siegel(5);
julia> fgsm(τ)
source
fgsm(R, chars=fgsm_chars())

Compute the FGSM relations in genus 5. Return the largest absolute value of the relation.

Arguments

  • R::RiemannMatrix: RiemannMatrix type containing matrix.
  • chars=fgsm_chars(): characteristics in the FGSM relations

Examples

julia> R = RiemannMatrix(random_siegel(5));
julia> fgsm(τ)
source
Theta.random_nonfgsmFunction.
random_nonfgsm(tol=0.1, trials=100)

Find a random genus 5 matrix in the Siegel upper half space which is not in the FGSM locus, such that the largest FGSM relation has absolute value at least tol, using input number of trials.

Arguments

  • tol::Real=0.1: tolerance for deciding whether absolute value of largest FGSM relation is nonzero.
  • trials::Integer=100: number of trials

Examples

julia> random_nonfgsm(tol=0.1, trials=10)
source
min_even_theta_char(R)

Compute the even theta characteristic for which the theta constant at the input matrix in R has the smallest absolute value.

Arguments

  • R::RiemannMatrix: RiemannMatrix type containing matrix.

Examples

julia> R = RiemannMatrix(random_siegel(5));
julia> min_even_theta_char(R)
source
Theta.schottky_nullFunction.
schottky_null(R, tol=1.0e-8)

Compute the even theta characteristic where the theta constant vanishes and the hessian matrix at the characteristic, as well as its nonzero eigenvalues, up to the input tolerance. Returns an array containing the characteristic, the hessian and the nonzero eigenvalues. If there is no such characteristic, return the absolute value of the theta constant at the smallest even characteristic.

Arguments

  • R::RiemannMatrix: RiemannMatrix type containing matrix.
  • tol::Real=1.0e-8: tolerance for deciding whether eigenvalues are nonzero.

Examples

julia> schottky_null(RiemannMatrix([1+im -1; -1 1+im]))
source
schottky_null(τ, tol=1.0e-8)

Compute the even theta characteristic where the theta constant vanishes and the hessian matrix at the characteristic, as well as its nonzero eigenvalues, up to the input tolerance. Returns an array containing the characteristic, the hessian and the nonzero eigenvalues. If there is no such characteristic, return the absolute value of the theta constant at the smallest even characteristic.

Arguments

  • τ::Array{<:Number}: 2-dimensional array of complex numbers
  • tol::Real=1.0e-8: tolerance for deciding whether eigenvalues are nonzero.

Examples

julia> schottky_null([1+im -1; -1 1+im])
source
Theta.hessianFunction.
hessian(z, R, char=[])

Compute the hessian matrix of the theta function at input z and input matrix in R, with optional characteristics.

Arguments

  • z::Array{<:Number}: Array of size g, each entry is a complex number.
  • R::RiemannMatrix: RiemannMatrix type containing matrix.
  • char::Array{}=[]: 2-element array containing 2 arrays of size g, where entries are either 0 or 1.

Examples

julia> hessian([0,0], RiemannMatrix([1+im -1; -1 1+im]), char=[[1,0],[1,1]])
source
hessian(z, τ, char=[], siegel=true, ϵ=1.0e-12)

Compute the hessian matrix of the theta function at inputs z and τ, with optional characteristics, and optional inputs specifying whether to compute the Siegel transformation and the error in the value of the theta function.

Arguments

  • z::Array{<:Number}: Array of size g, each entry is a complex number.
  • τ::Array{<:Number}: 2-dimensional array of complex numbers
  • char::Array{}=[]: 2-element array containing 2 arrays of size g, where entries are either 0 or 1.
  • siegel::Bool=false: do a siegel transformation on τ if true, and use the input matrix otherwise
  • ϵ::Real=1.0e-12: absolute error in value of theta function or its derivatives

Examples

julia> hessian([0,0], [1+im -1; -1 1+im], char=[[1,0],[1,1]])
source

Index