Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Inverse Laplacian #141

Closed
maximilian-gelbrecht opened this issue Sep 9, 2022 · 7 comments · Fixed by #142
Closed

Inverse Laplacian #141

maximilian-gelbrecht opened this issue Sep 9, 2022 · 7 comments · Fixed by #142
Labels
enhancement 🔍 New feature or request
Milestone

Comments

@maximilian-gelbrecht
Copy link
Member

maximilian-gelbrecht commented Sep 9, 2022

I’d need the inverse laplacian to convert vorticity to streamfunction. I’ve seen that this is available in the code, but commented out in spectral_gradients! as ∇⁻²!. Why is that so? Could we make it functional again? SpectralTransform still has the eigenvalues⁻¹ field that I also only used by ∇⁻²!

UV_from_vor! also somehow inverts the laplacian but I don't understand how. Anyway though, I need to compute the stream function not the wind field.

@milankl
Copy link
Member

milankl commented Sep 9, 2022

Yes, we can make that functional again. UV_from_vor! or UV_from_vordiv! invert the Laplacian and take the gradients (nabla in the former, rotated nabla + nabla in the latter) in one go (which is faster and more memory efficient). In the following diagram you basically jump directly from vor (,div) to U,V without explicitly calculating stream function or velocity potential.

image

I don't think it's necessary to use precalculated eigenvalues, as they are so simple, so I'd probably drop ::SpectralTransform as an argument. Which methods do you need

  • ∇⁻²(::AbstractMatrix)
  • ∇⁻²(::LowerTriangularMatrix)
  • ∇⁻²!(::AbstractMatrix,::AbstractMatrix)
  • ∇⁻²!(::LowerTriangularMatrix,::LowerTriangularMatrix)
  • ...?

@milankl
Copy link
Member

milankl commented Sep 9, 2022

I know the documentation is far from up to date, but at least that is already described

https://milankl.github.io/SpeedyWeather.jl/dev/spectral_transform/#U,V-from-vorticity-and-divergence

@maximilian-gelbrecht
Copy link
Member Author

maximilian-gelbrecht commented Sep 9, 2022

it should also work to just extend ∇²! shouldn't it?

function ∇²!(   ∇²alms::LowerTriangularMatrix{Complex{NF}}, # Output: Laplacian of alms
                alms::LowerTriangularMatrix{Complex{NF}},   # Input: spectral coefficients
                S::SpectralTransform{NF};                   # precomputed eigenvalues
                add::Bool=false,                            # add to output array or overwrite
                flipsign::Bool=false,                       # -∇² or ∇²
                inverse::Bool=false                        # ∇⁻² or ∇²
                ) where {NF<:AbstractFloat}

    @boundscheck size(alms) == size(∇²alms) || throw(BoundsError)
    lmax,mmax = size(alms) .- (2,1)     # degree l, order m of the Legendre polynomials
    
    @unpack eigenvalues = S
    eigenvalues = inverse ? inv.(eigenvalues) : eigenvalues 
    @boundscheck length(eigenvalues) >= lmax+1 || throw(BoundsError)

    kernel(o,a) = flipsign ? (add ? (o-a) : -a) :
                             (add ? (o+a) :  a)    

    lm = 0
    @inbounds for m in 1:mmax+1     # order m = 0:mmax but 1-based
        for l in m:lmax+1           # degree l = m:lmax but 1-based
            lm += 1
            ∇²alms[lm] = kernel(∇²alms[lm],alms[lm]*eigenvalues[l])
        end
        lm += 1
    end

    return ∇²alms
end

I didn't test it, though

@milankl
Copy link
Member

milankl commented Sep 9, 2022

Yes, good idea, I didn't realise I had already written ∇²! with a kernel. With #142 I decided to externalise the kernel such that we can also define ∇⁻² which might be nice to avoid confusion? And I didn't like the inv.(eigenvalues) as it would allocate memory ;)

@milankl
Copy link
Member

milankl commented Sep 9, 2022

Oh we could also leave the kernel inside and define ∇⁻²!(...) = ∇²!(....;inverse=true) true yeah I don't know whether it makes a difference

@milankl
Copy link
Member

milankl commented Sep 9, 2022

Internalised kernel is committed now, interestingly for 3x ?-clauses in a function you need to @inline them otherwise you get a 10x performance penalty...

    @inline kernel(o,a,b) = inverse ?   (flipsign ? (add ? (o - a/b) : -a/b)  :
                                                    (add ? (o + a/b) :  a/b)) :
                                        (flipsign ? (add ? (o - a*b) : -a*b)  :
                                                    (add ? (o + a*b) :  a*b))

@milankl milankl added the enhancement 🔍 New feature or request label Sep 12, 2022
@milankl milankl added this to the v0.4 milestone Sep 12, 2022
@milankl
Copy link
Member

milankl commented Sep 14, 2022

For the record, writing a kernel like this

    # use eigenvalues⁻¹/eigenvalues for ∇⁻²/∇² based but name both eigenvalues
    eigenvalues = inverse ? S.eigenvalues⁻¹ : S.eigenvalues
    @boundscheck length(eigenvalues) >= lmax+1 || throw(BoundsError)

    @inline kernel(o,a) = flipsign ? (add ? (o-a) : -a)  :
                                     (add ? (o+a) :  a)

    lm = 0
    @inbounds for m in 1:mmax+1     # order m = 0:mmax but 1-based
        for l in m:lmax+1           # degree l = m:lmax but 1-based
            lm += 1
            ∇²alms[lm] = kernel(∇²alms[lm],alms[lm]*eigenvalues[l])
        end
    end

Gives the same speed as hardcoding +=, -=, = +, = - (and same for inverse). The @inline is technically not necessary, but using 3x ? will not inline that function automatically. I suspect as the speed is the same, this is actually directly compiled to the same code -- at least on my hardware.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
enhancement 🔍 New feature or request
Projects
None yet
Development

Successfully merging a pull request may close this issue.

2 participants