Riemann Mapping#

AUTHORS:

  • Ethan Van Andel (2009-2011): initial version and upgrades

  • Robert Bradshaw (2009): his “complex_plot” was adapted for plot_colored

Development supported by NSF award No. 0702939.

class sage.calculus.riemann.Riemann_Map#

Bases: object

The Riemann_Map class computes an interior or exterior Riemann map, or an Ahlfors map of a region given by the supplied boundary curve(s) and center point. The class also provides various methods to evaluate, visualize, or extract data from the map.

A Riemann map conformally maps a simply connected region in the complex plane to the unit disc. The Ahlfors map does the same thing for multiply connected regions.

Note that all the methods are numerical. As a result all answers have some imprecision. Moreover, maps computed with small number of collocation points, or for unusually shaped regions, may be very inaccurate. Error computations for the ellipse can be found in the documentation for analytic_boundary() and analytic_interior().

[BSV2010] provides an overview of the Riemann map and discusses the research that lead to the creation of this module.

INPUT:

  • fs – A list of the boundaries of the region, given as complex-valued functions with domain \(0\) to \(2*pi\). Note that the outer boundary must be parameterized counter clockwise (i.e. e^(I*t)) while the inner boundaries must be clockwise (i.e. e^(-I*t)).

  • fprimes – A list of the derivatives of the boundary functions. Must be in the same order as fs.

  • a – Complex, the center of the Riemann map. Will be mapped to the origin of the unit disc. Note that a MUST be within the region in order for the results to be mathematically valid.

The following inputs may be passed in as named parameters:

  • N – integer (default: 500), the number of collocation points used to compute the map. More points will give more accurate results, especially near the boundaries, but will take longer to compute.

  • exterior – boolean (default: False), if set to True, the exterior map will be computed, mapping the exterior of the region to the exterior of the unit circle.

The following inputs may be passed as named parameters in unusual circumstances:

  • ncorners – integer (default: 4), if mapping a figure with (equally t-spaced) corners – corners that make a significant change in the direction of the boundary – better results may be sometimes obtained by accurately giving this parameter. Used to add the proper constant to the theta correspondence function.

  • opp – boolean (default: False), set to True in very rare cases where the theta correspondence function is off by pi, that is, if red is mapped left of the origin in the color plot.

EXAMPLES:

The unit circle identity map:

sage: f(t) = e^(I*t)
sage: fprime(t) = I*e^(I*t)
sage: m = Riemann_Map([f], [fprime], 0)  # long time (4 sec)
sage: m.plot_colored() + m.plot_spiderweb()  # long time
Graphics object consisting of 22 graphics primitives

The exterior map for the unit circle:

sage: m = Riemann_Map([f], [fprime], 0, exterior=True)  # long time (4 sec)
sage: #spiderwebs are not supported for exterior maps
sage: m.plot_colored() # long time
Graphics object consisting of 1 graphics primitive

The unit circle with a small hole:

sage: f(t) = e^(I*t)
sage: fprime(t) = I*e^(I*t)
sage: hf(t) = 0.5*e^(-I*t)
sage: hfprime(t) = 0.5*-I*e^(-I*t)
sage: m = Riemann_Map([f, hf], [fprime, hfprime], 0.5 + 0.5*I)
sage: #spiderweb and color plots cannot be added for multiply
sage: #connected regions. Instead we do this.
sage: m.plot_spiderweb(withcolor = True)  # long time
Graphics object consisting of 3 graphics primitives

A square:

sage: ps = polygon_spline([(-1, -1), (1, -1), (1, 1), (-1, 1)])
sage: f = lambda t: ps.value(real(t))
sage: fprime = lambda t: ps.derivative(real(t))
sage: m = Riemann_Map([f], [fprime], 0.25, ncorners=4)
sage: m.plot_colored() + m.plot_spiderweb()  # long time
Graphics object consisting of 22 graphics primitives

Compute rough error for this map:

