Riemann Mapping

AUTHORS:

  • Ethan Van Andel (2009): initial version
  • Robert Bradshaw (2009): his “complex_plot” was adapted for plot_colored

Development supported by NSF award No. 0702939.

class sage.calculus.riemann.ColorPlot

Bases: sage.plot.primitive.GraphicPrimitive

The GraphicsPrimitive to display complex functions in using the domain
coloring method

INPUT:

  • z_values – An array of complex values to be plotted.
  • x_range – A minimum and maximum x value for the plot.
  • y_range – A minimum and maximum y value for the plot.
EXAMPLES:
:: sage: p = complex_plot(lambda z: z^2-1, (-2, 2), (-2, 2))
get_minmax_data()

Returns a dictionary with the bounding box data.

EXAMPLES:

sage: p = complex_plot(lambda z: z, (-1, 2), (-3, 4))
sage: sorted(p.get_minmax_data().items())
[('xmax', 2.0), ('xmin', -1.0), ('ymax', 4.0), ('ymin', -3.0)]
class sage.calculus.riemann.Riemann_Map

Bases: object

The Riemann_Map class computes a Riemann or Ahlfors map from supplied data. It also has various methods to provide information about 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 numeric rather than analytic, for unusual regions or insufficient collocation points may give very inaccurate results.

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.

The following inputs must all 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.
  • ncorners – integer (default: 4), if mapping a figure with (equally t-spaced) corners, better results may be obtained by accurately giving this parameter. Used to add the proper constant to the theta correspondance function.
  • opp – boolean (default: False), set to True in very rare cases where the theta correspondance 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: m = Riemann_Map([lambda t: e^(I*t)], [lambda t: I*e^(I*t)], 0)  # long time (4 sec)

The unit circle with a small hole:

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

A square:

sage: ps = polygon_spline([(-1, -1), (1, -1), (1, 1), (-1, 1)])  # long time
sage: f = lambda t: ps.value(t)  # long time
sage: fprime = lambda t: ps.derivative(t)  # long time
sage: m = Riemann_Map([f], [fprime], 0.25, ncorners=4)  # long time
sage: m.plot_colored() + m.plot_spiderweb()  # long time

Compute rough error for this map:

sage: x = 0.75  # long time
sage: print "error =", m.inverse_riemann_map(m.riemann_map(x)) - x  # long time
error = (-0.0001211863...+0.001606139...j)

ALGORITHM:

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

REFERENCES:

[KT]N. Kerzman and M. R. Trummer. “Numerical Conformal Mapping via the Szego kernel”. Journal of Computational and Applied Mathematics, 14(1-2): 111–123, 1986.
get_szego(boundary=-1, absolute_value=False)

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

INPUT:

The following inputs 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.
  • 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: m = Riemann_Map([lambda t: e^(I*t) - 0.5*e^(-I*t)], [lambda t: I*e^(I*t) + 0.5*I*e^(-I*t)], 0)  # long time (4 sec)
sage: sz = m.get_szego(boundary=0)  # long time
sage: points = m.get_szego(absolute_value=True)  # long time
sage: list_plot(points)  # long time

Extending the points by a spline:

sage: s = spline(points)  # long time
sage: s(3*pi / 4)  # long time
0.00121587378429...

The unit circle with a small hole:

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

Getting the szego for a specifc boundary:

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

Returns 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. 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: m = Riemann_Map([lambda t: e^(I*t) - 0.5*e^(-I*t)], [lambda t: I*e^(I*t) + 0.5*I*e^(-I*t)], 0)  # long time (4 sec)
sage: points = m.get_theta_points()  # long time
sage: list_plot(points)  # long time

Extending the points by a spline:

sage: s = spline(points)  # long time
sage: s(3*pi / 4)  # long time
1.62766037996...

The unit circle with a small hole:

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

Getting the szego for a specifc boundary:

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

Returns the inverse Riemann mapping of a point. That is, given pt on the interior of the unit disc, inverse_reimann_map() will return the point on the original region that would be Riemann mapped to pt.

