Compute Taylor coefficients for conformal map

This example demonstrates how to reproduce the the Taylor coefficients for the conformal map $Z \to W$ and also of the inverse map, $W \to Z$, where $Z = z^4$ and $W = w^3$.

The algorithm to obtain the coefficients follows the procedure described in the Appendix of the paper by Rančić et al. [1]. An elaborate explanation of the algorithm is documented online. Furthermore, the repository github.com/kburns/LightningCubedSphere.jl provides an independent algorithm using least squares fitting and rational functions ("lightning" rational approximations) to compute the forward and inverse maps.

This example, follows the steps described in the notes and showcases how to obtain the coefficients $A_k$ of the Taylor series

\[W(Z) = \sum_{k=1}^\infty A_k Z^k\]

and also coefficients $B_k$ of the inverse Taylor series

\[Z(W) = \sum_{k=1}^\infty B_k Z^k\]

We start by importing some useful packages and define some utility methods.

using FFTW, SpecialFunctions, TaylorSeries

function find_angles(φ)
    φ⁻ = -φ
    φ⁺ = +φ

    w⁻ = cis(φ⁻)
    w⁺ = cis(φ⁺)

    w′⁻ = (1 - w⁻) / (1 + w⁻/2)
    w′⁺ = (1 - w⁺) / (1 + w⁺/2)

    φ′⁻ = angle(w′⁻)
    φ′⁺ = angle(w′⁺)

    return φ′⁻, φ′⁺
end
find_angles (generic function with 1 method)

Next we define the cubic roots and be careful to choose the appropriate branch. We use cbrt to go from $W$ to $w$ and cbrt′ to go from $W'$ to $w'$.

function Base.cbrt(z::Complex)
    r = abs(z)
    φ = angle(z) # ∈ [-π, +π]
    θ = φ / 3
    return cbrt(r) * cis(θ)
end

function cbrt′(z::Complex)
    φ′⁻, φ′⁺ = find_angles(π/3)
    r = abs(z)
    φ = angle(z) # ∈ [-π, +π]

    θ = φ / 3

    if 0 < θ ≤ φ′⁻
        θ -= 2π/3
    elseif φ′⁺ ≤ θ < 0
        θ += 2π/3
    end
    return r^(1/3) * cis(θ)
end
cbrt′ (generic function with 1 method)

A few more utility methods and the actual iterative algorithm.

"""
    find_N(r; decimals=15)

Return the required number of points we need to consider around the circle of
radius `r` to compute the conformal map series coefficients up to `decimals`
points. The number of points is computed based on the estimate of eq. (B9) in
the paper by [Rancic-etal-1996](@citet). That is, is the smallest integer ``N`` (which
is chosen to be a power of 2 so that the FFTs are efficient) for which
```math
N - \frac{7}{12} \frac{\\mathrm{log}_{10}(N)}{\\mathrm{log}_{10}(r)} - \frac{r + \\mathrm{log}_{10}(A₁ / C)}{-4 \\mathrm{log}_{10}(r)} > 0
```
where ``r`` is the number of `decimals` we are aiming for and
```math
C = \frac{\\sqrt{3} \\Gamma(1/3) A_1^{1/3}}{256^{1/3} π}
```
with ``A_1`` an estimate of the ``Z^1`` Taylor series coefficient of ``W(Z)``.

For ``A_1 \\approx 1.4771`` we get ``C \\approx 0.265``.

#References

* [Rancic-etal-1996](@cite) Rančić et al., *Q. J. R. Meteorol.*, (1996).
"""
function find_N(r; decimals=15)
    A₁ = 1.4771 # an approximation of the first coefficient
    C = sqrt(3) * gamma(1/3) * A₁^(1/3) / (256^(1/3) * π)

    N = 2
    while N + 7/12 * log10(N) / (-log10(r)) - (decimals + log10(A₁ / C)) / (-4 * log10(r)) < 0
        N *= 2
    end

    return N
end

function _update_coefficients!(A, r, Nφ)
    Ncoefficients = length(A)

    Lφ = π/2
    dφ = Lφ / Nφ
    φ = range(-Lφ/2 + dφ/2, stop=Lφ/2 - dφ/2, length=Nφ)

    z = @. r * cis(φ)

    W̃′ = 0z
    for k = Ncoefficients:-1:1
        @. W̃′ += A[k] * (1 - z)^(4k)
    end

    w̃′ = @. cbrt′(W̃′)
    w̃  = @. (1 - w̃′) / (1 + w̃′/2)
    W̃  = @. w̃^3

    k = collect(fftfreq(Nφ, Nφ))
    g̃ = fft(W̃) ./ (Nφ * cis.(k * 4φ[1])) # divide with Nϕ * exp(-ikπ) to account for FFT's normalization
    g̃ = g̃[2:Ncoefficients+1]             # drop coefficient for 0-th power

    A .= [real(g̃[k] / r^(4k)) for k in 1:Ncoefficients]

    return nothing
end

"""
    find_taylor_coefficients(r = 1 - 1e-7;
                             Niterations = 30,
                             maximum_coefficients = 256,
                             Nevaluations = find_N(r; decimals=15))

Return the Taylor coefficients for the conformal map ``Z \to W`` and also of the
inverse map, ``W \to Z``, where ``Z = z^4`` and ``W = w^3``. In particular, it returns the
coefficients ``A_k`` of the Taylor series

```math
W(Z) = \\sum_{k=1}^\\infty A_k Z^k
```

and also coefficients ``B_k`` the inverse Taylor series

```math
Z(W) = \\sum_{k=1}^\\infty B_k Z^k
```

The algorithm to obtain the coefficients follows the procedure described in the
Appendix of the paper by [Rancic-etal-1996](@citet)

Arguments
=========

* `r` (positional): the radius about the center and the edge of the cube used in the
  algorithm described by [Rancic-etal-1996](@citet). `r` must be less than 1; default: 1 - 10``^{-7}``.

* `maximum_coefficients` (keyword): the truncation for the Taylor series; default: 256.

* `Niterations` (keyword): the number of update iterations we perform on the
  Taylor coefficients ``A_k``; default: 30.

* `Nevaluations` (keyword): the number of function evaluations in over the circle of radius `r`;
  default `find_N(r; decimals=15)`; see [`find_N`](@ref).

Example
=======

```@example
julia> A, B = find_taylor_coefficients(1 - 1e-4);
[ Info: Computing the first 256 coefficients of the Taylor serieses
[ Info: using 32768 function evaluations on a circle with radius 0.9999.
[ Info: Iteration 1 out of 30
[ Info: Iteration 2 out of 30
[ Info: Iteration 3 out of 30
[ Info: Iteration 4 out of 30
[ Info: Iteration 5 out of 30
[ Info: Iteration 6 out of 30
[ Info: Iteration 7 out of 30
[ Info: Iteration 8 out of 30
[ Info: Iteration 9 out of 30
[ Info: Iteration 10 out of 30
[ Info: Iteration 11 out of 30
[ Info: Iteration 12 out of 30
[ Info: Iteration 13 out of 30
[ Info: Iteration 14 out of 30
[ Info: Iteration 15 out of 30
[ Info: Iteration 16 out of 30
[ Info: Iteration 17 out of 30
[ Info: Iteration 18 out of 30
[ Info: Iteration 19 out of 30
[ Info: Iteration 20 out of 30
[ Info: Iteration 21 out of 30
[ Info: Iteration 22 out of 30
[ Info: Iteration 23 out of 30
[ Info: Iteration 24 out of 30
[ Info: Algorithm converged after 24 iterations

julia> A[1:10]
10-element Vector{Float64}:
  1.4771306289227293
 -0.3818351018795475
 -0.05573057838030261
 -0.008958833150428962
 -0.007913155711663374
 -0.004866251689037038
 -0.003292515242976284
 -0.0023548122712604494
 -0.0017587029515141275
 -0.0013568087584722149
```

!!! info "Reproducing coefficient table by Rančić et al. (1996)"
    To reproduce the coefficients tabulated by [Rancic-etal-1996](@citet) use
    the default values, i.e., ``r = 1 - 10^{-7}``.

#References

* [Rancic-etal-1996](@cite) Rančić et al., *Q. J. R. Meteorol.*, (1996).
"""
function find_taylor_coefficients(r = 1 - 1e-7;
                                  Niterations = 30,
                                  maximum_coefficients = 256,
                                  Nevaluations = find_N(r; decimals=15))

    (r < 0 || r ≥ 1) && error("r needs to be within 0 < r < 1")

    Ncoefficients = Int(Nevaluations/2) - 2 > maximum_coefficients ? maximum_coefficients : Int(Nevaluations/2) - 2

    @info "Computing the first $Ncoefficients coefficients of the Taylor serieses"
    @info "using $Nevaluations function evaluations on a circle with radius $r."

    # initialize coefficients
    A_coefficients = rand(Ncoefficients)

    # we can use Rančić as initial guess; this makes convergence faster
    # A_coefficients[1:min(maximum_coefficients, 30)] = CubedSphere.A_Rancic[2:min(maximum_coefficients, 30)+1]

    A_coefficients_old = deepcopy(A_coefficients)

    for iteration in 1:Niterations
        @info "Iteration $iteration out of $Niterations"

        _update_coefficients!(A_coefficients, r, Nevaluations)

        rel_error = (abs.(A_coefficients - A_coefficients_old) ./ abs.(A_coefficients))[1:maximum_coefficients]

        A_coefficients_old .= A_coefficients

        if all(rel_error .< 1e-15)
            iteration_str = iteration == 1 ? "iteration" : "iterations"
            @info "Algorithm converged after $iteration $iteration_str"
            break
        end
    end

    # convert to Taylor series; add coefficient for 0-th power
    A_series = Taylor1([0; A_coefficients])

    B_series = inverse(A_series) # This is the inverse Taylor series
    B_series.coeffs[1] !== 0.0 && error("coefficient that corresponds to W^0 is non-zero; something went wrong")
    B_coefficients = B_series.coeffs[2:end] # don't return coefficient for 0-th power

    return A_coefficients, B_coefficients