sage: x = 0.75  # long time
sage: print("error = {}".format(m.inverse_riemann_map(m.riemann_map(x)) - x))  # long time
error = (-0.000...+0.0016...j)

A fun, complex region for demonstration purposes:

sage: f(t) = e^(I*t)
sage: fp(t) = I*e^(I*t)
sage: ef1(t) = .2*e^(-I*t) +.4+.4*I
sage: ef1p(t) = -I*.2*e^(-I*t)
sage: ef2(t) = .2*e^(-I*t) -.4+.4*I
sage: ef2p(t) = -I*.2*e^(-I*t)
sage: pts = [(-.5,-.15-20/1000),(-.6,-.27-10/1000),(-.45,-.45),(0,-.65+10/1000),(.45,-.45),(.6,-.27-10/1000),(.5,-.15-10/1000),(0,-.43+10/1000)]
sage: pts.reverse()
sage: cs = complex_cubic_spline(pts)
sage: mf = lambda x:cs.value(x)
sage: mfprime = lambda x: cs.derivative(x)
sage: m = Riemann_Map([f,ef1,ef2,mf],[fp,ef1p,ef2p,mfprime],0,N = 400) # long time
sage: p = m.plot_colored(plot_points = 400) # long time

ALGORITHM:

This class computes the Riemann Map via the Szego kernel using an adaptation of the method described by [KT1986].

compute_on_grid(plot_range, x_points)#

Compute the Riemann map on a grid of points.

Note that these points are complex of the form z = x + y*i.

INPUT:

  • plot_range – a tuple of the form [xmin, xmax, ymin, ymax]. If the value is [], the default plotting window of the map will be used.

  • x_points – int, the size of the grid in the x direction The number of points in the y_direction is scaled accordingly

OUTPUT:

  • a tuple containing [z_values, xmin, xmax, ymin, ymax] where z_values is the evaluation of the map on the specified grid.

EXAMPLES:

General usage:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)
sage: data = m.compute_on_grid([],5)
sage: data[0][8,1]
(-0.0879...+0.9709...j)
get_szego(boundary=-1, absolute_value=False)#

Return a discretized version of the Szego kernel for each boundary function.

INPUT:

The following inputs may be passed in as named parameters:

  • boundary – integer (default: -1) if < 0, get_theta_points() will return the points for all boundaries. If >= 0, get_theta_points() will return only the points for the boundary specified.

  • absolute_value – boolean (default: False) if True, will return the absolute value of the (complex valued) Szego kernel instead of the kernel itself. Useful for plotting.

OUTPUT:

A list of points of the form [t value, value of the Szego kernel at that t].

EXAMPLES:

Generic use:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)
sage: sz = m.get_szego(boundary=0)
sage: points = m.get_szego(absolute_value=True)
sage: list_plot(points)
Graphics object consisting of 1 graphics primitive

Extending the points by a spline:

sage: s = spline(points)
sage: s(3*pi / 4)
0.0012158...
sage: plot(s,0,2*pi) # plot the kernel
Graphics object consisting of 1 graphics primitive

The unit circle with a small hole:

sage: f(t) = e^(I*t)
sage: fprime(t) = I*e^(I*t)
sage: hf(t) = 0.5*e^(-I*t)
sage: hfprime(t) = 0.5*-I*e^(-I*t)
sage: m = Riemann_Map([f, hf], [fprime, hfprime], 0.5 + 0.5*I)

Getting the szego for a specific boundary:

sage: sz0 = m.get_szego(boundary=0)
sage: sz1 = m.get_szego(boundary=1)
get_theta_points(boundary=-1)#

Return an array of points of the form [t value, theta in e^(I*theta)], that is, a discretized version of the theta/boundary correspondence function. In other words, a point in this array [t1, t2] represents that the boundary point given by f(t1) is mapped to a point on the boundary of the unit circle given by e^(I*t2).

For multiply connected domains, get_theta_points will list the points for each boundary in the order that they were supplied.

