Calculus plots with Makie

The Makie.jl webpage says

From the Jpanese word Maki-e, which is a technique to sprinkle lacquer with gold and silver powder. Data is basically the gold and silver of our age, so let's spread it out beautifully on the screen!

Makie itself is a metapackage for a rich ecosystem. We show how to use the interface provided by AbstractPlotting and the GLMakie backend to produce the familiar graphics of calculus. We do not discuss the MakieLayout package which provides a means to layout multiple graphics and add widgets, such as sliders and buttons, to a layout. We do not discuss MakieRecipes. For Plots, there are "recipes" that make some of the plots more straightforward. We do not discuss the AlgebraOfGraphics which presents an interface for the familiar graphics of statistics. The MakieGallery shows many exmaples of the use of Makie.

Scenes

Makie draws graphics onto a canvas termed a "scene" in the Makie documentation. There are GLMakie, WGLMakie, and CairoMakie backends for different types of canvases. In the following, we have used GLMakie. WGLMakie is useful for incorporating Makie plots into web-based technologies.

We begin by loading our two packages:

using AbstractPlotting
using GLMakie
#using WGLMakie; WGLMakie.activate!()
#AbstractPlotting.set_theme!(scale_figure=false, resolution = (480, 400))

The Makie developers have workarounds for the delayed time to first plot, but without utilizing these the time to load the package is lengthy.

A scene is produced with Scene() or through a plotting primitive:

scene = Scene()

We see next how to move beyond the blank canvas.

Points (scatter)

The task of plotting the points, say $(1,2)$, $(2,3)$, $(3,2)$ can be done different ways. Most plotting packages, and Makie is no exception, allow the following: form vectors of the $x$ and $y$ values then plot those with scatter:

xs = [1,2,3]
ys = [2,3,2]
scatter(xs, ys)

The scatter function creates and returns a Scene object, which when displayed shows the plot.

The more generic plot function can also be used for this task.

Point2, Point3

When learning about points on the Cartesian plane, a "t"-chart is often produced:

x | y
-----
1 | 2
2 | 3
3 | 2

The scatter usage above used the columns. The rows are associated with the points, and these too can be used to produce the same graphic. Rather than make vectors of $x$ and $y$ (and optionally $z$) coordinates, it is more idiomatic to create a vector of "points." Makie utilizes a Point type to store a 2 or 3 dimensional point. The Point2 and Point3 constructors will be utilized.

Makie uses a GPU, when present, to accelerate the graphic rendering. GPUs employ 32-bit numbers. Julia uses an f0 to indicate 32-bit floating points. Hence the alternate types Point2f0 to store 2D points as 32-bit numbers and Points3f0 to store 3D points as 32-bit numbers are seen in the documentation for Makie.

We can plot vector of points in as direct manner as vectors of their coordinates:

pts = [Point2(1,2), Point2(2,3), Point2(3,2)]
scatter(pts)

A typical usage is to generate points from some vector-valued function. Say we have a parameterized function r taking $R$ into $R^2$ defined by:

r(t) = [sin(t), cos(t)]
r (generic function with 1 method)

Then broadcasting values gives a vector of vectors, each identified with a point:

ts = [1,2,3]
r.(ts)
3-element Array{Array{Float64,1},1}:
 [0.8414709848078965, 0.5403023058681398]
 [0.9092974268256817, -0.4161468365471424]
 [0.1411200080598672, -0.9899924966004454]

We can broadcast Point2 over this to create a vector of Point objects:

pts = Point2.(r.(ts))
3-element Array{GeometryBasics.Point{2,Float64},1}:
 [0.8414709848078965, 0.5403023058681398]
 [0.9092974268256817, -0.4161468365471424]
 [0.1411200080598672, -0.9899924966004454]

These then can be plotted directly:

scatter(pts)

The ploting of points in three dimesions is essentially the same, save the use of Point3 instead of Point2.

r(t) = [sin(t), cos(t), t]
ts = range(0, 4pi, length=100)
pts = Point3.(r.(ts))
scatter(pts)

To plot points generated in terms of vectors of coordinates, the component vectors must be created. The "t"-table shows how, simply loop over each column and add the corresponding $x$ or $y$ (or $z$) value. This utility function does exactly that, returning the vectors in a tuple.

unzip(vs) = Tuple([vs[j][i] for j in eachindex(vs)] for i in eachindex(vs[1]))
unzip (generic function with 1 method)

(The functionality is essentially a reverse of the zip function, hence the name.)

We might have then:

scatter(unzip(r.(ts))...)

where splatting is used to specify the xs, ys, and zs to scatter.

(Compare to scatter(Point3.(r.(ts))) or scatter(Point3∘r).(ts)).)

