In [1]:
using Plots
using Blink,Interact

Unable to load WebIO. Please make sure WebIO works for your Jupyter client. For troubleshooting, please see the WebIO/IJulia documentation.

In [12]:
"""
Generate the coefficients of the Legendre polynomials in
three-term recurrence form.

From Gautschi (Numerical Analysis) Table 3.1 (p.178)
"""
function legendre_coeffs(n::Int)
    a = zeros(n)
    b = zeros(n)
    b[1] = 2.
    b[2:end] = 1 ./ (4 .- (1.:(n-1)).^(-2))
    return a, b
end

"""
Evaluate a set of orthogonal polynomials given
the coefficients of the 3-term recurrence
and a set of points. The result
is a matrix where column j is the
polynomial pi_{j-1} evaluated at points ts.
"""
function ortho_eval(a, b, ts)
    n = length(a)
    N = length(ts)
    
    P = zeros(N,n)
    P[:,1] .= 1.
    for j=2:n
        if j==2 # annoying start condition
            P[:,j] = (ts .- a[j-1]).*P[:,j-1]
        else
            P[:,j] = (ts .- a[j-1]).*P[:,j-1] - b[j-1].*P[:,j-2]
        end
    end
    return P
end

ts = collect(range(-1,stop=1,length=1000))
ui = @manipulate for n = 2:25
    a,b = legendre_coeffs(n)
    plot(ts,ortho_eval(a,b,ts))
    xlabel!("t")