INPUT:

The following input must all be passed in as named parameters:

  • boundary – integer (default: -1) if < 0, get_theta_points() will return the points for all boundaries. If >= 0, get_theta_points() will return only the points for the boundary specified.

OUTPUT:

A list of points of the form [t value, theta in e^(I*theta)].

EXAMPLES:

Getting the list of points, extending it via a spline, getting the points for only the outside of a multiply connected domain:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)
sage: points = m.get_theta_points()
sage: list_plot(points)
Graphics object consisting of 1 graphics primitive

Extending the points by a spline:

sage: s = spline(points)
sage: s(3*pi / 4)
1.627660...

The unit circle with a small hole:

sage: f(t) = e^(I*t)
sage: fprime(t) = I*e^(I*t)
sage: hf(t) = 0.5*e^(-I*t)
sage: hfprime(t) = 0.5*-I*e^(-I*t)
sage: m = Riemann_Map([f, hf], [hf, hfprime], 0.5 + 0.5*I)

Getting the boundary correspondence for a specific boundary:

sage: tp0 = m.get_theta_points(boundary=0)
sage: tp1 = m.get_theta_points(boundary=1)
inverse_riemann_map(pt)#

Return the inverse Riemann mapping of a point.

That is, given pt on the interior of the unit disc, inverse_riemann_map() will return the point on the original region that would be Riemann mapped to pt. Note that this method does not work for multiply connected domains.

INPUT:

  • pt – A complex number (usually with absolute value <= 1) representing the point to be inverse mapped.

OUTPUT:

The point on the region that Riemann maps to the input point.

EXAMPLES:

Can work for different types of complex numbers:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)
sage: m.inverse_riemann_map(0.5 + sqrt(-0.5))
(0.406880...+0.3614702...j)
sage: m.inverse_riemann_map(0.95)
(0.486319...-4.90019052...j)
sage: m.inverse_riemann_map(0.25 - 0.3*I)
(0.1653244...-0.180936...j)
sage: m.inverse_riemann_map(complex(-0.2, 0.5))
(-0.156280...+0.321819...j)
plot_boundaries(plotjoined=True, rgbcolor=[0, 0, 0], thickness=1)#

Plots the boundaries of the region for the Riemann map. Note that this method DOES work for multiply connected domains.

INPUT:

The following inputs may be passed in as named parameters:

  • plotjoined – boolean (default: True) If False, discrete points will be drawn; otherwise they will be connected by lines. In this case, if plotjoined=False, the points shown will be the original collocation points used to generate the Riemann map.

  • rgbcolor – float array (default: [0,0,0]) the red-green-blue color of the boundary.

  • thickness – positive float (default: 1) the thickness of the lines or points in the boundary.

EXAMPLES:

General usage:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)

Default plot:

sage: m.plot_boundaries()
Graphics object consisting of 1 graphics primitive

Big blue collocation points:

sage: m.plot_boundaries(plotjoined=False, rgbcolor=[0,0,1], thickness=6)
Graphics object consisting of 1 graphics primitive
plot_colored(plot_range=[], plot_points=100, interpolation='catrom', **options)#

Generates a colored plot of the Riemann map. A red point on the colored plot corresponds to a red point on the unit disc.

INPUT:

The following inputs may be passed in as named parameters:

  • plot_range – (default: []) list of 4 values (xmin, xmax, ymin, ymax). Declare if you do not want the plot to use the default range for the figure.

  • plot_points – integer (default: 100), number of points to plot in the x direction. Points in the y direction are scaled accordingly. Note that very large values can cause this function to run slowly.

EXAMPLES:

Given a Riemann map m, general usage:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)
sage: m.plot_colored()
Graphics object consisting of 1 graphics primitive

Plot zoomed in on a specific spot:

sage: m.plot_colored(plot_range=[0,1,.25,.75])
Graphics object consisting of 1 graphics primitive

High resolution plot:

sage: m.plot_colored(plot_points=1000)  # long time (29s on sage.math, 2012)
Graphics object consisting of 1 graphics primitive

To generate the unit circle map, it’s helpful to see what the colors correspond to:

sage: f(t) = e^(I*t)
sage: fprime(t) = I*e^(I*t)
sage: m = Riemann_Map([f], [fprime], 0, 1000)
sage: m.plot_colored()
Graphics object consisting of 1 graphics primitive
plot_spiderweb(spokes=16, circles=4, pts=32, linescale=0.99, rgbcolor=[0, 0, 0], thickness=1, plotjoined=True, withcolor=False, plot_points=200, min_mag=0.001, interpolation='catrom', **options)#

Generate a traditional “spiderweb plot” of the Riemann map.

This shows what concentric circles and radial lines map to. The radial lines may exhibit erratic behavior near the boundary; if this occurs, decreasing linescale may mitigate the problem.

For multiply connected domains the spiderweb is by necessity generated using the forward mapping. This method is more computationally intensive. In addition, these spiderwebs cannot be added to color plots. Instead the withcolor option must be used.

In addition, spiderweb plots are not currently supported for exterior maps.

INPUT:

The following inputs may be passed in as named parameters:

  • spokes – integer (default: 16) the number of equally spaced radial lines to plot.

  • circles – integer (default: 4) the number of equally spaced circles about the center to plot.

  • pts – integer (default: 32) the number of points to plot. Each radial line is made by 1*pts points, each circle has 2*pts points. Note that high values may cause erratic behavior of the radial lines near the boundaries. - only for simply connected domains

  • linescale – float between 0 and 1. Shrinks the radial lines away from the boundary to reduce erratic behavior. - only for simply connected domains

  • rgbcolor – float array (default: [0,0,0]) the red-green-blue color of the spiderweb.

  • thickness – positive float (default: 1) the thickness of the lines or points in the spiderweb.

  • plotjoined – boolean (default: True) If False, discrete points will be drawn; otherwise they will be connected by lines. - only for simply connected domains

  • withcolor – boolean (default: False) If True, The spiderweb will be overlaid on the basic color plot.

  • plot_points – integer (default: 200) the size of the grid in the x direction The number of points in the y_direction is scaled accordingly. Note that very large values can cause this function to run slowly. - only for multiply connected domains

  • min_mag – float (default: 0.001) The magnitude cutoff below which spiderweb points are not drawn. This only applies to multiply connected domains and is designed to prevent “fuzz” at the edge of the domain. Some complicated multiply connected domains (particularly those with corners) may require a larger value to look clean outside.

EXAMPLES:

General usage:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)

Default plot:

sage: m.plot_spiderweb()
Graphics object consisting of 21 graphics primitives

Simplified plot with many discrete points:

sage: m.plot_spiderweb(spokes=4, circles=1, pts=400, linescale=0.95, plotjoined=False)
Graphics object consisting of 6 graphics primitives

Plot with thick, red lines:

sage: m.plot_spiderweb(rgbcolor=[1,0,0], thickness=3)
Graphics object consisting of 21 graphics primitives

To generate the unit circle map, it’s helpful to see what the original spiderweb looks like:

sage: f(t) = e^(I*t)
sage: fprime(t) = I*e^(I*t)
sage: m = Riemann_Map([f], [fprime], 0, 1000)
sage: m.plot_spiderweb()
Graphics object consisting of 21 graphics primitives

A multiply connected region with corners. We set min_mag higher to remove “fuzz” outside the domain:

sage: ps = polygon_spline([(-4,-2),(4,-2),(4,2),(-4,2)])
sage: z1 = lambda t: ps.value(t); z1p = lambda t: ps.derivative(t)
sage: z2(t) = -2+exp(-I*t); z2p(t) = -I*exp(-I*t)
sage: z3(t) = 2+exp(-I*t); z3p(t) = -I*exp(-I*t)
sage: m = Riemann_Map([z1,z2,z3],[z1p,z2p,z3p],0,ncorners=4) # long time
sage: p = m.plot_spiderweb(withcolor=True,plot_points=500, thickness = 2.0, min_mag=0.1) # long time
riemann_map(pt)#