Attributes

A point is drawn with a "marker" with a certain size and color. These attributes can be adjusted, as in the following:

scatter(xs, ys, marker=[:x,:cross, :circle], markersize=25, color=:blue)

Marker attributes include

Text (text)

Text can be placed at a point, as a marker is. To place text the desired text and a position need to be specified.

For example:

pts = Point2.(1:5, 1:5)
scene = scatter(pts)
[text!(scene, "text", position=pt, textsize=1/i, rotation=2pi/i) for (i,pt) in enumerate(pts)]
scene

The graphic shows that position positions the text, textsize adjusts the displayed size, and rotation adjusts the orientation.

Attributes for text include:

Curves

Plots of univariate functions

The basic plot of univariate calculus is the graph of a function $f$ over an interval $[a,b]$. This is implemented using a familiar strategy: produce a series of representative values between $a$ and $b$; produce the corresponding $f(x)$ values; plot these as points and connect the points with straight lines. The lines function of AbstractPlotting will do the last step.

By taking a sufficient number of points within $[a,b]$ the connect-the-dot figure will appear curved, when the function is.

To create regular values between a and b either the range function, the related LinRange function, or the range operator (a:h:b) are employed.

For example:

f(x) = sin(x)
a, b = 0, 2pi
xs = range(a, b, length=250)
lines(xs, f.(xs))

Or

f(x) = cos(x)
a, b = -pi, pi
xs = a:pi/100:b
lines(xs, f.(xs))

As with scatter, lines returns a Scene object that produces a graphic when displayed.

As with scatter, lines can can also be drawn using a vector of points:

lines([Point2(x, fx) for (x,fx) in zip(xs, f.(xs))])