Note

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: m = Riemann_Map([lambda t: e^(I*t) - 0.5*e^(-I*t)], [lambda t: I*e^(I*t) + 0.5*I*e^(-I*t)], 0)  # long time (4 sec)
sage: m.inverse_riemann_map(0.5 + sqrt(-0.5))  # long time
(0.406880548363...+0.361470279816...j)
sage: m.inverse_riemann_map(0.95)  # long time
(0.486319431795...-4.90019052...j)
sage: m.inverse_riemann_map(0.25 - 0.3*I)  # long time
(0.165324498558...-0.180936785500...j)
sage: import numpy as np  # long time
sage: m.inverse_riemann_map(np.complex(-0.2, 0.5))  # long time
(-0.156280570579...+0.321819151891...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 must all 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: m = Riemann_Map([lambda t: e^(I*t) - 0.5*e^(-I*t)], [lambda t: I*e^(I*t) + 0.5*I*e^(-I*t)], 0)  # long time (4 sec)

Default plot:

sage: m.plot_boundaries()  # long time

Big blue collocation points:

sage: m.plot_boundaries(plotjoined=False, rgbcolor=[0,0,1], thickness=6)  # long time
plot_colored(plot_range=[], plot_points=100)

Draws a colored plot of the Riemann map. A red point on the colored plot corresponds to a red point on the unit disc. Note that this method DOES work for multiply connected domains.

INPUT:

The following inputs must all 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 each direction of the grid. Note that very large values can cause this function to run slowly.

EXAMPLES:

Given a Riemann map m, general usage:

sage: m = Riemann_Map([lambda t: e^(I*t) - 0.5*e^(-I*t)], [lambda t: I*e^(I*t) + 0.5*I*e^(-I*t)], 0)  # long time (4 sec)
sage: m.plot_colored() #long time

Plot zoomed in on a specific spot:

sage: m.plot_colored(plot_range=[-1,1,2,3])  # long time

High resolution plot:

sage: m.plot_colored(plot_points=1000)  # long time (40 sec)

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

sage: m = Riemann_Map([lambda t: e^(I*t)], [lambda t: I*e^(I*t)], 0, 1000)  # long time (10 sec)
sage: m.plot_colored()  # long time
plot_spiderweb(spokes=16, circles=4, pts=32, linescale=0.98999999999999999, rgbcolor=[, 0, 0, 0], thickness=1, plotjoined=True)

Generates a traditional “spiderweb plot” of the Riemann map. Shows what concentric circles and radial lines map to. Note that this method DOES NOT work for multiply connected domains.

INPUT:

The following inputs must all 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.
  • linescale – float between 0 and 1. Shrinks the radial lines away from the boundary to reduce erratic behavior.
  • 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.

EXAMPLES:

General usage:

sage: m = Riemann_Map([lambda t: e^(I*t) - 0.5*e^(-I*t)], [lambda t: I*e^(I*t) + 0.5*I*e^(-I*t)], 0)  # long time (4 sec)

Default plot:

sage: m.plot_spiderweb()  # long time

Simplified plot with many discrete points:

sage: m.plot_spiderweb(spokes=4, circles=1, pts=400, linescale=0.95, plotjoined=False)  # long time

Plot with thick, red lines:

sage: m.plot_spiderweb(rgbcolor=[1,0,0], thickness=3)  # long time

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

sage: m = Riemann_Map([lambda t: e^(I*t)], [lambda t: I*e^(I*t)], 0, 1000)  # long time (8 sec)
sage: m.plot_spiderweb()  # long time
riemann_map(pt)

Returns 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; it breaks down very close to the boundary. To get boundary corrospondance, 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: m = Riemann_Map([lambda t: e^(I*t) - 0.5*e^(-I*t)], [lambda t: I*e^(I*t) + 0.5*I*e^(-I*t)], 0)  # long time (4 sec)
sage: m.riemann_map(0.25 + sqrt(-0.5))  # long time
(0.137514137885...+0.876696023004...j)
sage: m.riemann_map(1.3*I)  # long time
(-1.561029396...+0.989694535737...j)
sage: I = CDF.gen()  # long time
sage: m.riemann_map(0.4)  # long time
(0.733242677182...+3.50767714...j)
sage: import numpy as np  # long time
sage: m.riemann_map(np.complex(-3, 0.0001))  # long time
(1.405757135...+1.05106...j)

Previous topic

Solving ODE numericaly by GSL

Next topic

Complex Interpolation

This Page