Return the Riemann mapping of a point.

That is, given pt on the interior of the mapped region, riemann_map will return the point on the unit disk that pt maps to. Note that this method only works for interior points; accuracy breaks down very close to the boundary. To get boundary correspondence, use get_theta_points().

INPUT:

  • pt – A complex number representing the point to be inverse mapped.

OUTPUT:

A complex number representing the point on the unit circle that the input point maps to.

EXAMPLES:

Can work for different types of complex numbers:

sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)
sage: m.riemann_map(0.25 + sqrt(-0.5))
(0.137514...+0.876696...j)
sage: I = CDF.gen()
sage: m.riemann_map(1.3*I)
(-1.56...e-05+0.989694...j)
sage: m.riemann_map(0.4)
(0.73324...+3.2...e-06j)
sage: m.riemann_map(complex(-3, 0.0001))
(1.405757...e-05+8.06...e-10j)
sage.calculus.riemann.analytic_boundary(t, n, epsilon)#

Provides an exact (for n = infinity) Riemann boundary correspondence for the ellipse with axes 1 + epsilon and 1 - epsilon. The boundary is therefore given by e^(I*t)+epsilon*e^(-I*t). It is primarily useful for testing the accuracy of the numerical Riemann_Map.

INPUT:

  • t – The boundary parameter, from 0 to 2*pi

  • n – integer - the number of terms to include. 10 is fairly accurate, 20 is very accurate.

  • epsilon – float - the skew of the ellipse (0 is circular)

OUTPUT:

A theta value from 0 to 2*pi, corresponding to the point on the circle e^(I*theta)

sage.calculus.riemann.analytic_interior(z, n, epsilon)#

Provides a nearly exact computation of the Riemann Map of an interior point of the ellipse with axes 1 + epsilon and 1 - epsilon. It is primarily useful for testing the accuracy of the numerical Riemann Map.

INPUT:

  • z – complex - the point to be mapped.

  • n – integer - the number of terms to include. 10 is fairly accurate, 20 is very accurate.

sage.calculus.riemann.cauchy_kernel(t, args)#

Intermediate function for the integration in analytic_interior().

INPUT:

  • t – The boundary parameter, meant to be integrated over

  • args – a tuple containing:

    • epsilon – float - the skew of the ellipse (0 is circular)

    • z – complex - the point to be mapped.

    • n – integer - the number of terms to include. 10 is fairly accurate, 20 is very accurate.

    • part – will return the real (‘r’), imaginary (‘i’) or complex (‘c’) value of the kernel

sage.calculus.riemann.complex_to_rgb(z_values)#

Convert from a (Numpy) array of complex numbers to its corresponding matrix of RGB values. For internal use of plot_colored() only.

INPUT:

  • z_values – A Numpy array of complex numbers.

OUTPUT:

An \(N \times M \times 3\) floating point Numpy array X, where X[i,j] is an (r,g,b) tuple.

EXAMPLES:

sage: from sage.calculus.riemann import complex_to_rgb
sage: import numpy
sage: complex_to_rgb(numpy.array([[0, 1, 1000]], dtype=numpy.complex128))
array([[[1.        , 1.        , 1.        ],
        [1.        , 0.05558355, 0.05558355],
        [0.17301243, 0.        , 0.        ]]])

sage: complex_to_rgb(numpy.array([[0, 1j, 1000j]], dtype=numpy.complex128))
array([[[1.        , 1.        , 1.        ],
        [0.52779177, 1.        , 0.05558355],
        [0.08650622, 0.17301243, 0.        ]]])