(Though the advantage isn't clear here, this will be useful when the points are more naturally generated.)

When a y value is NaN or infinite, the connecting lines are not drawn:

xs = 1:5
ys = [1,2,NaN, 4, 5]
lines(xs, ys)

As with other plotting packages, this is useful to represent discontinuous functions, such as what occurs at a vertical asymptote.

Adding to a scene (lines!, scatter!, ...)

To add or modify a scene can be done using a mutating version of a plotting primitive, such as lines! or scatter!. The names follow Julia's convention of using an ! to indicate that a function modifies an argument, in this case the scene.

Here is one way to show two plots at once:

xs = range(0, 2pi, length=100)
scene = lines(xs, sin.(xs))
lines!(scene, xs, cos.(xs))

We will see soon how to modify the line attributes so that the curves can be distinguished.

The following shows the construction details in the graphic, and that the initial scene argument is implicitly assumed:

xs = range(0, 2pi, length=10)
lines(xs, sin.(xs))
scatter!(xs, sin.(xs), markersize=10)

The current scene will have data limits that can be of interest. The following indicates how they can be manipulated to get the limits of the displayed x values.

xs = range(0, 2pi, length=200)
scene = plot(xs, sin.(xs))
rect = scene.data_limits[] # get limits for g from f
a, b = rect.origin[1],  rect.origin[1] + rect.widths[1]
(-0.9633175f0, 6.2831855f0)

In the output it can be discerned that the values are 32-bit floating point numbers and yield a slightly larger interval than specified in xs.

As an example, this shows how to add the tangent line to a graph. The slope of the tangent line being computed by ForwardDiff.derivative.

using ForwardDiff
f(x) = x^x
a, b= 0, 2
c = 0.5
xs = range(a, b, length=200)

tl(x) = f(c) + ForwardDiff.derivative(f, c) * (x-c)

scene = lines(xs, f.(xs))
lines!(scene, xs, tl.(xs), color=:blue)

Attributes

In the last example, we added the argument color=:blue to the lines! call. This set an attribute for the line being drawn. Lines have other attributes that allow different ones to be distinguished, as above where colors indicate the different graphs.

Other attributes can be seen from the help page for lines, and include:

A legend can also be used to help identify different curves on the same graphic, though this is not illustrated. There are examples in the Makie gallery.

Scene attributes

Attributes of the scene include any titles and labels, the limits that define the coordinates being displayed, the presentation of tick marks, etc.

The title function can be used to add a title to a scene. The calling syntax is title(scene, text).

To set the labels of the graph, there are "shorthand" functions xlabel!, ylabel!, and zlabel!. The calling pattern would follow xlabel!(scene, "x-axis").

The plotting ticks and their labels are returned by the unexported functions tickranges and ticklabels. The unexported xtickrange, ytickrange, and ztickrange; and xticklabels, yticklabels, and zticklabels return these for the indicated axes.

These can be dynamically adjusted using xticks!, yticks!, or zticks!.

pts = [Point2(1,2), Point2(2,3), Point2(3,2)]
scene = scatter(pts)
title(scene, "3 points")
ylabel!(scene, "y values")
xticks!(scene, xtickrange=[1,2,3], xticklabels=["a", "b", "c"])

To set the limits of the graph there are shorthand functions xlims!, ylims!, and zlims!. This might prove useful if vertical asymptotes are encountered, as in this example:

f(x) = 1/x
a,b = -1, 1
xs = range(-1, 1, length=200)
scene = lines(xs, f.(xs))
ylims!(scene, (-10, 10))
center!(scene)

Plots of parametric functions

A space curve is a plot of a function $f:R^2 \rightarrow R$ or $f:R^3 \rightarrow R$.

To construct a curve from a set of points, we have a similar pattern in both $2$ and $3$ dimensions:

r(t) = [sin(2t), cos(3t)]
ts = range(0, 2pi, length=200)
pts = Point2.(r.(ts))  # or (Point2∘r).(ts)
lines(pts)

Or

r(t) = [sin(2t), cos(3t), t]
ts = range(0, 2pi, length=200)
pts = Point3.(r.(ts))
lines(pts)

Alternatively, vectors of the $x$, $y$, and $z$ components can be produced and then plotted using the pattern lines(xs, ys) or lines(xs, ys, zs). For example, using unzip, as above, we might have done the prior example with:

xs, ys, zs = unzip(r.(ts))
lines(xs, ys, zs)

Tangent vectors (arrows)

A tangent vector along a curve can be drawn quite easily using the arrows function. There are different interfaces for arrows, but we show the one which uses a vector of positions and a vector of "vectors". For the latter, we utilize the derivative function from ForwardDiff:

using ForwardDiff
r(t) = [sin(t), cos(t)] # vector, not tuple
ts = range(0, 4pi, length=200)
scene = Scene()
lines!(scene, Point2.(r.(ts)))

nts = 0:pi/4:2pi
us = r.(nts)
dus = ForwardDiff.derivative.(r, nts)

arrows!(scene, Point2.(us), Point2.(dus))

In 3 dimensions the differences are minor:

r(t) = [sin(t), cos(t), t] # vector, not tuple
ts = range(0, 4pi, length=200)
scene = Scene()
lines!(scene, Point3.(r.(ts)))

nts = pi:pi/4:3pi
us = r.(nts)
dus = ForwardDiff.derivative.(r, nts)

arrows!(scene, Point3.(us), Point3.(dus))
Attributes

Attributes for arrows include

Implicit equations (2D)

The graph of an equation is the collection of all $(x,y)$ values satisfying the equation. This is more general than the graph of a function, which can be viewed as the graph of the equation $y=f(x)$. An equation in $x$-$y$ can be graphed if the set of solutions to a related equation $f(x,y)=0$ can be identified, as one can move all terms to one side of an equation and define $f$ as the rule of the side with the terms.

The MDBM (Multi-Dimensional Bisection Method) package can be used for the task of characterizing when $f(x,y)=0$. (Also IntervalConstraintProgramming can be used.) We first wrap its interface and then define a "plot" recipe (through method overloading, not through MakieRecipes).

using MDBM
function implicit_equation(f, axes...; iteration::Int=4, constraint=nothing)

    axes = [axes...]

    if constraint == nothing
        prob = MDBM_Problem(f, axes)
    else
        prob = MDBM_Problem(f, axes, constraint=constraint)
    end

    solve!(prob, iteration)

    prob
end
implicit_equation (generic function with 1 method)

The implicit_equation function is just a simplified wrapper for the MDBM_Problem interface. It creates an object to be plotted in a manner akin to:

f(x,y) = x^3 + x^2 + x + 1 - x*y        # solve x^3 + x^2 + x + 1 = x*y
ie = implicit_equation(f, -5:5, -10:10)
MDBM.MDBM_Problem{MDBM.MemF{typeof(Main.##WeaveSandBox#253.f),MDBM.var"#17#
19",Float64,Bool,Tuple{Float64,Float64}},2,1,1,StaticArrays.SArray{Tuple{4}
,StaticArrays.SArray{Tuple{2},Bool,1,2},1,4},StaticArrays.SArray{Tuple{3,4}
,Float64,2,12},Int64,Float64,Tuple{MDBM.Axis{Float64},MDBM.Axis{Float64}}}(
MDBM.MemF{typeof(Main.##WeaveSandBox#253.f),MDBM.var"#17#19",Float64,Bool,T
uple{Float64,Float64}}(Main.##WeaveSandBox#253.f, MDBM.var"#17#19"(), MDBM.
MDBMcontainer{Float64,Bool,Tuple{Float64,Float64}}[MDBM.MDBMcontainer{Float
64,Bool,Tuple{Float64,Float64}}(-154.0, true, (-5.0, -10.0)), MDBM.MDBMcont
ainer{Float64,Bool,Tuple{Float64,Float64}}(-149.0, true, (-5.0, -9.0)), MDB
M.MDBMcontainer{Float64,Bool,Tuple{Float64,Float64}}(-144.0, true, (-5.0, -
8.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{Float64,Float64}}(-139.0, true
, (-5.0, -7.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{Float64,Float64}}(-1
34.0, true, (-5.0, -6.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{Float64,Fl
oat64}}(-129.0, true, (-5.0, -5.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{
Float64,Float64}}(-124.0, true, (-5.0, -4.0)), MDBM.MDBMcontainer{Float64,B
ool,Tuple{Float64,Float64}}(-119.0, true, (-5.0, -3.0)), MDBM.MDBMcontainer
{Float64,Bool,Tuple{Float64,Float64}}(-114.0, true, (-5.0, -2.0)), MDBM.MDB
Mcontainer{Float64,Bool,Tuple{Float64,Float64}}(-109.0, true, (-5.0, -1.0))
  …  MDBM.MDBMcontainer{Float64,Bool,Tuple{Float64,Float64}}(151.0, true, (
5.0, 1.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{Float64,Float64}}(146.0, 
true, (5.0, 2.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{Float64,Float64}}(
141.0, true, (5.0, 3.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{Float64,Flo
at64}}(136.0, true, (5.0, 4.0)), MDBM.MDBMcontainer{Float64,Bool,Tuple{Floa
t64,Float64}}(131.0, true, (5.0, 5.0)), MDBM.MDBMcontainer{Float64,Bool,Tup
le{Float64,Float64}}(126.0, true, (5.0, 6.0)), MDBM.MDBMcontainer{Float64,B
ool,Tuple{Float64,Float64}}(121.0, true, (5.0, 7.0)), MDBM.MDBMcontainer{Fl
oat64,Bool,Tuple{Float64,Float64}}(116.0, true, (5.0, 8.0)), MDBM.MDBMconta
iner{Float64,Bool,Tuple{Float64,Float64}}(111.0, true, (5.0, 9.0)), MDBM.MD
BMcontainer{Float64,Bool,Tuple{Float64,Float64}}(106.0, true, (5.0, 10.0))]
, [18976]), MDBM.Axes{2,Tuple{MDBM.Axis{Float64},MDBM.Axis{Float64}}}(([-5.
0, -4.9375, -4.875, -4.8125, -4.75, -4.6875, -4.625, -4.5625, -4.5, -4.4375
  …  4.4375, 4.5, 4.5625, 4.625, 4.6875, 4.75, 4.8125, 4.875, 4.9375, 5.0],
 [-10.0, -9.9375, -9.875, -9.8125, -9.75, -9.6875, -9.625, -9.5625, -9.5, -
9.4375  …  9.4375, 9.5, 9.5625, 9.625, 9.6875, 9.75, 9.8125, 9.875, 9.9375,
 10.0])), MDBM.NCube{Int64,Float64,2}[MDBM.NCube{Int64,Float64,2}([23, 316]
, [1, 1], [1.6063732892989169, 0.25275584507525284], true), MDBM.NCube{Int6
4,Float64,2}([23, 317], [1, 1], [1.3026770110598742, 0.205532987510733], tr
ue), MDBM.NCube{Int64,Float64,2}([23, 318], [1, 1], [0.9973904105883838, 0.
15779865481363942], true), MDBM.NCube{Int64,Float64,2}([23, 319], [1, 1], [
0.6905013026761199, 0.10954668019693573], true), MDBM.NCube{Int64,Float64,2
}([23, 320], [1, 1], [0.3819973811799878, 0.06077080938709818], true), MDBM
.NCube{Int64,Float64,2}([24, 310], [1, 1], [1.5730148223057674, 0.252441931
36240256], true), MDBM.NCube{Int64,Float64,2}([24, 311], [1, 1], [1.2635163
443540753, 0.2033503249592459], true), MDBM.NCube{Int64,Float64,2}([24, 312
], [1, 1], [0.9523386316264015, 0.15370721710235083], true), MDBM.NCube{Int
64,Float64,2}([24, 313], [1, 1], [0.6394683680980245, 0.1035057125802353], 
true), MDBM.NCube{Int64,Float64,2}([24, 314], [1, 1], [0.324892101175531, 0
.05273881477153529], true)  …  MDBM.NCube{Int64,Float64,2}([119, 316], [1, 
1], [0.5737847603815756, -0.10253884896868308], true), MDBM.NCube{Int64,Flo
at64,2}([119, 317], [1, 1], [0.9241595306873377, -0.16592313346062373], tru
e), MDBM.NCube{Int64,Float64,2}([119, 318], [1, 1], [1.2776102439518204, -0
.23045619627359337], true), MDBM.NCube{Int64,Float64,2}([119, 319], [1, 1],
 [1.6341761101488126, -0.2961613662405615], true), MDBM.NCube{Int64,Float64
,2}([120, 315], [1, 1], [-1.6639239167933426, 0.28200587702155877], true), 
MDBM.NCube{Int64,Float64,2}([120, 316], [1, 1], [-1.3398508634660045, 0.228
05972144102202], true), MDBM.NCube{Int64,Float64,2}([120, 317], [1, 1], [-1
.013129945313339, 0.17319397496295114], true), MDBM.NCube{Int64,Float64,2}(
[120, 318], [1, 1], [-0.6837296229134533, 0.11739123110059835], true), MDBM
.NCube{Int64,Float64,2}([120, 319], [1, 1], [-0.35161787289425867, 0.060633
695953389216], true), MDBM.NCube{Int64,Float64,2}([120, 320], [1, 1], [-0.0
1676217906067428, 0.002903178176581569], true)], StaticArrays.SArray{Tuple{
2},Bool,1,2}[[0, 0], [1, 0], [0, 1], [1, 1]], [-0.25 0.25 -0.25 0.25; -0.25
 -0.25 0.25 0.25; -0.25 -0.25 -0.25 -0.25])

The function definition is straightforward. The limits for x and y are specified in the above using ranges. This specifies the initial grid of points for the apdaptive algorithm used by MDBM to identify solutions.

To visualize the output, we make a new method for plot and plot!. There is a distinction between 2 and 3 dimensions. Below in two dimensions curve(s) are drawn. In three dimensions, scaled cubes are used to indicate the surface.

AbstractPlotting.plot(m::MDBM_Problem; kwargs...) = plot!(Scene(), m; kwargs...)
AbstractPlotting.plot!(m::MDBM_Problem; kwargs...) = plot!(AbstractPlotting.current_scene(), m; kwargs...)
AbstractPlotting.plot!(scene::AbstractPlotting.Scene, m::MDBM_Problem; kwargs...) =
    plot!(Val(_dim(m)), scene, m; kwargs...)

_dim(m::MDBM.MDBM_Problem{a,N,b,c,d,e,f,g,h}) where {a,N,b,c,d,e,f,g,h} = N
_dim (generic function with 1 method)

Dispatch is used for the two different dimesions, identified through _dim, defined above.

# 2D plot
function AbstractPlotting.plot!(::Val{2}, scene::AbstractPlotting.Scene,
    m::MDBM_Problem; color=:black, kwargs...)

    mdt=MDBM.connect(m)
    for i in 1:length(mdt)
        dt=mdt[i]
        P1=getinterpolatedsolution(m.ncubes[dt[1]], m)
        P2=getinterpolatedsolution(m.ncubes[dt[2]], m)
        lines!(scene, [P1[1],P2[1]],[P1[2],P2[2]], color=color, kwargs...)
    end

    scene
end
# 3D plot
function AbstractPlotting.plot!(::Val{3}, scene::AbstractPlotting.Scene,
    m::MDBM_Problem; color=:black, kwargs...)

    positions = Point{3, Float32}[]
    scales = Vec3[]
    
    mdt=MDBM.connect(m)
    for i in 1:length(mdt)
        dt=mdt[i]
        P1=getinterpolatedsolution(m.ncubes[dt[1]], m)
        P2=getinterpolatedsolution(m.ncubes[dt[2]], m)

        a, b = Vec3(P1), Vec3(P2)
        push!(positions, Point3(P1))
        push!(scales, b-a)
    end

    cube = Rect{3, Float32}(Vec3(-0.5, -0.5, -0.5), Vec3(1, 1, 1))
    meshscatter!(scene, positions, marker=cube, scale = scales, color=color, transparency=true, kwargs...)
    
    scene
end

We see that the equation ie has two pieces. (This is known as Newton's trident, as Newton was interested in this form of equation.)

plot(ie)

Surfaces

Plots of surfaces in 3 dimensions are useful to help understand the behavior of multivariate functions.

Surfaces defined through $z=f(x,y)$

The "peaks" function generates the logo for MATLAB. Here we see how it can be plotted over the region $[-5,5]\times[-5,5]$.

peaks(x,y) = 3*(1-x)^2*exp(-x^2 - (y+1)^2) - 10(x/5-x^3-y^5)*exp(-x^2-y^2)- 1/3*exp(-(x+1)^2-y^2)
xs = ys = range(-5, 5, length=25)
surface(xs, ys, peaks)

The calling pattern surface(xs, ys, f) implies a rectangular grid over the $x$-$y$ plane defined by xs and ys with $z$ values given by $f(x,y)$.

Alternatively a "matrix" of $z$ values can be specified. For a function f, this is conveniently generated by the pattern f.(xs, ys'), the ' being important to get a matrix of all $x$-$y$ pairs through Julia's broadcasting syntax.

zs = peaks.(xs, ys')
surface(xs, ys, zs)

To see how this graph is constructed, the points $(x,y,f(x,y))$ are plotted over the grid and displayed.

Here we downsample to illutrate

xs = ys = range(-5, 5, length=5)
pts = [Point3(x, y, peaks(x,y)) for x in xs for y in ys]
scatter(pts, markersize=25)

These points are connected. The wireframe function illustrates just the frame

wireframe(xs, ys, peaks.(xs, ys'), linewidth=5)

The surface call triangulates the frame and fills in the shading:

surface!(xs, ys, peaks.(xs, ys'))

Implicitly defined surfaces, $F(x,y,z)=0$

The set of points $(x,y,z)$ satisfying $F(x,y,z) = 0$ will form a surface that can be visualized using the MDBM package. We illustrate showing two nested surfaces.

r₂(x,y,z) = x^2 + y^2 + z^2 - 5/4 # a sphere
r₄(x,y,z) = x^4 + y^4 + z^4 - 1
xs = ys = zs = -2:2
m2,m4 = implicit_equation(r₂, xs, ys, zs), implicit_equation(r₄, xs, ys, zs)

plot(m4, color=:yellow)
plot!(m2, color=:red)

Parametrically defined surfaces

A surface may be parametrically defined through a function $r(u,v) = (x(u,v), y(u,v), z(u,v))$. For example, the surface generated by $z=f(x,y)$ is of the form with $r(u,v) = (u,v,f(u,v))$.

The surface function and the wireframe function can be used to display such surfaces. In previous usages, the x and y values were vectors from which a 2-dimensional grid is formed. For parametric surfaces, a grid for the x and y values must be generated. This function will do so:

function parametric_grid(us, vs, r)
    n,m = length(us), length(vs)
    xs, ys, zs = zeros(n,m), zeros(n,m), zeros(n,m)
    for (i, uᵢ) in enumerate(us)
        for (j, vⱼ) in enumerate(vs)
            x,y,z = r(uᵢ, vⱼ)
            xs[i,j] = x
            ys[i,j] = y
            zs[i,j] = z
        end
    end
    (xs, ys, zs)
end
parametric_grid (generic function with 1 method)

With the data suitably massaged, we can directly plot either a surface or wireframe plot.

For example, a sphere can be parameterized by $r(u,v) = (\sin(u)\cos(v), \sin(u)\sin(v), \cos(u))$ and visualized through:

r(u,v) = [sin(u)*cos(v), sin(u)*sin(v), cos(u)]
us = range(0, pi, length=25)
vs = range(0, pi/2, length=25)
xs, ys, zs = parametric_grid(us, vs, r)

scene = Scene()
surface!(scene, xs, ys, zs)
wireframe!(scene, xs, ys, zs)

A surface of revolution for $g(u)$ revolved about the $z$ axis can be visualized through:

g(u) = u^2 * exp(-u)
r(u,v) = (g(u)*sin(v), g(u)*cos(v), u)
us = range(0, 3, length=10)
vs = range(0, 2pi, length=10)
xs, ys, zs = parametric_grid(us, vs, r)

scene = Scene()
surface!(scene, xs, ys, zs)
wireframe!(scene, xs, ys, zs)

A torus with big radius $2$ and inner radius $1/2$ can be visualized as follows

r1, r2 = 2, 1/2
r(u,v) = ((r1 + r2*cos(v))*cos(u), (r1 + r2*cos(v))*sin(u), r2*sin(v))
us = vs = range(0, 2pi, length=25)
xs, ys, zs = parametric_grid(us, vs, r)

scene = Scene()
surface!(scene, xs, ys, zs)
wireframe!(scene, xs, ys, zs)

A Möbius strip can be produced with:

ws = range(-1/4, 1/4, length=8)
thetas = range(0, 2pi, length=30)
r(w, θ) = ((1+w*cos(θ/2))*cos(θ), (1+w*cos(θ/2))*sin(θ), w*sin(θ/2))
xs, ys, zs = parametric_grid(ws, thetas, r)

scene = Scene()
surface!(scene, xs, ys, zs)
wireframe!(scene, xs, ys, zs)

Contour plots (contour, heatmap)

For a function $z = f(x,y)$ an alternative to a surface plot, is a contour plot. That is, for different values of $c$ the level curves $f(x,y)=c$ are drawn.

For a function $f(x,y)$, the syntax for generating a contour plot follows that for surface.

For example, using the peaks function, previously defined, we have a contour plot over the region $[-5,5]\times[-5,5]$ is generated through:

xs = ys = range(-5, 5, length=100)
contour(xs, ys, peaks)

The default of $5$ levels can be adjusted using the levels keyword:

contour(xs, ys, peaks.(xs, ys'), levels = 20)

(As a reminder, the above also shows how to generate values "zs" to pass to contour instead of a function.)

The contour graph makes identification of peaks and valleys easy as the limits of patterns of nested contour lines.

An alternative visualzation using color to replace contour lines is a heatmap. The heatmap function produces these. The calling syntax is similar to contour and surface:

heatmap(xs, ys, peaks)

This graph shows peaks and valleys through "hotspots" on the graph.

The MakieGallery package includes an example of a surface plot with both a wireframe and 2D contour graph added. It is replicated here using the peaks function scaled by $5$.

The function and domain to plot are described by:

xs = ys = range(-5, 5, length=51)
zs = peaks.(xs, ys') / 5;
51×51 Array{Float64,2}:
 3.24008e-17  1.52592e-16  6.61687e-16  …  5.8987e-18   1.05344e-18
 2.14596e-16  1.01055e-15  4.38226e-15     4.38874e-17  7.76195e-18
 1.30819e-15  6.15854e-15  2.67006e-14     2.96905e-16  5.21764e-17
 7.33909e-15  3.45339e-14  1.49654e-13     1.83635e-15  3.21308e-16
 3.78835e-14  1.78148e-13  7.71488e-13     1.04182e-14  1.81732e-15
 1.79883e-13  8.45241e-13  3.65718e-12  …  5.43295e-14  9.45621e-15
 7.8547e-13   3.68731e-12  1.59369e-11     2.60777e-13  4.53148e-14
 3.15296e-12  1.47845e-11  6.38161e-11     1.15312e-12  2.00128e-13
 1.16297e-11  5.44604e-11  2.34705e-10     4.70007e-12  8.1493e-13
 3.93975e-11  1.84202e-10  7.9236e-10      1.76656e-11  3.06067e-12
 ⋮                                      ⋱               ⋮
 2.89253e-12  1.29362e-11  5.26501e-11     4.8661e-12   8.38318e-13
 8.72615e-13  3.94936e-12  1.63263e-11     1.20336e-12  2.07204e-13
 2.37814e-13  1.08573e-12  4.53856e-12     2.74793e-13  4.72891e-14
 5.8816e-14   2.70298e-13  1.1393e-12   …  5.79456e-14  9.96556e-15
 1.32429e-14  6.11735e-14  2.59494e-13     1.12835e-14  1.93921e-15
 2.72096e-15  1.26206e-14  5.38057e-14     2.02902e-15  3.48448e-16
 5.11068e-16  2.37838e-15  1.0181e-14      3.36938e-16  5.78154e-17
 8.78729e-17  4.10061e-16  1.76118e-15     5.16704e-17  8.85825e-18
 1.3846e-17   6.47614e-17  2.78916e-16  …  7.31756e-18  1.2533e-18

The zs were generated, as wireframe does not provide the interface for passing a function.

The surface and wireframe are produced as follows:

scene = surface(xs, ys, zs)
wireframe!(scene, xs, ys, zs, overdraw = true, transparency = true, color = (:black, 0.1))

To add the contour, a simple call via contour!(scene, xs, ys, zs) will place the contour at the $z=0$ level which will make it hard to read. Rather, placing at the "bottom" of the scene is desirable. To identify that the "scene limits" are queried and the argument transformation = (:xy, zmin) is passed to contour!:

xmin, ymin, zmin = minimum(scene_limits(scene))
contour!(scene, xs, ys, zs, levels = 15, linewidth = 2, transformation = (:xy, zmin))
center!(scene)

The transformation plot attribute sets the "plane" (one of :xy, :yz, or :xz) at a location, in this example zmin.

Three dimensional contour plots

The contour function can also plot $3$-dimensional contour plots. Concentric spheres, contours of $x^2 + y^2 + z^2 = c$ for $c > 0$ are presented by the following:

f(x,y,z) = x^2 + y^2 + z^2
xs = ys = zs = range(-3, 3, length=100)
contour(xs, ys, zs, f)

Vector fields. Visualizations of $f:R^2 \rightarrow R^2$

The vector field $f(x,y) = (y, -x)$ can be visualized as a set of vectors, $f(x,y)$, positioned at a grid. These can be produced with the arrows function. Below we pass a vector of points for the anchors and a vector of points representing the vectors.

We can generate these on a regular grid through:

f(x, y) = [y, -x]
xs = ys = -5:5
pts = vec(Point2.(xs, ys'))
dus = vec(Point2.(f.(xs, ys')))
121-element Array{GeometryBasics.Point{2,Int64},1}:
 [-5, 5]
 [-5, 4]
 [-5, 3]
 [-5, 2]
 [-5, 1]
 [-5, 0]
 [-5, -1]
 [-5, -2]
 [-5, -3]
 [-5, -4]
 ⋮
 [5, 3]
 [5, 2]
 [5, 1]
 [5, 0]
 [5, -1]
 [5, -2]
 [5, -3]
 [5, -4]
 [5, -5]

Broadcasting over (xs, ys') ensures each pair of possible values is encountered. The vec call reshapes an array into a vector.

Calling arrows on the prepared data produces the graphic:

arrows(pts, dus)

The grid seems rotated at first glance. This is due to the length of the vectors as the $(x,y)$ values get farther from the origin. Plotting the normalized values (each will have length $1$) can be done easily using norm (which requires LinearAlgebra to be loaded):

using LinearAlgebra
dvs = dus ./ norm.(dus)
arrows(pts, dvs)

The rotational pattern becomes quite clear now.

The streamplot function also illustrates this phenomenon. This implements an "algorithm [that] puts an arrow somewhere and extends the streamline in both directions from there. Then, it chooses a new position (from the remaining ones), repeating the the exercise until the streamline gets blocked, from which on a new starting point, the process repeats."

The streamplot function expects a point not a pair of values, so we adjust f slightly and call the function using the pattern streamplot(f, xs, ys):

g(x,y) = Point2(f(x,y))
streamplot(g, xs, ys)
ERROR: MethodError: no method matching ^(::GeometryBasics.Point{2,Float64}, ::Int64)
Closest candidates are:
  ^(!Matched::IntervalArithmetic.Interval{Float64}, ::Integer) at /Users/verzani/.julia/packages/IntervalArithmetic/sRFlx/src/intervals/functions.jl:10
  ^(!Matched::Irrational{:ℯ}, ::Integer) at mathconstants.jl:91
  ^(!Matched::Irrational{:ℯ}, ::Number) at mathconstants.jl:91
  ...

(The viewing range could also be adjusted with the -5..5 notation from the IntervalSets package which is brought in when AbstractPlotting is loaded.)

Interacting with a scene

Interaction with a scene is very much integrated into Makie, as the design has a "sophisticated referencing system" which allows sharing of attributes. Adjusting one attribute, can then propogate to others.

In Makie, a Node is a structure that allows its value to be updated, similar to an array. Nodes are Observables, which when changed can trigger an event. Nodes can rely on other nodes, so events can be cascaded.

A simple example is a means to dynamically adjust a label for a scene.

xs, = 1:5, rand(5)
scene = scatter(xs, ys)

We can create a "Node" through:

x = Node("x values")

The value of the node is retrieved by x[], though the function call to_value(x) is recommened, as it will be defined even when x is not a node. This stored value could be used to set the $x$-label in our scene:

xlabel!(scene, x[])

We now set up an observer to update the label whenever the value of x is updated:

on(x) do val
   xlabel!(scen, val)
end

Now setting the value of x will also update the label:

x[] = "The x axis"

A node can be more complicated. This shows how a node of $x$ value can be used to define dependent $y$ values. A scatter plot will update when the $x$ values are updated:

xs = Node(1:10)
ys = lift(a -> f.(a), xs)

The lift function lifts the values of xs to the values of ys.

These can be plotted directly:

scene =  lines(xs, ys)

Changes to the xs values will be reflected in the scene.

xs[] = 2:9

We can incoporporate the two:

lab = lift(val -> "Showing from $(val.start) to $(val.stop)", xs)
on(lab) do val
  xlabel!(scene, val)
  udpate!(scene)
end

The update! call redraws the scene to adjust to increased or decreased range of $x$ values.

The mouse position can be recorded. An example in the gallery of examples shows how.

Here is a hint:

scene = lines(1:5, rand(5))
pos = lift(scene.events.mouseposition) do mpos
    @show AbstractPlotting.to_world(scene, Point2f0(mpos))
end

This will display the coordinates of the mouse in the terminal, as the mouse is moved around.