Spherical harmonic addition theorem
This example confirms numerically that
is actually a degree-$(n-1)$ polynomial on $\mathbb{S}^2$, where $P_n$ is the degree-$n$ Legendre polynomial, and $x,y,z \in \mathbb{S}^2$. To verify, we sample the function on a $N\times M$ equiangular grid defined by:
we convert the function samples to Fourier coefficients using plan_sph_analysis
; and finally, we transform the Fourier coefficients to spherical harmonic coefficients using plan_sph2fourier
.
In the basis of spherical harmonics, it is plain to see the addition theorem in action, since $P_n(x\cdot y)$ should only consist of exact-degree-$n$ harmonics.
For the storage pattern of the arrays, please consult the documentation.
function threshold!(A::AbstractArray, ϵ)
for i in eachindex(A)
if abs(A[i]) < ϵ A[i] = 0 end
end
A
end
using FastTransforms, LinearAlgebra, Plots
const GENFIGS = joinpath(pkgdir(FastTransforms), "docs/src/generated")
!isdir(GENFIGS) && mkdir(GENFIGS)
plotlyjs()
Plots.PlotlyJSBackend()
The colatitudinal grid (mod $\pi$):
N = 15
θ = (0.5:N-0.5)/N
0.03333333333333333:0.06666666666666667:0.9666666666666667
The longitudinal grid (mod $\pi$):
M = 2*N-1
φ = (0:M-1)*2/M
0.0:0.06896551724137931:1.9310344827586206
Arbitrarily, we place $x$ at the North pole:
x = [0,0,1]
3-element Vector{Int64}:
0
0
1
Another vector is completely free:
y = normalize([.123,.456,.789])
3-element Vector{Float64}:
0.13375998748853216
0.4958906853233388
0.8580213831581455
Thus $z \in \mathbb{S}^2$ is our variable vector, parameterized in spherical coordinates:
z = (θ,φ) -> [sinpi(θ)*cospi(φ), sinpi(θ)*sinpi(φ), cospi(θ)]
#1 (generic function with 1 method)
On the tensor product grid, the Legendre polynomial $P_n(z\cdot y)$ is:
A = [(2k+1)/(k+1) for k in 0:N-1]
B = zeros(N)
C = [k/(k+1) for k in 0:N]
c = zeros(N); c[N] = 1
pts = vec([z(θ, φ)⋅y for θ in θ, φ in φ])
phi0 = ones(N*M)
F = reshape(FastTransforms.clenshaw!(c, A, B, C, pts, phi0, zeros(N*M)), N, M)
15×29 Matrix{Float64}:
0.264259 0.299171 0.304985 0.288068 0.261038 0.237692 0.228241 0.236528 0.259192 0.286391 0.304437 0.300462 0.26755 0.207601 0.130208 0.0483744 -0.0262845 -0.0865116 -0.129804 -0.156955 -0.169985 -0.170403 -0.158265 -0.132139 -0.0899937 -0.0308858 0.0429655 0.124641 0.202756
0.212068 0.303332 0.129215 -0.186607 -0.387155 -0.383993 -0.338172 -0.379747 -0.392602 -0.206438 0.109443 0.299961 0.225231 -0.0101802 -0.207053 -0.258628 -0.187516 -0.070652 0.0323131 0.0988479 0.12977 0.130741 0.102005 0.0380582 -0.0627707 -0.180155 -0.256899 -0.215501 -0.0267606
-0.115553 0.236712 0.221776 -0.251073 -0.32857 0.421619 0.986223 0.482429 -0.299585 -0.280602 0.198511 0.253762 -0.0930651 -0.258938 -0.0793421 0.160653 0.233501 0.154045 0.0336174 -0.0560314 -0.0990238 -0.100375 -0.0604041 0.0261339 0.145985 0.231751 0.1719 -0.0610864 -0.255787
-0.18454 -0.186901 0.216999 0.207436 -0.278408 -0.36511 -0.189358 -0.35156 -0.304334 0.180746 0.238643 -0.165767 -0.202157 0.128252 0.217098 -0.00943667 -0.200512 -0.194373 -0.0821865 0.0208264 0.0732008 0.0748654 0.0260895 -0.074028 -0.188494 -0.206184 -0.027257 0.20837 0.147067
0.233218 -0.0443128 -0.25662 0.0443223 0.302577 0.209245 0.107407 0.198264 0.305263 0.0716125 -0.25045 -0.0709297 0.228915 0.0560634 -0.209167 -0.109098 0.135993 0.209367 0.121014 0.0107106 -0.049922 -0.0518826 0.00472316 0.112882 0.206906 0.148044 -0.0921222 -0.214775 0.0327421
-0.14959 0.176398 0.159536 -0.157586 -0.253853 -0.135785 -0.0630475 -0.127292 -0.249419 -0.174706 0.141057 0.191837 -0.130845 -0.171771 0.128955 0.181407 -0.058914 -0.206946 -0.154207 -0.0410986 0.0276329 0.029906 -0.034464 -0.146673 -0.208687 -0.0754426 0.17029 0.145612 -0.156134
0.03232 -0.216991 -0.0542418 0.201432 0.206916 0.0892457 0.0316764 0.0822717 0.200041 0.210302 -0.0334057 -0.219505 0.00906962 0.211107 -0.0222883 -0.209099 -0.0231776 0.189332 0.183999 0.0724419 -0.00503206 -0.00766855 0.065182 0.177672 0.19593 -0.0044746 -0.206821 -0.0443208 0.20777
0.0772106 0.200813 -0.0345837 -0.214416 -0.166372 -0.0543381 -0.0063074 -0.0483886 -0.158601 -0.217382 -0.0533452 0.193042 0.0976693 -0.187906 -0.085199 0.193101 0.105896 -0.154519 -0.211246 -0.107108 -0.019323 -0.0162213 -0.09921 -0.206977 -0.166608 0.0879047 0.201099 -0.0639311 -0.196895
-0.161616 -0.14919 0.105731 0.2116 0.130202 0.0251362 -0.0160611 0.0199531 0.122174 0.210305 0.12046 -0.134582 -0.17451 0.114917 0.174042 -0.130385 -0.182374 0.0957622 0.234753 0.148573 0.0475436 0.0437788 0.140028 0.233975 0.113813 -0.1691 -0.148371 0.159712 0.133414
0.208908 0.0728989 -0.161008 -0.198486 -0.0961065 0.00135321 0.037119 0.00590617 -0.0881639 -0.194067 -0.170766 0.0549875 0.211431 -0.00338118 -0.21986 0.016122 0.234375 0.00021735 -0.247537 -0.203094 -0.0836434 -0.0787903 -0.19406 -0.253043 -0.023216 0.231655 0.0407554 -0.218225 -0.0261102
-0.203495 0.0220031 0.200371 0.176174 0.0619327 -0.0270702 -0.058089 -0.0310551 0.0543296 0.169516 0.204695 0.0395021 -0.194934 -0.128043 0.183785 0.141369 -0.211961 -0.150004 0.2205 0.283345 0.137864 0.130789 0.275142 0.23809 -0.126501 -0.226893 0.11958 0.198142 -0.109095
0.12233 -0.127714 -0.219526 -0.142772 -0.0253159 0.0536809 0.080125 0.0571016 -0.0183211 -0.134724 -0.218347 -0.140698 0.10537 0.230635 -0.0169741 -0.260275 0.00868827 0.306062 -0.0323035 -0.400665 -0.252906 -0.23893 -0.40434 -0.0688019 0.303902 0.0393269 -0.259455 -0.0416906 0.225524
0.0533261 0.221676 0.20588 0.0928744 -0.0169863 -0.083121 -0.10464 -0.0859165 -0.0229991 0.0845139 0.199903 0.226054 0.0700661 -0.185075 -0.235303 0.0671548 0.306034 0.0030983 -0.399433 -0.020757 0.815719 0.855083 0.0350084 -0.404293 -0.0323252 0.304047 0.0929697 -0.224541 -0.198835
-0.251975 -0.233798 -0.13019 -0.0152352 0.0703864 0.11835 0.133806 0.120358 0.0748091 -0.00825682 -0.122061 -0.22853 -0.255028 -0.136944 0.097903 0.288106 0.247787 -0.0273998 -0.313828 -0.4047 -0.353264 -0.349855 -0.403561 -0.327222 -0.0497156 0.234981 0.294129 0.11465 -0.123547
0.0890404 0.00957584 -0.0584588 -0.110278 -0.145276 -0.165114 -0.171728 -0.165964 -0.14708 -0.113177 -0.062519 0.00451035 0.0834519 0.16497 0.23634 0.285407 0.305513 0.298672 0.274919 0.248139 0.230823 0.230209 0.246547 0.273036 0.297446 0.305825 0.287718 0.240496 0.170292
We superpose a surface plot of $f$ on top of the grid:
X = [sinpi(θ)*cospi(φ) for θ in θ, φ in φ]
Y = [sinpi(θ)*sinpi(φ) for θ in θ, φ in φ]
Z = [cospi(θ) for θ in θ, φ in φ]
scatter3d(vec(X), vec(Y), vec(Z); markersize=1.25, markercolor=:violetred)
surface!(X, Y, Z; surfacecolor=F, legend=false, xlabel="x", ylabel="y", zlabel="f")
savefig(joinpath(GENFIGS, "sphere1.html"))
"/home/runner/work/FastTransforms.jl/FastTransforms.jl/docs/src/generated/sphere1.html"
We show the cut in the surface to help illustrate the definition of the grid. In particular, we do not sample the poles.
We precompute a spherical harmonic–Fourier plan:
P = plan_sph2fourier(F)
FastTransforms Spherical harmonic--Fourier plan for 15×29-element array of Float64
And an FFTW Fourier analysis plan on $\mathbb{S}^2$:
PA = plan_sph_analysis(F)
FastTransforms plan for FFTW Fourier analysis on the sphere for 15×29-element array of Float64
Its spherical harmonic coefficients demonstrate that it is exact-degree-$n$:
V = PA*F
U = threshold!(P\V, 400*eps())
15×29 Matrix{Float64}:
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 -1.66397e-5 2.73274e-5
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 -0.000271538 -7.91038e-5 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 3.29805e-5 -0.00164832 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.00703526 -0.00174747 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0122275 0.0220184 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 -0.0503776 0.0488881 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 -0.134224 -0.0799001 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0729185 -0.259962 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.337609 0.00337719 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0610941 0.235853 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0289877 -0.0164793 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.190197 0.192109 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 -0.0788841 0.135585 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.188313 0.0507949 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.142187 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
The $L^2(\mathbb{S}^2)$ norm of the function is:
nrm1 = norm(U)
0.6582728344942347
Similarly, on the tensor product grid, our function samples are:
Pnxy = FastTransforms.clenshaw!(c, A, B, C, [x⋅y], [1.0], [0.0])[1]
F = [(F[n, m] - Pnxy)/(z(θ[n], φ[m])⋅y - x⋅y) for n in 1:N, m in 1:M]
15×29 Matrix{Float64}:
5.19948 4.13834 2.99067 1.90658 1.02691 0.455272 0.249877 0.429314 0.977414 1.83893 2.91328 4.06148 5.13287 6.0026 6.60156 6.92544 7.02254 6.96965 6.84793 6.72603 6.65164 6.64905 6.71906 6.83899 6.96269 7.02276 6.93845 6.63197 6.05227
5.95491 2.78879 -1.44975 -4.81316 -5.93907 -5.32122 -4.74802 -5.26014 -5.93537 -4.96628 -1.72861 2.51536 5.80241 7.01649 6.35996 4.76275 3.10172 1.83264 1.03553 0.610092 0.43191 0.426508 0.591355 0.996377 1.76414 3.00054 4.64424 6.26964 7.0226
6.89687 5.65236 0.11585 -5.28724 -4.64181 1.51819 5.42995 1.95544 -4.33608 -5.4796 -0.331651 5.39124 6.95616 4.78949 1.90996 0.265776 -0.0679633 0.206459 0.546026 0.760489 0.851318 0.854026 0.770114 0.565133 0.231314 -0.0604333 0.208258 1.74968 4.58324
3.06027 6.5357 5.89997 -0.165743 -5.46592 -5.06636 -3.28816 -4.90175 -5.62054 -0.644121 5.6253 6.677 3.31846 0.440384 -0.00401181 0.652985 1.01303 0.876591 0.580633 0.35795 0.254199 0.250991 0.347269 0.561893 0.857659 1.01648 0.69514 0.0273999 0.337402
-0.070137 1.60923 5.43093 6.93528 3.91137 -0.131309 -1.75656 -0.336142 3.61343 6.85843 5.66343 1.83509 -0.053904 0.487406 1.01741 0.644213 0.136921 0.0101545 0.134135 0.275013 0.346961 0.349218 0.282296 0.144972 0.0138328 0.115298 0.603776 1.01615 0.548272
0.945078 0.135483 0.27158 2.71505 5.55048 6.82449 7.00499 6.85725 5.69314 2.92844 0.369819 0.0844295 0.912316 0.807551 0.149683 0.0509116 0.357273 0.499174 0.407141 0.269902 0.193177 0.1907 0.262329 0.397299 0.498435 0.37596 0.0666419 0.119637 0.764012
0.3347 0.96741 0.756526 0.0516758 0.0403323 0.667862 1.0399 0.710559 0.0718486 0.0205398 0.70939 0.987443 0.382027 0.0075432 0.314675 0.493074 0.249535 0.0255986 0.0288774 0.124344 0.187529 0.189637 0.130366 0.0344678 0.0191719 0.228563 0.486554 0.340488 0.0125441
0.19163 0.0244645 0.474228 0.953439 0.971777 0.756739 0.645454 0.743468 0.95977 0.969839 0.515354 0.0373859 0.165042 0.488731 0.32132 0.0219213 0.0963195 0.301654 0.329627 0.240493 0.171995 0.169627 0.234193 0.325354 0.310183 0.111448 0.0141702 0.296273 0.495085
0.41699 0.454012 0.154439 0.00689009 0.147415 0.349199 0.434542 0.359728 0.162073 0.00898106 0.134876 0.43942 0.434649 0.100528 0.0376636 0.284172 0.302882 0.085971 -0.0128005 0.0446154 0.109712 0.112106 0.0501911 -0.0122374 0.0727995 0.291434 0.297253 0.0501979 0.0815347
0.00653828 0.144506 0.416175 0.496049 0.398793 0.286378 0.242449 0.28086 0.390091 0.493228 0.429442 0.163601 0.00423886 0.186207 0.340647 0.145106 -0.0125214 0.139641 0.288402 0.253913 0.179097 0.176118 0.248108 0.291182 0.154331 -0.0106262 0.126633 0.337571 0.204303
0.334489 0.165749 0.014265 0.0385386 0.156172 0.253936 0.289478 0.258454 0.164306 0.0451436 0.010362 0.151498 0.32921 0.257363 0.0225829 0.0493168 0.268651 0.220231 -0.0026208 -0.0383697 0.0440181 0.047986 -0.0336588 -0.0128432 0.205569 0.277144 0.0634765 0.012466 0.242082
0.0666146 0.257095 0.341563 0.293211 0.203674 0.139862 0.117949 0.137043 0.198126 0.287354 0.341672 0.267702 0.0789497 -0.00992474 0.150427 0.293905 0.12293 -0.051658 0.138809 0.33873 0.255252 0.247571 0.340446 0.15898 -0.0503184 0.104503 0.292551 0.165851 -0.00643622
0.106032 -0.00382806 0.00704658 0.088086 0.17021 0.221457 0.238521 0.223661 0.174806 0.0942214 0.0112304 -0.00679647 0.0953468 0.252997 0.275056 0.0878466 -0.0516354 0.119191 0.338121 0.128404 -0.323093 -0.344224 0.0980961 0.340427 0.138821 -0.0504052 0.0724608 0.267888 0.261055
0.286625 0.280994 0.220196 0.149338 0.0951419 0.0642734 0.0542416 0.0629728 0.0923113 0.144964 0.215269 0.278062 0.288893 0.21207 0.0695276 -0.0416382 -0.0180382 0.136094 0.292814 0.340282 0.310784 0.308881 0.339519 0.300016 0.148439 -0.0107592 -0.0450601 0.0595874 0.203734
0.0747987 0.122395 0.163672 0.195512 0.217272 0.229726 0.233902 0.230262 0.2184 0.197306 0.166153 0.125451 0.0781262 0.0298766 -0.0118277 -0.0400996 -0.0514099 -0.0472393 -0.0335309 -0.0182388 -0.00839924 -0.0080514 -0.0173325 -0.0324517 -0.0465243 -0.0515701 -0.0414172 -0.0142377 0.0267486
We superpose a surface plot of $f$ on top of the grid:
scatter3d(vec(X), vec(Y), vec(Z); markersize=1.25, markercolor=:violetred)
surface!(X, Y, Z; surfacecolor=F, legend=false, xlabel="x", ylabel="y", zlabel="f")
savefig(joinpath(GENFIGS, "sphere2.html"))
"/home/runner/work/FastTransforms.jl/FastTransforms.jl/docs/src/generated/sphere2.html"
Its spherical harmonic coefficients demonstrate that it is degree-$(n-1)$:
V = PA*F
U = threshold!(P\V, 400*eps())
15×29 Matrix{Float64}:
2.07907 0.87453 0.235893 0.069594 -0.119618 0.0581001 0.0586842 0.0859934 -0.0488866 -0.015775 -0.060899 -0.0291623 -0.000291717 -0.00253358 0.00903247 0.000694039 0.000413142 -0.00117608 0.0011413 0.000780662 0.00140576 0.000942806 -0.000234181 8.63587e-6 -0.00043161 -0.000121732 -3.54626e-5 0.0 0.0
1.51317 0.44642 0.120416 -0.169934 0.292081 0.328839 0.332145 0.288911 -0.164244 -0.0423205 -0.163378 -0.0586349 -0.000586539 -0.00145765 0.00519666 -0.00998075 -0.00594126 -0.00862106 0.00836616 0.00369554 0.00665466 0.00342882 -0.000851678 2.16057e-5 -0.00107982 0.0 0.0 0.0 0.0
0.36606 -0.627659 -0.169303 -0.592046 1.0176 0.703541 0.710614 0.503834 -0.286426 -0.056082 -0.216504 -0.0224552 -0.000224625 0.0140577 -0.0501173 -0.0493515 -0.0293776 -0.0276576 0.0268398 0.00914327 0.0164645 0.00585456 -0.0014542 0.0 0.0 0.0 0.0 0.0 0.0
-0.225238 -1.3923 -0.375555 -0.884141 1.51965 0.893241 0.902221 0.498016 -0.283118 -0.0162747 -0.0628284 0.165928 0.00165982 0.0536982 -0.191439 -0.123087 -0.0732701 -0.0534663 0.0518855 0.0122482 0.0220556 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0279478 -1.24976 -0.337105 -0.80053 1.37594 0.669715 0.676448 0.113165 -0.0643333 0.096013 0.370657 0.512707 0.00512874 0.109365 -0.389897 -0.195697 -0.116493 -0.0592355 0.0574841 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.642706 -0.460978 -0.124343 -0.408306 0.701792 0.115506 0.116667 -0.53107 0.30191 0.2422 0.935011 0.866631 0.00866913 0.145962 -0.520371 -0.183527 -0.109249 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.924445 0.171516 0.0462641 -0.0394459 0.0677993 -0.395183 -0.399156 -1.05835 0.601662 0.335991 1.29709 0.970497 0.00970813 0.116526 -0.415426 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.619857 0.0934262 0.0252005 0.0128822 -0.0221418 -0.505099 -0.510177 -1.1197 0.636541 0.30573 1.18027 0.65018 0.00650392 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0816312 -0.512221 -0.138165 -0.239451 0.411566 -0.194804 -0.196763 -0.706203 0.401471 0.16128 0.622619 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
-0.15003 -0.975789 -0.263206 -0.507765 0.872741 0.197111 0.199092 -0.186337 0.105931 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.105958 -0.834521 -0.225101 -0.515565 0.886148 0.287159 0.290046 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.530802 -0.266039 -0.0717606 -0.261012 0.448626 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.668154 0.140073 0.0377829 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.390394 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
Finally, the Legendre polynomial $P_n(z\cdot x)$ is aligned with the grid:
pts = vec([z(θ, φ)⋅x for θ in θ, φ in φ])
F = reshape(FastTransforms.clenshaw!(c, A, B, C, pts, phi0, zeros(N*M)), N, M)
15×29 Matrix{Float64}:
0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808
-0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968
0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489
-0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835
0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901
-0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532
0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677
-0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473 -0.209473
0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677 0.210677
-0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532 -0.214532
0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901 0.221901
-0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835 -0.234835
0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489 0.258489
-0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968 -0.30968
0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808 0.501808
We superpose a surface plot of $f$ on top of the grid:
scatter3d(vec(X), vec(Y), vec(Z); markersize=1.25, markercolor=:violetred)
surface!(X, Y, Z; surfacecolor=F, legend=false, xlabel="x", ylabel="y", zlabel="f")
savefig(joinpath(GENFIGS, "sphere3.html"))
"/home/runner/work/FastTransforms.jl/FastTransforms.jl/docs/src/generated/sphere3.html"
It only has one nonnegligible spherical harmonic coefficient. Can you spot it?
V = PA*F
U = threshold!(P\V, 400*eps())
15×29 Matrix{Float64}:
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
0.658273 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
That nonnegligible coefficient should be
ret = eval("√(2π/($(N-1)+1/2))")
"√(2π/(14+1/2))"
which is approximately
eval(Meta.parse(ret))
0.6582728344942353
since the convention in this library is to orthonormalize.
nrm2 = norm(U)
0.6582728344942352
Note that the integrals of both functions $P_n(z\cdot y)$ and $P_n(z\cdot x)$ and their $L^2(\mathbb{S}^2)$ norms are the same because of rotational invariance. The integral of either is perhaps not interesting as it is mathematically zero, but the norms of either should be approximately the same.
nrm1 ≈ nrm2
true
This page was generated using Literate.jl.