Title: | Jacobi Theta Functions and Related Functions |
---|---|
Description: | Evaluation of the Jacobi theta functions and related functions: Weierstrass elliptic function, Weierstrass sigma function, Weierstrass zeta function, Klein j-function, Dedekind eta function, lambda modular function, Jacobi elliptic functions, Neville theta functions, Eisenstein series, lemniscate elliptic functions, elliptic alpha function, Rogers-Ramanujan continued fractions, and Dixon elliptic functions. Complex values of the variable are supported. |
Authors: | Stéphane Laurent [aut, cre], Mikael Fremling [aut] (author of the original Fortran code for the theta functions) |
Maintainer: | Stéphane Laurent <[email protected]> |
License: | GPL-3 |
Version: | 3.1.1 |
Built: | 2024-11-02 04:43:12 UTC |
Source: | https://github.com/stla/jacobi |
Evaluation of the arithmetic-geometric mean of two complex numbers.
agm(x, y)
agm(x, y)
x , y
|
complex numbers |
A complex number, the arithmetic-geometric mean of x
and y
.
agm(1, sqrt(2)) 2*pi^(3/2)*sqrt(2) / gamma(1/4)^2
agm(1, sqrt(2)) 2*pi^(3/2)*sqrt(2) / gamma(1/4)^2
Evaluation of the amplitude function.
am(u, m)
am(u, m)
u |
complex number |
m |
square of elliptic modulus, a complex number |
A complex number.
library(Carlson) phi <- 1 + 1i m <- 2 u <- elliptic_F(phi, m) am(u, m) # should be phi
library(Carlson) phi <- 1 + 1i m <- 2 u <- elliptic_F(phi, m) am(u, m) # should be phi
Computes a mesh of the Costa surface.
CostaMesh(nu = 50L, nv = 50L)
CostaMesh(nu = 50L, nv = 50L)
nu , nv
|
numbers of subdivisions |
A triangle rgl mesh (object of class mesh3d
).
library(jacobi) library(rgl) mesh <- CostaMesh(nu = 250, nv = 250) open3d(windowRect = c(50, 50, 562, 562), zoom = 0.9) bg3d("#15191E") shade3d(mesh, color = "darkred", back = "cull") shade3d(mesh, color = "orange", front = "cull")
library(jacobi) library(rgl) mesh <- CostaMesh(nu = 250, nv = 250) open3d(windowRect = c(50, 50, 562, 562), zoom = 0.9) bg3d("#15191E") shade3d(mesh, color = "darkred", back = "cull") shade3d(mesh, color = "orange", front = "cull")
Conformal map from the unit disk to the upper half-plane. The function is vectorized.
disk2H(z)
disk2H(z)
z |
a complex number in the unit disk |
A complex number in the upper half-plane.
# map the disk to H and calculate kleinj f <- function(x, y) { z <- complex(real = x, imaginary = y) K <- rep(NA_complex_, length(x)) inDisk <- Mod(z) < 1 K[inDisk] <- kleinj(disk2H(z[inDisk])) K } n <- 1024L x <- y <- seq(-1, 1, length.out = n) Grid <- expand.grid(X = x, Y = y) K <- f(Grid$X, Grid$Y) dim(K) <- c(n, n) # plot if(require("RcppColors")) { img <- colorMap5(K) } else { img <- as.raster(1 - abs(Im(K))/Mod(K)) } opar <- par(mar = c(0, 0, 0, 0)) plot(NULL, xlim = c(0, 1), ylim = c(0, 1), asp = 1, axes = FALSE, xaxs = "i", yaxs = "i", xlab = NA, ylab = NA) rasterImage(img, 0, 0, 1, 1) par(opar)
# map the disk to H and calculate kleinj f <- function(x, y) { z <- complex(real = x, imaginary = y) K <- rep(NA_complex_, length(x)) inDisk <- Mod(z) < 1 K[inDisk] <- kleinj(disk2H(z[inDisk])) K } n <- 1024L x <- y <- seq(-1, 1, length.out = n) Grid <- expand.grid(X = x, Y = y) K <- f(Grid$X, Grid$Y) dim(K) <- c(n, n) # plot if(require("RcppColors")) { img <- colorMap5(K) } else { img <- as.raster(1 - abs(Im(K))/Mod(K)) } opar <- par(mar = c(0, 0, 0, 0)) plot(NULL, xlim = c(0, 1), ylim = c(0, 1), asp = 1, axes = FALSE, xaxs = "i", yaxs = "i", xlab = NA, ylab = NA) rasterImage(img, 0, 0, 1, 1) par(opar)
Conformal map from the unit disk to the square
. The function is vectorized.
disk2square(z)
disk2square(z)
z |
a complex number in the unit disk |
A complex number in the square .
n <- 70L r <- seq(0, 1, length.out = n) theta <- seq(0, 2*pi, length.out = n+1L)[-1L] Grid <- transform( expand.grid(R = r, Theta = theta), Z = R*exp(1i*Theta) ) s <- vapply(Grid$Z, disk2square, complex(1L)) plot(Re(s), Im(s), pch = ".", asp = 1, cex = 2) # # a more insightful plot #### r_ <- seq(0, 1, length.out = 10L) theta_ <- seq(0, 2*pi, length.out = 33)[-1L] plot( NULL, xlim = c(-1, 1), ylim = c(-1, 1), asp = 1, xlab = "x", ylab = "y" ) for(r in r_) { theta <- sort( c(seq(0, 2, length.out = 200L), c(1/4, 3/4, 5/4, 7/4)) ) z <- r*(cospi(theta) + 1i*sinpi(theta)) s <- vapply(z, disk2square, complex(1L)) lines(Re(s), Im(s), col = "blue", lwd = 2) } for(theta in theta_) { r <- seq(0, 1, length.out = 30L) z <- r*exp(1i*theta) s <- vapply(z, disk2square, complex(1L)) lines(Re(s), Im(s), col = "green", lwd = 2) }
n <- 70L r <- seq(0, 1, length.out = n) theta <- seq(0, 2*pi, length.out = n+1L)[-1L] Grid <- transform( expand.grid(R = r, Theta = theta), Z = R*exp(1i*Theta) ) s <- vapply(Grid$Z, disk2square, complex(1L)) plot(Re(s), Im(s), pch = ".", asp = 1, cex = 2) # # a more insightful plot #### r_ <- seq(0, 1, length.out = 10L) theta_ <- seq(0, 2*pi, length.out = 33)[-1L] plot( NULL, xlim = c(-1, 1), ylim = c(-1, 1), asp = 1, xlab = "x", ylab = "y" ) for(r in r_) { theta <- sort( c(seq(0, 2, length.out = 200L), c(1/4, 3/4, 5/4, 7/4)) ) z <- r*(cospi(theta) + 1i*sinpi(theta)) s <- vapply(z, disk2square, complex(1L)) lines(Re(s), Im(s), col = "blue", lwd = 2) } for(theta in theta_) { r <- seq(0, 1, length.out = 30L) z <- r*exp(1i*theta) s <- vapply(z, disk2square, complex(1L)) lines(Re(s), Im(s), col = "green", lwd = 2) }
The Dixon elliptic functions.
sm(z) cm(z)
sm(z) cm(z)
z |
a real or complex number |
A complex number.
# cubic Fermat curve x^3+y^3=1 pi3 <- beta(1/3, 1/3) epsilon <- 0.7 t_ <- seq(-pi3/3 + epsilon, 2*pi3/3 - epsilon, length.out = 100) pts <- t(vapply(t_, function(t) { c(Re(cm(t)), Re(sm(t))) }, FUN.VALUE = numeric(2L))) plot(pts, type = "l", asp = 1)
# cubic Fermat curve x^3+y^3=1 pi3 <- beta(1/3, 1/3) epsilon <- 0.7 t_ <- seq(-pi3/3 + epsilon, 2*pi3/3 - epsilon, length.out = 100) pts <- t(vapply(t_, function(t) { c(Re(cm(t)), Re(sm(t))) }, FUN.VALUE = numeric(2L))) plot(pts, type = "l", asp = 1)
Evaluation of Eisenstein series with weight 2, 4 or 6.
EisensteinE(n, q)
EisensteinE(n, q)
n |
the weight, can be 2, 4 or 6 |
q |
nome, complex number with modulus smaller than one |
A complex number, the value of the Eisenstein series.
Evaluates the elliptic alpha function.
ellipticAlpha(z)
ellipticAlpha(z)
z |
a complex number |
A complex number.
Weisstein, Eric W. "Elliptic Alpha Function".
Elliptic invariants from half-periods.
ellipticInvariants(omega1omega2)
ellipticInvariants(omega1omega2)
omega1omega2 |
the half-periods, a vector of two complex numbers |
The elliptic invariants, a vector of two complex numbers.
Evaluation of the Dedekind eta function.
eta(tau)
eta(tau)
tau |
a vector of complex numbers with strictly positive imaginary parts |
A vector of complex numbers.
eta(2i) gamma(1/4) / 2^(11/8) / pi^(3/4)
eta(2i) gamma(1/4) / 2^(11/8) / pi^(3/4)
Half-periods from elliptic invariants.
halfPeriods(g2g3)
halfPeriods(g2g3)
g2g3 |
the elliptic invariants, a vector of two complex numbers |
The half-periods, a vector of two complex numbers.
Evaluation of the Jacobi elliptic functions.
jellip(kind, u, tau = NULL, m = NULL)
jellip(kind, u, tau = NULL, m = NULL)
kind |
a string with two characters among |
u |
a complex number, vector or matrix |
tau |
complex number with strictly positive imaginary part; it is
related to |
m |
the "parameter", square of the elliptic modulus; it is related to
|
A complex number, vector or matrix.
u <- 2 + 2i tau <- 1i jellip("cn", u, tau)^2 + jellip("sn", u, tau)^2 # should be 1
u <- 2 + 2i tau <- 1i jellip("cn", u, tau)^2 + jellip("sn", u, tau)^2 # should be 1
Evaluates the Jacobi theta function with characteristics.
jtheta_ab(a, b, z, tau = NULL, q = NULL)
jtheta_ab(a, b, z, tau = NULL, q = NULL)
a , b
|
the characteristics, two complex numbers |
z |
complex number, vector, or matrix |
tau |
lattice parameter, a complex number with strictly positive
imaginary part; the two complex numbers |
q |
the nome, a complex number whose modulus is strictly less than one, but not zero |
The Jacobi theta function with characteristics generalizes the four Jacobi
theta functions. It is denoted by
𝜃[a,b](z|τ).
One gets the four Jacobi theta functions when a
and b
take the
values 0
or 0.5
:
a=b=0.5
then one gets -𝜗1(z|τ)
a=0.5
and b=0
then one gets 𝜗2(z|τ)
a=b=0
then one gets 𝜗3(z|τ)
a=0
and b=0.5
then one gets 𝜗4(z|τ)
Both 𝜃[a,b](z+π|τ) and 𝜃[a,b](z+π×τ|τ) are equal to 𝜃[a,b](z|τ) up to a factor - see the examples for the details.
A complex number, vector or matrix, like z
.
Different conventions are used in the book cited as reference.
Hershel M. Farkas, Irwin Kra. Theta Constants, Riemann Surfaces and the Modular Group: An Introduction with Applications to Uniformization Theorems, Partition Identities and Combinatorial Number Theory. Graduate Studies in Mathematics, volume 37, 2001.
a <- 2 + 0.3i b <- 1 - 0.6i z <- 0.1 + 0.4i tau <- 0.2 + 0.3i jab <- jtheta_ab(a, b, z, tau) # first property #### jtheta_ab(a, b, z + pi, tau) # is equal to: jab * exp(2i*pi*a) # second property #### jtheta_ab(a, b, z + pi*tau, tau) # is equal to: jab * exp(-1i*(pi*tau + 2*z + 2*pi*b))
a <- 2 + 0.3i b <- 1 - 0.6i z <- 0.1 + 0.4i tau <- 0.2 + 0.3i jab <- jtheta_ab(a, b, z, tau) # first property #### jtheta_ab(a, b, z + pi, tau) # is equal to: jab * exp(2i*pi*a) # second property #### jtheta_ab(a, b, z + pi*tau, tau) # is equal to: jab * exp(-1i*(pi*tau + 2*z + 2*pi*b))
Evaluates the first Jacobi theta function.
jtheta1(z, tau = NULL, q = NULL) ljtheta1(z, tau = NULL, q = NULL)
jtheta1(z, tau = NULL, q = NULL) ljtheta1(z, tau = NULL, q = NULL)
z |
complex number, vector, or matrix |
tau |
lattice parameter, a complex number with strictly positive
imaginary part; the two complex numbers |
q |
the nome, a complex number whose modulus is strictly less than one, but not zero |
A complex number, vector or matrix; jtheta1
evaluates the
first Jacobi theta function and ljtheta1
evaluates its logarithm.
jtheta1(1 + 1i, q = exp(-pi/2))
jtheta1(1 + 1i, q = exp(-pi/2))
Evaluates the second Jacobi theta function.
jtheta2(z, tau = NULL, q = NULL) ljtheta2(z, tau = NULL, q = NULL)
jtheta2(z, tau = NULL, q = NULL) ljtheta2(z, tau = NULL, q = NULL)
z |
complex number, vector, or matrix |
tau |
lattice parameter, a complex number with strictly positive
imaginary part; the two complex numbers |
q |
the nome, a complex number whose modulus is strictly less than one, but not zero |
A complex number, vector or matrix; jtheta2
evaluates the
second Jacobi theta function and ljtheta2
evaluates its logarithm.
jtheta2(1 + 1i, q = exp(-pi/2))
jtheta2(1 + 1i, q = exp(-pi/2))
Evaluates the third Jacobi theta function.
jtheta3(z, tau = NULL, q = NULL) ljtheta3(z, tau = NULL, q = NULL)
jtheta3(z, tau = NULL, q = NULL) ljtheta3(z, tau = NULL, q = NULL)
z |
complex number, vector, or matrix |
tau |
lattice parameter, a complex number with strictly positive
imaginary part; the two complex numbers |
q |
the nome, a complex number whose modulus is strictly less than one, but not zero |
A complex number, vector or matrix; jtheta3
evaluates the
third Jacobi theta function and ljtheta3
evaluates its logarithm.
jtheta3(1 + 1i, q = exp(-pi/2))
jtheta3(1 + 1i, q = exp(-pi/2))
Evaluates the fourth Jacobi theta function.
jtheta4(z, tau = NULL, q = NULL) ljtheta4(z, tau = NULL, q = NULL)
jtheta4(z, tau = NULL, q = NULL) ljtheta4(z, tau = NULL, q = NULL)
z |
complex number, vector, or matrix |
tau |
lattice parameter, a complex number with strictly positive
imaginary part; the two complex numbers |
q |
the nome, a complex number whose modulus is strictly less than one, but not zero |
A complex number, vector or matrix; jtheta4
evaluates the
fourth Jacobi theta function and ljtheta4
evaluates its logarithm.
jtheta4(1 + 1i, q = exp(-pi/2))
jtheta4(1 + 1i, q = exp(-pi/2))
Evaluation of the Klein j-invariant function and its inverse.
kleinj(tau, transfo = FALSE) kleinjinv(j)
kleinj(tau, transfo = FALSE) kleinjinv(j)
tau |
a complex number with strictly positive imaginary part, or a vector or matrix of such complex numbers; missing values allowed |
transfo |
Boolean, whether to use a transformation of the values
of |
j |
a complex number |
A complex number, vector or matrix.
The Klein-j function is the one with the factor 1728.
( j <- kleinj(2i) ) 66^3 kleinjinv(j)
( j <- kleinj(2i) ) 66^3 kleinjinv(j)
Evaluation of the lambda modular function.
lambda(tau, transfo = FALSE)
lambda(tau, transfo = FALSE)
tau |
a complex number with strictly positive imaginary part, or a vector or matrix of such complex numbers; missing values allowed |
transfo |
Boolean, whether to use a transformation of the values
of |
A complex number, vector or matrix.
The lambda function is the square of the elliptic modulus.
x <- 2 lambda(1i*sqrt(x)) + lambda(1i*sqrt(1/x)) # should be one
x <- 2 lambda(1i*sqrt(x)) + lambda(1i*sqrt(1/x)) # should be one
Lemniscate sine, cosine, arcsine, arccosine, hyperbolic sine, and hyperbolic cosine functions.
sl(z) cl(z) asl(z) acl(z) slh(z) clh(z)
sl(z) cl(z) asl(z) acl(z) slh(z) clh(z)
z |
a real number or a complex number |
A complex number.
sl(1+1i) * cl(1+1i) # should be 1 ## | the lemniscate #### # lemniscate parameterization p <- Vectorize(function(s) { a <- Re(cl(s)) b <- Re(sl(s)) c(a, a * b) / sqrt(1 + b*b) }) # lemnniscate constant ombar <- 2.622 # gamma(1/4)^2 / (2 * sqrt(2*pi)) # plot s_ <- seq(0, ombar, length.out = 100) lemniscate <- t(p(s_)) plot(lemniscate, type = "l", col = "blue", lwd = 3) lines(cbind(lemniscate[, 1L], -lemniscate[, 2L]), col="red", lwd = 3)
sl(1+1i) * cl(1+1i) # should be 1 ## | the lemniscate #### # lemniscate parameterization p <- Vectorize(function(s) { a <- Re(cl(s)) b <- Re(sl(s)) c(a, a * b) / sqrt(1 + b*b) }) # lemnniscate constant ombar <- 2.622 # gamma(1/4)^2 / (2 * sqrt(2*pi)) # plot s_ <- seq(0, ombar, length.out = 100) lemniscate <- t(p(s_)) plot(lemniscate, type = "l", col = "blue", lwd = 3) lines(cbind(lemniscate[, 1L], -lemniscate[, 2L]), col="red", lwd = 3)
The nome in function of the parameter .
nome(m)
nome(m)
m |
the parameter, square of elliptic modulus, real or complex number |
A complex number.
nome(-2)
nome(-2)
Evaluates the Rogers-Ramanujan continued fraction.
RR(q)
RR(q)
q |
the nome, a complex number whose modulus is strictly less than one, and which is not zero |
A complex number
This function is sometimes denoted by .
Evaluates the alternating Rogers-Ramanujan continued fraction.
RRa(q)
RRa(q)
q |
the nome, a complex number whose modulus is strictly less than one, and which is not zero |
A complex number
This function is sometimes denoted by .
Conformal map from the unit square to the unit disk. The function is vectorized.
square2disk(z)
square2disk(z)
z |
a complex number in the unit square |
A complex number in the unit disk.
x <- y <- seq(0, 1, length.out = 25L) Grid <- transform( expand.grid(X = x, Y = y), Z = complex(real = X, imaginary = Y) ) u <- square2disk(Grid$Z) plot(u, pch = 19, asp = 1)
x <- y <- seq(0, 1, length.out = 25L) Grid <- transform( expand.grid(X = x, Y = y), Z = complex(real = X, imaginary = Y) ) u <- square2disk(Grid$Z) plot(u, pch = 19, asp = 1)
Conformal map from the unit square to the upper half-plane. The function is vectorized.
square2H(z)
square2H(z)
z |
a complex number in the unit square |
A complex number in the upper half-plane.
n <- 1024L x <- y <- seq(0.0001, 0.9999, length.out = n) Grid <- transform( expand.grid(X = x, Y = y), Z = complex(real = X, imaginary = Y) ) K <- kleinj(square2H(Grid$Z)) dim(K) <- c(n, n) # plot if(require("RcppColors")) { img <- colorMap5(K) } else { img <- as.raster((Arg(K) + pi)/(2*pi)) } opar <- par(mar = c(0, 0, 0, 0)) plot(NULL, xlim = c(0, 1), ylim = c(0, 1), asp = 1, axes = FALSE, xaxs = "i", yaxs = "i", xlab = NA, ylab = NA) rasterImage(img, 0, 0, 1, 1) par(opar)
n <- 1024L x <- y <- seq(0.0001, 0.9999, length.out = n) Grid <- transform( expand.grid(X = x, Y = y), Z = complex(real = X, imaginary = Y) ) K <- kleinj(square2H(Grid$Z)) dim(K) <- c(n, n) # plot if(require("RcppColors")) { img <- colorMap5(K) } else { img <- as.raster((Arg(K) + pi)/(2*pi)) } opar <- par(mar = c(0, 0, 0, 0)) plot(NULL, xlim = c(0, 1), ylim = c(0, 1), asp = 1, axes = FALSE, xaxs = "i", yaxs = "i", xlab = NA, ylab = NA) rasterImage(img, 0, 0, 1, 1) par(opar)
Evaluation of the Neville theta functions.
theta.s(z, tau = NULL, m = NULL) theta.c(z, tau = NULL, m = NULL) theta.n(z, tau = NULL, m = NULL) theta.d(z, tau = NULL, m = NULL)
theta.s(z, tau = NULL, m = NULL) theta.c(z, tau = NULL, m = NULL) theta.n(z, tau = NULL, m = NULL) theta.d(z, tau = NULL, m = NULL)
z |
a complex number, vector, or matrix |
tau |
complex number with strictly positive imaginary part; it is
related to |
m |
the "parameter", square of the elliptic modulus; it is related to
|
A complex number, vector or matrix.
Evaluation of the Weierstrass elliptic function and its derivatives.
wp(z, g = NULL, omega = NULL, tau = NULL, derivative = 0L)
wp(z, g = NULL, omega = NULL, tau = NULL, derivative = 0L)
z |
complex number, vector or matrix |
g |
the elliptic invariants, a vector of two complex numbers; only
one of |
omega |
the half-periods, a vector of two complex numbers; only
one of |
tau |
the half-periods ratio; supplying |
derivative |
differentiation order, an integer between 0 and 3 |
A complex number, vector or matrix.
omega1 <- 1.4 - 1i omega2 <- 1.6 + 0.5i omega <- c(omega1, omega2) e1 <- wp(omega1, omega = omega) e2 <- wp(omega2, omega = omega) e3 <- wp(-omega1-omega2, omega = omega) e1 + e2 + e3 # should be 0
omega1 <- 1.4 - 1i omega2 <- 1.6 + 0.5i omega <- c(omega1, omega2) e1 <- wp(omega1, omega = omega) e2 <- wp(omega2, omega = omega) e3 <- wp(-omega1-omega2, omega = omega) e1 + e2 + e3 # should be 0
Evaluation of the inverse of the Weierstrass elliptic function.
wpinv(w, g = NULL, omega = NULL, tau = NULL)
wpinv(w, g = NULL, omega = NULL, tau = NULL)
w |
complex number |
g |
the elliptic invariants, a vector of two complex numbers; only
one of |
omega |
the half-periods, a vector of two complex numbers; only
one of |
tau |
the half-periods ratio; supplying |
A complex number.
library(jacobi) omega <- c(1.4 - 1i, 1.6 + 0.5i) w <- 1 + 1i z <- wpinv(w, omega = omega) wp(z, omega = omega) # should be w
library(jacobi) omega <- c(1.4 - 1i, 1.6 + 0.5i) w <- 1 + 1i z <- wpinv(w, omega = omega) wp(z, omega = omega) # should be w
Evaluation of the Weierstrass sigma function.
wsigma(z, g = NULL, omega = NULL, tau = NULL)
wsigma(z, g = NULL, omega = NULL, tau = NULL)
z |
a complex number, vector or matrix |
g |
the elliptic invariants, a vector of two complex numbers; only
one of |
omega |
the half-periods, a vector of two complex numbers; only
one of |
tau |
the half-periods ratio; supplying |
A complex number, vector or matrix.
wsigma(1, g = c(12, -8)) # should be equal to: sin(1i*sqrt(3))/(1i*sqrt(3)) / sqrt(exp(1))
wsigma(1, g = c(12, -8)) # should be equal to: sin(1i*sqrt(3))/(1i*sqrt(3)) / sqrt(exp(1))
Evaluation of the Weierstrass zeta function.
wzeta(z, g = NULL, omega = NULL, tau = NULL)
wzeta(z, g = NULL, omega = NULL, tau = NULL)
z |
complex number, vector or matrix |
g |
the elliptic invariants, a vector of two complex numbers; only
one of |
omega |
the half-periods, a vector of two complex numbers; only
one of |
tau |
the half-periods ratio; supplying |
A complex number, vector or matrix.
# Mirror symmetry property: z <- 1 + 1i g <- c(1i, 1+2i) wzeta(Conj(z), Conj(g)) Conj(wzeta(z, g))
# Mirror symmetry property: z <- 1 + 1i g <- c(1i, 1+2i) wzeta(Conj(z), Conj(g)) Conj(wzeta(z, g))