Mercurial > repos > public > sbplib_julia
changeset 2055:274f4c1ce4b5 feature/sbp_operators/laplace_curvilinear tip
Add test for SBP property of laplace in mapped grid
| author | Jonatan Werpers <jonatan@werpers.com> |
|---|---|
| date | Sun, 08 Feb 2026 00:05:40 +0100 |
| parents | c08bc343d1cd |
| children | |
| files | test/SbpOperators/volumeops/laplace/laplace_test.jl |
| diffstat | 1 files changed, 33 insertions(+), 0 deletions(-) [+] |
line wrap: on
line diff
--- a/test/SbpOperators/volumeops/laplace/laplace_test.jl Sat Feb 07 22:54:36 2026 +0100 +++ b/test/SbpOperators/volumeops/laplace/laplace_test.jl Sun Feb 08 00:05:40 2026 +0100 @@ -5,6 +5,9 @@ using Diffinitive.LazyTensors using StaticArrays +using SparseArrays +using Tokens +using LinearAlgebra @testset "Laplace" begin # Default stencils (4th order) @@ -107,6 +110,36 @@ @test collect(Δ*gf) isa Array{<:Any,2} @test Δ*gf ≈ map(Δf, g) rtol=2e-2 + + + @testset "SBP property" begin + g = equidistant_grid(c, 20,20) + Δ = laplace(g, stencil_set) + H = inner_product(g, stencil_set) + es = map(boundary_identifiers(g)) do id + boundary_restriction(g, stencil_set, id) + end + ds = map(boundary_identifiers(g)) do id + normal_derivative(g, stencil_set, id) + end + Hᵧs = map(boundary_identifiers(g)) do id + inner_product(boundary_grid(g, id), stencil_set) + end + + BT = mapreduce(+, es, ds,Hᵧs) do e, d, Hᵧ + e'∘Hᵧ∘d + end + M = -H∘Δ + BT + + M = sparse(M) + @test M ≈ M' + + function issemiposdef(A, tol=1e-8) + return isposdef(A+tol*I) + end + + @test issemiposdef(Symmetric(M)) + end end end