end
w = Window()
body!(w,ui)
Out[12]:
Page(5, WebSocket(server, CONNECTED), Dict{String,Any}("webio" => Blink.AtomShell.var"#24#25"{Blink.AtomShell.WebIOBlinkComm}(Blink.AtomShell.WebIOBlinkComm(Window(5, Electron(Process(`/Users/dgleich/.julia/packages/Blink/u1xcH/deps/Julia.app/Contents/MacOS/Julia /Users/dgleich/.julia/packages/Blink/u1xcH/src/AtomShell/main.js port 9218`, ProcessRunning), Sockets.TCPSocket(RawFD(0x0000003a) active, 0 bytes waiting), Dict{String,Any}("callback" => Blink.var"#1#2"())), Page(#= circular reference @-5 =#), Task (done) @0x000000013185af50))),"callback" => Blink.var"#1#2"()), Distributed.Future(1, 1, 5, Some(true)))
In [3]:
sqrt(1/3)
Out[3]:
0.5773502691896257
In [13]:
# Use the eigenvalues to get the zeros!
using LinearAlgebra
"""
Build the Jacobi matrix for a given set of orthogonal polynomials
"""
function jacobi_matrix(a,b)
    n = length(a)
   # J = diagm(a) +diagm(b[2:end],-1) + diagm(ones(n-1),1)
    J = diagm(0=>a) +diagm(-1=>b[2:end]) + diagm(1=>ones(n-1))
end
jacobi_matrix(ab) = jacobi_matrix(ab[1],ab[2]) # little Julia helper

function jacobi_matrix_sym(a,b)
    n = length(a)
    #J = diagm(a) +diagm(sqrt(b[2:end]),-1) + diagm(sqrt(b[2:end]),1)
    J = diagm(0=>a) +diagm(-1=>sqrt.(b[2:end])) + diagm(1=>sqrt.(b[2:end]))
end
jacobi_matrix_sym(ab) = jacobi_matrix_sym(ab[1],ab[2]) # little Julia helper

display(jacobi_matrix(legendre_coeffs(5)))
display(jacobi_matrix_sym(legendre_coeffs(5)))
5×5 Array{Float64,2}:
 0.0       1.0       0.0       0.0       0.0
 0.333333  0.0       1.0       0.0       0.0
 0.0       0.266667  0.0       1.0       0.0
 0.0       0.0       0.257143  0.0       1.0
 0.0       0.0       0.0       0.253968  0.0
5×5 Array{Float64,2}:
 0.0      0.57735   0.0       0.0       0.0
 0.57735  0.0       0.516398  0.0       0.0
 0.0      0.516398  0.0       0.507093  0.0
 0.0      0.0       0.507093  0.0       0.503953
 0.0      0.0       0.0       0.503953  0.0
In [14]:
using Blink,Interact
ts = collect(range(-1,stop=1,length=1000))
ui = @manipulate for n = 2:10
    a,b = legendre_coeffs(n)
    plot(ts,ortho_eval(a,b,ts),framestyle=:origin)
    lams = eigvals(jacobi_matrix_sym(a[1:n-1],b[1:n-1]))
    scatter!(lams, zeros(n-1))
    xlabel!("t")
end
w = Window()
body!(w,ui)
Out[14]:
Page(6, WebSocket(server, CONNECTED), Dict{String,Any}("webio" => Blink.AtomShell.var"#24#25"{Blink.AtomShell.WebIOBlinkComm}(Blink.AtomShell.WebIOBlinkComm(Window(6, Electron(Process(`/Users/dgleich/.julia/packages/Blink/u1xcH/deps/Julia.app/Contents/MacOS/Julia /Users/dgleich/.julia/packages/Blink/u1xcH/src/AtomShell/main.js port 9218`, ProcessRunning), Sockets.TCPSocket(RawFD(0x0000003a) active, 0 bytes waiting), Dict{String,Any}("callback" => Blink.var"#1#2"())), Page(#= circular reference @-5 =#), Task (done) @0x0000000119f2d210))),"callback" => Blink.var"#1#2"()), Distributed.Future(1, 1, 6, Some(true)))
In [18]:
""" Generate the n-point Gauss-Legendre quadrature rule. """

using LinearAlgebra

function gauss_legendre_quadrature(n)
    a,b = legendre_coeffs(n)
    Jsym = jacobi_matrix_sym(a,b)
    x,V = eigen(Symmetric(Jsym)) # force the symmetric solver
    w = V[1,:].^2*b[1]
    p = sortperm(x) # get them in sorted order
    return x[p],w[p]
end

gauss_legendre_quadrature(100)
Out[18]:
([-0.9997137267734407, -0.9984919506395958, -0.9962951347331248, -0.9931249370374435, -0.9889843952429915, -0.9838775407060569, -0.9778093584869183, -0.970785775763706, -0.9628136542558154, -0.9539007829254913  …  0.9539007829254917, 0.9628136542558156, 0.9707857757637063, 0.9778093584869183, 0.983877540706057, 0.9889843952429918, 0.9931249370374434, 0.9962951347331251, 0.9984919506395958, 0.9997137267734413], [0.0007346344905055609, 0.001709392653517898, 0.0026839253715532914, 0.003655961201325863, 0.004624450063423359, 0.0055884280038653485, 0.006546948450845261, 0.007499073255463434, 0.008443871469668752, 0.009380419653694568  …  0.009380419653694504, 0.008443871469669016, 0.007499073255464638, 0.006546948450845328, 0.00558842800386554, 0.0046244500634221326, 0.003655961201326404, 0.0026839253715535264, 0.0017093926535180707, 0.0007346344905056346])
Unexpected end of input
 ...when parsing byte with value '0'
Stacktrace:
 [1] error(::String) at ./error.jl:33
 [2] _error(::String, ::JSON.Parser.StreamingParserState{Sockets.TCPSocket}) at /Users/dgleich/.julia/packages/JSON/d89fA/src/Parser.jl:150
 [3] byteat(::JSON.Parser.StreamingParserState{Sockets.TCPSocket}) at /Users/dgleich/.julia/packages/JSON/d89fA/src/Parser.jl:59
 [4] current at /Users/dgleich/.julia/packages/JSON/d89fA/src/Parser.jl:72 [inlined]
 [5] chomp_space! at /Users/dgleich/.julia/packages/JSON/d89fA/src/Parser.jl:117 [inlined]
 [6] parse_value(::JSON.Parser.ParserContext{Dict{String,Any},Int64,true,nothing}, ::JSON.Parser.StreamingParserState{Sockets.TCPSocket}) at /Users/dgleich/.julia/packages/JSON/d89fA/src/Parser.jl:160
 [7] parse(::Sockets.TCPSocket; dicttype::Type{T} where T, inttype::Type{Int64}, allownan::Bool, null::Nothing) at /Users/dgleich/.julia/packages/JSON/d89fA/src/Parser.jl:494
 [8] parse at /Users/dgleich/.julia/packages/JSON/d89fA/src/Parser.jl:492 [inlined]
 [9] macro expansion at /Users/dgleich/.julia/packages/Lazy/9Xnd3/src/macros.jl:268 [inlined]
 [10] macro expansion at /Users/dgleich/.julia/packages/Blink/u1xcH/src/AtomShell/process.jl:111 [inlined]
 [11] (::Blink.AtomShell.var"#9#10"{Electron})() at ./task.jl:356
In [ ]: