elliptic-package {elliptic}R Documentation

Weierstrass and Jacobi Elliptic Functions

Description

A suite of elliptic and related functions including Weierstrass and Jacobi forms. Also includes various tools for manipulating and visualizing complex functions.

Details

The DESCRIPTION file:

Package: elliptic
Version: 1.4-0
Title: Weierstrass and Jacobi Elliptic Functions
Authors@R: person(given=c("Robin", "K. S."), family="Hankin", role = c("aut","cre"), email="hankin.robin@gmail.com", comment = c(ORCID = "0000-0001-5982-0415"))
Depends: R (>= 2.5.0)
Imports: MASS
Suggests: emulator, calibrator (>= 1.2-8)
SystemRequirements: pari/gp
Description: A suite of elliptic and related functions including Weierstrass and Jacobi forms. Also includes various tools for manipulating and visualizing complex functions.
Maintainer: Robin K. S. Hankin <hankin.robin@gmail.com>
License: GPL-2
URL: https://github.com/RobinHankin/elliptic.git
BugReports: https://github.com/RobinHankin/elliptic/issues
Author: Robin K. S. Hankin [aut, cre] (<https://orcid.org/0000-0001-5982-0415>)

Index of help topics:

Im<-                    Manipulate real or imaginary components of an
                        object
J                       Various modular functions
K.fun                   quarter period K
P.laurent               Laurent series for elliptic and related
                        functions
WeierstrassP            Weierstrass P and related functions
amn                     matrix a on page 637
as.primitive            Converts basic periods to a primitive pair
ck                      Coefficients of Laurent expansion of
                        Weierstrass P function
congruence              Solves mx+by=1 for x and y
coqueraux               Fast, conceptually simple, iterative scheme for
                        Weierstrass P functions
divisor                 Number theoretic functions
e16.28.1                Numerical verification of equations 16.28.1 to
                        16.28.5
e18.10.9                Numerical checks of equations 18.10.9-11, page
                        650
e1e2e3                  Calculate e1, e2, e3 from the invariants
elliptic-package        Weierstrass and Jacobi Elliptic Functions
equianharmonic          Special cases of the Weierstrass elliptic
                        function
eta                     Dedekind's eta function
farey                   Farey sequences
fpp                     Fundamental period parallelogram
g.fun                   Calculates the invariants g2 and g3
half.periods            Calculates half periods in terms of e
latplot                 Plots a lattice of periods on the complex plane
lattice                 Lattice of complex numbers
limit                   Limit the magnitude of elements of a vector
massage                 Massages numbers near the real line to be real
mob                     Moebius transformations
myintegrate             Complex integration
near.match              Are two vectors close to one another?
newton_raphson          Newton Raphson iteration to find roots of
                        equations
nome                    Nome in terms of m or k
p1.tau                  Does the right thing when calling g2.fun() and
                        g3.fun()
parameters              Parameters for Weierstrass's P function
pari                    Wrappers for PARI functions
sn                      Jacobi form of the elliptic functions
sqrti                   Generalized square root
theta                   Jacobi theta functions 1-4
theta.neville           Neville's form for the theta functions
theta1.dash.zero        Derivative of theta1
theta1dash              Derivatives of theta functions
unimodular              Unimodular matrices
view                    Visualization of complex functions

The primary function in package elliptic is P(): this calculates the Weierstrass \wp function, and may take named arguments that specify either the invariants g or half periods Omega. The derivative is given by function Pdash and the Weierstrass sigma and zeta functions are given by functions sigma() and zeta() respectively; these are documented in ?P. Jacobi forms are documented under ?sn and modular forms under ?J.

Notation follows Abramowitz and Stegun (1965) where possible, although there only real invariants are considered; ?e1e2e3 and ?parameters give a more detailed discussion. Various equations from AMS-55 are implemented (for fun); the functions are named after their equation numbers in AMS-55; all references are to this work unless otherwise indicated.

The package uses Jacobi's theta functions (?theta and ?theta.neville) where possible: they converge very quickly.

Various number-theoretic functions that are required for (eg) converting a period pair to primitive form (?as.primitive) are implemented; see ?divisor for a list.

The package also provides some tools for numerical verification of complex analysis such as contour integration (?myintegrate) and Newton-Raphson iteration for complex functions (?newton_raphson).

Complex functions may be visualized using view(); this is customizable but has an extensive set of built-in colourmaps.

Author(s)

NA

Maintainer: Robin K. S. Hankin <hankin.robin@gmail.com>

References

Examples

     ## Example 8, p666, RHS:
 P(z=0.07 + 0.1i, g=c(10,2)) 

     ## Now a nice little plot of the zeta function:
 x <- seq(from=-4,to=4,len=100)
 z <- outer(x,1i*x,"+")
 par(pty="s")
 view(x,x,limit(zeta(z,c(1+1i,2-3i))),nlevels=3,scheme=1)
 view(x,x,P(z*3,params=equianharmonic()),real=FALSE)

     ## Some number theory:
 mobius(1:10)
 plot(divisor(1:300,k=1),type="s",xlab="n",ylab="divisor(n,1)")

    ## Primitive periods:
 as.primitive(c(3+4.01i , 7+10i))
 as.primitive(c(3+4.01i , 7+10i),n=10)   # Note difference

    ## Now some contour integration:
 f <- function(z){1/z}
 u <- function(x){exp(2i*pi*x)}
 udash <- function(x){2i*pi*exp(2i*pi*x)}
 integrate.contour(f,u,udash) - 2*pi*1i


 x <- seq(from=-10,to=10,len=200)
 z <- outer(x,1i*x,"+")
 view(x,x,P(z,params=lemniscatic()),real=FALSE)
 view(x,x,P(z,params=pseudolemniscatic()),real=FALSE)
 view(x,x,P(z,params=equianharmonic()),real=FALSE)

[Package elliptic version 1.4-0 Index]