sage.calculus.riemann.complex_to_spiderweb(z_values, dr, dtheta, spokes, circles, rgbcolor, thickness, withcolor, min_mag)#

Converts a grid of complex numbers into a matrix containing rgb data for the Riemann spiderweb plot.

INPUT:

  • z_values – A grid of complex numbers, as a list of lists.

  • dr – grid of floats, the r derivative of z_values. Used to determine precision.

  • dtheta – grid of floats, the theta derivative of z_values. Used to determine precision.

  • spokes – integer - the number of equally spaced radial lines to plot.

  • circles – integer - the number of equally spaced circles about the center to plot.

  • rgbcolor – float array - the red-green-blue color of the lines of the spiderweb.

  • thickness – positive float - the thickness of the lines or points in the spiderweb.

  • withcolor – boolean - If True the spiderweb will be overlaid on the basic color plot.

  • min_mag – float - The magnitude cutoff below which spiderweb points are not drawn. This only applies to multiply connected domains and is designed to prevent “fuzz” at the edge of the domain. Some complicated multiply connected domains (particularly those with corners) may require a larger value to look clean outside.

OUTPUT:

An \(N x M x 3\) floating point Numpy array X, where X[i,j] is an (r,g,b) tuple.

EXAMPLES:

sage: from sage.calculus.riemann import complex_to_spiderweb
sage: import numpy
sage: zval = numpy.array([[0,1,1000], [.2+.3j,1,-.3j], [0,0,0]],
....:                    dtype=numpy.complex128)
sage: deriv = numpy.array([[.1]],dtype = numpy.float64)
sage: complex_to_spiderweb(zval, deriv, deriv, 4, 4, [0,0,0], 1, False, 0.001)
array([[[1., 1., 1.],
        [1., 1., 1.],
        [1., 1., 1.]],

       [[1., 1., 1.],
        [0., 0., 0.],
        [1., 1., 1.]],

       [[1., 1., 1.],
        [1., 1., 1.],
        [1., 1., 1.]]])

sage: complex_to_spiderweb(zval, deriv, deriv, 4, 4, [0,0,0], 1, True, 0.001)
array([[[1.        , 1.        , 1.        ],
        [1.        , 0.05558355, 0.05558355],
        [0.17301243, 0.        , 0.        ]],

       [[1.        , 0.96804683, 0.48044583],
        [0.        , 0.        , 0.        ],
        [0.77351965, 0.5470393 , 1.        ]],

       [[1.        , 1.        , 1.        ],
        [1.        , 1.        , 1.        ],
        [1.        , 1.        , 1.        ]]])
sage.calculus.riemann.get_derivatives(z_values, xstep, ystep)#

Computes the r*e^(I*theta) form of derivatives from the grid of points. The derivatives are computed using quick-and-dirty taylor expansion and assuming analyticity. As such get_derivatives is primarily intended to be used for comparisons in plot_spiderweb and not for applications that require great precision.

INPUT:

  • z_values – The values for a complex function evaluated on a grid in the complex plane, usually from compute_on_grid.

  • xstep – float, the spacing of the grid points in the real direction

OUTPUT:

  • A tuple of arrays, [dr, dtheta], with each array 2 less in both dimensions than z_values

    • dr - the abs of the derivative of the function in the +r direction

    • dtheta - the rate of accumulation of angle in the +theta direction

EXAMPLES:

Standard usage with compute_on_grid:

sage: from sage.calculus.riemann import get_derivatives
sage: f(t) = e^(I*t) - 0.5*e^(-I*t)
sage: fprime(t) = I*e^(I*t) + 0.5*I*e^(-I*t)
sage: m = Riemann_Map([f], [fprime], 0)
sage: data = m.compute_on_grid([],19)
sage: xstep = (data[2]-data[1])/19
sage: ystep = (data[4]-data[3])/19
sage: dr, dtheta = get_derivatives(data[0],xstep,ystep)
sage: dr[8,8]
0.241...
sage: dtheta[5,5]
5.907...