end
Main.var"##225".find_taylor_coefficients

Now let's reproduce the results by Rančić et al. [1]. For that, we need to choose $r = 1 - 10^{-7}$.

r = 1 - 1e-7
A_coefficients, B_coefficients = find_taylor_coefficients(r)
([1.4771306260096393, -0.38183510510173474, -0.05573058001191046, -0.008958836068176124, -0.007913157852211802, -0.004866254377075463, -0.0032925175127914356, -0.002354814883250576, -0.0017587052747532216, -0.0013568113327830716, -0.001074598476991445, -0.0008694447594809179, -0.0007160711512062413, -0.0005986710009288835, -0.0005069906323874284, -0.00043415191278705083, -0.0003754100328585039, -0.00032741060099483667, -0.00028773091481459136, -0.00025458777518603393, -0.0002266464237139699, -0.00020289261022324973, -0.0001825451083023894, -0.00016499474459965587, -0.00014976117167335968, -0.00013646173947213006, -0.00012478875821918617, -0.00011449267278843712, -0.00010536946150265186, -9.72510937545457e-5, -8.999822956766894e-5, -8.349458099897983e-5, -7.764251836049124e-5, -7.235961803524566e-5, -6.757592905242176e-5, -6.323179302866274e-5, -5.927609356945838e-5, -5.5664841489802644e-5, -5.236002450166198e-5, -4.9328666575981884e-5, -4.654205459455576e-5, -4.397509927726449e-5, -4.1605804499664226e-5, -3.941482457921214e-5, -3.7385093323059846e-5, -3.550151190273581e-5, -3.375068517737713e-5, -3.212069809550705e-5, -3.060092539229762e-5, -2.9181869059567485e-5, -2.7855019071828996e-5, -2.6612733658798598e-5, -2.5448136065130857e-5, -2.435502526456832e-5, -2.332779852365339e-5, -2.2361384059483585e-5, -2.1451182322214663e-5, -2.059301466841553e-5, -1.9783078385925436e-5, -1.9017907191642877e-5, -1.8294336457795284e-5, -1.7609472533707903e-5, -1.696066562354835e-5, -1.6345485758904326e-5, -1.5761701471001904e-5, -1.5207260822981516e-5, -1.4680274510013218e-5, -1.4179000774514827e-5, -1.3701831918432124e-5, -1.3247282222904689e-5, -1.2813977110917616e-5, -1.240064340950306e-5, -1.2006100586490526e-5, -1.162925285244488e-5, -1.1269082031983413e-5, -1.0924641120429924e-5, -1.0595048451915383e-5, -1.0279482413770775e-5, -9.977176649846514e-6, -9.687415701957124e-6, -9.40953104456613e-6, -9.142897472854458e-6, -8.886929808817616e-6, -8.64107989395892e-6, -8.404833840552195e-6, -8.177709516519166e-6, -7.959254241596947e-6, -7.749042674869538e-6, -7.546674875796755e-6, -7.35177452271506e-6, -7.163987274447959e-6, -6.982979262099635e-6, -6.808435699400626e-6, -6.640059601122039e-6, -6.477570600128959e-6, -6.320703854499292e-6, -6.1692090370322125e-6, -6.022849400131344e-6, -5.88140090975986e-6, -5.744651442703552e-6, -5.612400041958636e-6, -5.484456225484151e-6, -5.360639344047111e-6, -5.240777984189024e-6, -5.124709412797815e-6, -5.012279059976641e-6, -4.903340037266291e-6, -4.797752688485152e-6, -4.695384170705692e-6, -4.596108063083692e-6, -4.49980400146898e-6, -4.406357336865874e-6, -4.315658816006367e-6, -4.227604282406696e-6, -4.14209439645335e-6, -4.05903437311164e-6, -3.978333736062199e-6, -3.8999060870459756e-6, -3.823668889416159e-6, -3.7495432648563915e-6, -3.6774538024125214e-6, -3.6073283789456564e-6, -3.5390979902983685e-6, -3.472696592386658e-6, -3.4080609516129257e-6, -3.3451305039539503e-6, -3.283847222172259e-6, -3.2241554906024696e-6, -3.166001987073713e-6, -3.1093355713865024e-6, -3.054107180197992e-6, -3.0002697274986304e-6, -2.9477780107521843e-6, -2.8965886220645523e-6, -2.8466598642074285e-6, -2.7979516711278777e-6, -2.7504255327365284e-6, -2.7040444236573376e-6, -2.6587727357360968e-6, -2.614576214082776e-6, -2.5714218964230907e-6, -2.5292780555851013e-6, -2.4881141449131714e-6, -2.4479007464727145e-6, -2.4086095218603583e-6, -2.370213165478379e-6, -2.332685360137418e-6, -2.2960007348507434e-6, -2.2601348247135713e-6, -2.225064032718586e-6, -2.1907655934424145e-6, -2.157217538465465e-6, -2.1243986634574264e-6, -2.092288496816828e-6, -2.0608672697999724e-6, -2.030115888040815e-6, -2.0000159044089475e-6, -1.9705494931200033e-6, -1.9416994250357375e-6, -1.913449044097987e-6, -1.8857822448291563e-6, -1.858683450855371e-6, -1.8321375943880216e-6, -1.8061300966244031e-6, -1.7806468490143896e-6, -1.7556741953593746e-6, -1.7311989146894875e-6, -1.7072082048861426e-6, -1.6836896670225374e-6, -1.6606312903681602e-6, -1.6380214380420707e-6, -1.615848833269593e-6, -1.594102546232839e-6, -1.5727719814556685e-6, -1.5518468657402953e-6, -1.5313172365803208e-6, -1.5111734310707121e-6, -1.4914060752721607e-6, -1.4720060739950955e-6, -1.4529646010213127e-6, -1.4342730897036566e-6, -1.4159232239529262e-6, -1.397906929582434e-6, -1.3802163659971953e-6, -1.3628439182081625e-6, -1.345782189170577e-6, -1.3290239924115658e-6, -1.3125623449522543e-6, -1.2963904604985074e-6, -1.2805017428962509e-6, -1.264889779842358e-6, -1.2495483368256495e-6, -1.2344713512945302e-6, -1.2196529271094517e-6, -1.205087329032529e-6, -1.1907689776553116e-6, -1.176692444338419e-6, -1.1628524464034445e-6, -1.1492438425049063e-6, -1.135861628157257e-6, -1.1227009314296771e-6, -1.109757008801272e-6, -1.0970252411546346e-6, -1.0845011299267541e-6, -1.072180293388257e-6, -1.0600584630590202e-6, -1.0481314802454338e-6, -1.0363952927018524e-6, -1.0248459514242113e-6, -1.013479607521337e-6, -1.002292509227999e-6, -9.912809990049161e-7, -9.80441510742953e-7, -9.697705670568836e-7, -9.592647766824409e-7, -9.48920831949032e-7, -9.387355063509303e-7, -9.287056521853151e-7, -9.188281982831935e-7, -9.091001478068767e-7, -8.995185761229876e-7, -8.900806287461875e-7, -8.80783519353514e-7, -8.716245278509693e-7, -8.626009985241213e-7, -8.53710338227849e-7, -8.449500146462649e-7, -8.363175546072889e-7, -8.278105424457603e-7, -8.194266184265126e-7, -8.111634772112687e-7, -8.030188663748191e-7, -7.949905849684276e-7, -7.870764821339853e-7, -7.792744557497762e-7, -7.715824511268408e-7, -7.639984597472419e-7, -7.565205180272074e-7, -7.49146706140871e-7, -7.418751468565845e-7, -7.347040044276617e-7, -7.276314834996564e-7, -7.206558280659872e-7, -7.137753204477404e-7, -7.069882802973736e-7, -7.002930636500063e-7, -6.93688061980744e-7, -6.871717013084617e-7, -6.807424413162179e-7, -6.743987744989821e-7, -6.681392253350889e-7, -6.619623494907355e-7, -6.55866733032405e-7, -6.498509916726902e-7, -6.439137700359888e-7, -6.38053740941834e-7], [0.6769881975173903, 0.11847293456554321, 0.05317178134667832, 0.02965810434051842, 0.019124473040276056, 0.01342565621117347, 0.009988733231804128, 0.007748689964055174, 0.00620346979888045, 0.0050901087488285975, 0.004259811843277695, 0.003623089560765813, 0.003123414689401724, 0.0027236094894184074, 0.0023983808655506857, 0.0021300190511767425, 0.0019058131613060752, 0.0017164415640357968, 0.001554937682548106, 0.0014160071520703192, 0.0012955659775392616, 0.0011904214022623186, 0.0010980471179021772, 0.0010164221662771601, 0.0009439136652232816, 0.0008791902122376714, 0.0008211571031080119, 0.0007689072877511454, 0.0007216838296902914, 0.0006788508775020999, 0.0006398710059121416, 0.0006042873734597352, 0.000571709557900082, 0.0005418022253450302, 0.0005142760014104456, 0.0004888800671092897, 0.0004653961157768634, 0.0004436333915709468, 0.00042342459316764943, 0.00040462247389303467, 0.0003870970057640944, 0.00037073300268894124, 0.0003554281195199249, 0.00034109116031953035, 0.00032764064223438744, 0.0003150035716301939, 0.00030311439725837593, 0.00029191411168540924, 0.0002813494773834344, 0.00027137235803536083, 0.00026193913896368636, 0.00025301022331564914, 0.0002445495928570399, 0.00023652442404379213, 0.00022890475153347685, 0.00022166317253051845, 0.00021477458637883951, 0.00020821596466318932, 0.00020196614778715904, 0.00019600566458719294, 0.0001903165720382067, 0.0001848823125242738, 0.00017968758650069427, 0.0001747182386725781, 0.0001699611560688461, 0.00016540417660663997, 0.0001610360069256186, 0.00015684614842951752, 0.0001528248306078242, 0.00014896295082693448, 0.0001452520198805875, 0.00014168411267614312, 0.0001382518235083838, 0.0001349482254376983, 0.00013176683334616382, 0.00012870157029439908, 0.00012574673684513643, 0.00012289698305712778, 0.00012014728288599328, 0.00011749291075757598, 0.00011492942010481859, 0.00011245262368159099, 0.00011005857548666311, 0.00010774355414848041, 0.00010550404763684836, 0.00010333673918132214, 0.00010123849428824722, 9.920634875919123e-5, 9.723749762311789e-5, 9.532928490321454e-5, 9.347919414692738e-5, 9.168483965458498e-5, 8.99439583481001e-5, 8.825440222671164e-5, 8.661413136163685e-5, 8.502120738591082e-5, 8.347378743965378e-5, 8.197011853457176e-5, 8.050853230471133e-5, 7.908744011338703e-5, 7.770532848881888e-5, 7.63607548633831e-5, 7.505234359352196e-5, 7.377878223929745e-5, 7.253881808433366e-5, 7.133125487848844e-5, 7.015494978704604e-5, 6.900881053154131e-5, 6.789179270852632e-5, 6.680289727368426e-5, 6.574116817969237e-5, 6.470569015714638e-5, 6.369558662868876e-5, 6.271001774724474e-5, 6.174817854996368e-5, 6.08092972201004e-5, 5.989263344965469e-5, 5.8997476896121244e-5, 5.812314572719385e-5, 5.726898524771792e-5, 5.643436660360022e-5, 5.5618685557765364e-5, 5.482136133359963e-5, 5.4041835521645644e-5, 5.3279571045609456e-5, 5.253405118401552e-5, 5.18047786440996e-5, 5.10912746847628e-5, 5.039307828562665e-5, 4.970974535942922e-5, 4.904084800518733e-5, 4.838597379972111e-5, 4.77447251252962e-5, 4.711671853128588e-5, 4.6501584127891924e-5, 4.589896501008968e-5, 4.5308516710080634e-5, 4.4729906676644663e-5, 4.416281377988616e-5, 4.360692783996248e-5, 4.3061949178471085e-5, 4.252758819125378e-5, 4.200356494145253e-5, 4.1489608771722616e-5, 4.0985457934575e-5, 4.0490859239881726e-5, 4.000556771863587e-5, 3.952934630211154e-5, 3.90619655156197e-5, 3.860320318610265e-5, 3.8152844162854226e-5, 3.77106800506938e-5, 3.727650895496072e-5, 3.68501352377322e-5, 3.643136928470099e-5, 3.602002728218158e-5, 3.5615931003742765e-5, 3.521890760599265e-5, 3.4828789433068347e-5, 3.4445413829407014e-5, 3.406862296039821e-5, 3.369826364053897e-5, 3.333418716873353e-5, 3.297624917039908e-5, 3.262430944605618e-5, 3.227823182610042e-5, 3.1937884031467146e-5, 3.1603137539916436e-5, 3.127386745767972e-5, 3.094995239622259e-5, 3.063127435389122e-5, 3.0317718602221373e-5, 3.0009173576700505e-5, 2.9705530771783866e-5, 2.940668463997542e-5, 2.911253249479415e-5, 2.8822974417454778e-5, 2.853791316710082e-5, 2.8257254094435407e-5, 2.7980905058603264e-5, 2.7708776347183984e-5, 2.7440780599163712e-5, 2.7176832730758693e-5, 2.6916849863970025e-5, 2.6660751257754955e-5, 2.6408458241705156e-5, 2.6159894152127843e-5, 2.5914984270430293e-5, 2.5673655763712937e-5, 2.5435837627480712e-5, 2.52014606303864e-5, 2.4970457260923553e-5, 2.4742761675990642e-5, 2.451830965125124e-5, 2.429703853321872e-5, 2.4078887192997004e-5, 2.3863795981611998e-5, 2.3651706686871294e-5, 2.3442562491692293e-5, 2.3236307933841822e-5, 2.3032888867032508e-5, 2.283225242332373e-5, 2.2634346976777212e-5, 2.243912210831937e-5, 2.2246528571764652e-5, 2.2056518260956134e-5, 2.1869044177981274e-5, 2.1684060402422726e-5, 2.1501522061605686e-5, 2.1321385301804827e-5, 2.1143607260375554e-5, 2.0968146038775592e-5, 2.0794960676444436e-5, 2.062401112550953e-5, 2.0455258226289225e-5, 2.028866368356389e-5, 2.0124190043587566e-5, 1.9961800671813885e-5, 1.9801459731310772e-5, 1.9643132161839664e-5, 1.9486783659575833e-5, 1.9332380657447373e-5, 1.9179890306071305e-5, 1.9029280455266026e-5, 1.8880519636120214e-5, 1.8733577043599097e-5, 1.8588422519669573e-5, 1.8445026536926604e-5, 1.8303360182703766e-5, 1.8163395143651688e-5, 1.8025103690768548e-5, 1.788845866486754e-5, 1.7753433462466706e-5, 1.762000202208708e-5, 1.7488138810945685e-5, 1.7357818812030372e-5, 1.7229017511543905e-5, 1.7101710886705345e-5, 1.6975875393896995e-5, 1.685148795714582e-5, 1.6728525956928485e-5, 1.6606967219289625e-5, 1.6486790005263343e-5, 1.6367973000588298e-5, 1.6250495305707004e-5, 1.6134336426040412e-5, 1.60194762625291e-5, 1.5905895102432707e-5, 1.5793573610379565e-5, 1.5682492819658715e-5, 1.5572634123746862e-5, 1.5463979268062938e-5, 1.535651034194339e-5, 1.525020977083131e-5, 1.5145060308673006e-5, 1.504104503051565e-5])

The first few coefficients are:

A_coefficients[1:10]
10-element Vector{Float64}:
  1.4771306260096393
 -0.38183510510173474
 -0.05573058001191046
 -0.008958836068176124
 -0.007913157852211802
 -0.004866254377075463
 -0.0032925175127914356
 -0.002354814883250576
 -0.0017587052747532216
 -0.0013568113327830716

And we can also plot them:

using GLMakie

fig = Figure(fontsize=20)
ax = Axis(fig[1, 1],
          xlabel = L"$k$ th coefficient",
          ylabel = L"|A_k|",
          xticks = 1:10,
          yscale = log10)

scatter!(ax, abs.(A_coefficients[1:10]))
current_figure()

This page was generated using Literate.jl.