Mercurial > repos > public > sbplib_julia
changeset 750:f88b2117dc69 feature/laplace_opset
Merge in default
author | Vidar Stiernström <vidar.stiernstrom@it.uu.se> |
---|---|
date | Fri, 19 Mar 2021 16:52:53 +0100 |
parents | c16abc564b82 (diff) 688767a6f3eb (current diff) |
children | f94feb005e7d |
files | Notes.md test/SbpOperators/volumeops/laplace/laplace_test.jl test/runtests.jl test/testSbpOperators.jl |
diffstat | 6 files changed, 231 insertions(+), 12 deletions(-) [+] |
line wrap: on
line diff
--- a/Notes.md Fri Mar 19 16:40:30 2021 +0100 +++ b/Notes.md Fri Mar 19 16:52:53 2021 +0100 @@ -130,6 +130,7 @@ - [ ] How do we handle mixes of periodic and non-periodic grids? Seems it should be supported on the grid level and on the 1d operator level. Between there it should be transparent. - [ ] Can we have a trait to tell if a TensorMapping is transposable? - [ ] Is it ok to have "Constructors" for abstract types which create subtypes? For example a Grids() functions that gives different kind of grids based on input? + - [ ] Figure out how to treat the borrowing parameters of operators. Include in into the struct? Expose via function dispatched on the operator type and grid? ## Regions and tensormappings - [ ] Use a trait to indicate if a TensorMapping uses indices with regions.
--- a/src/Grids/AbstractGrid.jl Fri Mar 19 16:40:30 2021 +0100 +++ b/src/Grids/AbstractGrid.jl Fri Mar 19 16:52:53 2021 +0100 @@ -7,7 +7,7 @@ """ abstract type AbstractGrid end - +export AbstractGrid function dimension end function points end export dimension, points
--- a/src/SbpOperators/volumeops/laplace/laplace.jl Fri Mar 19 16:40:30 2021 +0100 +++ b/src/SbpOperators/volumeops/laplace/laplace.jl Fri Mar 19 16:52:53 2021 +0100 @@ -1,19 +1,87 @@ """ - laplace(grid::EquidistantGrid{Dim}, inner_stencil, closure_stencils) + Laplace{T, Dim, ...} <: TensorMapping{T,Dim,Dim} + Laplace(grid::AbstractGrid, fn; order) + +Implements the Laplace operator, approximating ∑d²/xᵢ² , i = 1,...,`Dim` as a +`TensorMapping`. Additionally, `Laplace` stores the inner product and boundary +operators relevant for constructing a SBP finite difference scheme as `TensorMapping`s. + +Laplace(grid, fn; order) creates the Laplace operator defined on `grid`, +where the operators are read from TOML. The differential operator is created +using `laplace(grid::AbstractGrid,...)`. +""" +struct Laplace{T, Dim, Rb, TMdiffop<:TensorMapping{T,Dim,Dim}, # Differential operator + TMipop<:TensorMapping{T,Dim,Dim}, # Inner product operator + TMbop<:TensorMapping{T,Rb,Dim}, # Boundary operator + TMbqop<:TensorMapping{T,Rb,Rb}, # Boundary quadrature + BID<:BoundaryIdentifier} <: TensorMapping{T,Dim,Dim} + D::TMdiffop # Differential operator + H::TMipop # Inner product operator + H_inv::TMipop # Inverse inner product operator + e::Dict{BID,TMbop} # Boundary restriction operators + d::Dict{BID,TMbop} # Normal derivative operators + H_boundary::Dict{BID,TMbqop} # Boundary quadrature operators +end +export Laplace + +function Laplace(grid::AbstractGrid, fn; order) + # TODO: Removed once we can construct the volume and + # boundary operators by op(grid, fn; order,...). + # Read stencils + op = read_D2_operator(fn; order) + D_inner_stecil = op.innerStencil + D_closure_stencils = op.closureStencils + H_closure_stencils = op.quadratureClosure + e_closure_stencil = op.eClosure + d_closure_stencil = op.dClosure + + # Volume operators + Δ = laplace(grid, D_inner_stecil, D_closure_stencils) + H = inner_product(grid, H_closure_stencils) + H⁻¹ = inverse_inner_product(grid, H_closure_stencils) + + # Boundary operator - id pairs + ids = boundary_identifiers(grid) + n_ids = length(ids) + e_pairs = ntuple(i -> Pair(ids[i],boundary_restriction(grid,e_closure_stencil,ids[i])),n_ids) + d_pairs = ntuple(i -> Pair(ids[i],normal_derivative(grid,d_closure_stencil,ids[i])),n_ids) + Hᵧ_pairs = ntuple(i -> Pair(ids[i],inner_product(boundary_grid(grid,ids[i]),H_closure_stencils)),n_ids) + + return Laplace(Δ, H, H⁻¹, Dict(e_pairs), Dict(d_pairs), Dict(Hᵧ_pairs)) +end + +LazyTensors.range_size(L::Laplace) = LazyTensors.range_size(L.D) +LazyTensors.domain_size(L::Laplace) = LazyTensors.domain_size(L.D) +LazyTensors.apply(L::Laplace, v::AbstractArray, I...) = LazyTensors.apply(L.D,v,I...) + +inner_product(L::Laplace) = L.H +export inner_product +inverse_inner_product(L::Laplace) = L.H_inv +export inverse_inner_product +boundary_restriction(L::Laplace,bid::BoundaryIdentifier) = L.e[bid] +export boundary_restriction +normal_derivative(L::Laplace,bid::BoundaryIdentifier) = L.d[bid] +export normal_derivative +# TODO: boundary_inner_product? +boundary_quadrature(L::Laplace,bid::BoundaryIdentifier) = L.H_boundary[bid] +export boundary_quadrature + +""" + laplace(grid, inner_stencil, closure_stencils) Creates the Laplace operator operator `Δ` as a `TensorMapping` -`Δ` approximates the Laplace operator ∑d²/xᵢ² , i = 1,...,`Dim` on `grid`, using -the stencil `inner_stencil` in the interior and a set of stencils `closure_stencils` -for the points in the closure regions. +`Δ` approximates the Laplace operator ∑d²/xᵢ² , i = 1,...,N on the N-dimensional +`grid`, using the stencil `inner_stencil` in the interior and a set of stencils +`closure_stencils` for the points in the closure regions. On a one-dimensional `grid`, `Δ` is equivalent to `second_derivative`. On a multi-dimensional `grid`, `Δ` is the sum of multi-dimensional `second_derivative`s where the sum is carried out lazily. """ -function laplace(grid::EquidistantGrid{Dim}, inner_stencil, closure_stencils) where Dim +function laplace(grid::AbstractGrid, inner_stencil, closure_stencils) Δ = second_derivative(grid, inner_stencil, closure_stencils, 1) - for d = 2:Dim + for d = 2:dimension(grid) Δ += second_derivative(grid, inner_stencil, closure_stencils, d) end return Δ
--- a/test/SbpOperators/volumeops/laplace/laplace_test.jl Fri Mar 19 16:40:30 2021 +0100 +++ b/test/SbpOperators/volumeops/laplace/laplace_test.jl Fri Mar 19 16:52:53 2021 +0100 @@ -3,12 +3,95 @@ using Sbplib.SbpOperators using Sbplib.Grids using Sbplib.LazyTensors +using Sbplib.RegionIndices + +""" + cmp_fields(s1,s2) + +Compares the fields of two structs s1, s2, using the == operator. +""" +function cmp_fields(s1::T,s2::T) where T + f = fieldnames(T) + return getfield.(Ref(s1),f) == getfield.(Ref(s2),f) +end @testset "Laplace" begin g_1D = EquidistantGrid(101, 0.0, 1.) g_3D = EquidistantGrid((51,101,52), (0.0, -1.0, 0.0), (1., 1., 1.)) + op = read_D2_operator(sbp_operators_path()*"standard_diagonal.toml"; order=4) @testset "Constructors" begin - op = read_D2_operator(sbp_operators_path()*"standard_diagonal.toml"; order=4) + + @testset "1D" begin + # Create all tensor mappings included in Laplace + Δ = laplace(g_1D, op.innerStencil, op.closureStencils) + H = inner_product(g_1D, op.quadratureClosure) + Hi = inverse_inner_product(g_1D, op.quadratureClosure) + + (id_l, id_r) = boundary_identifiers(g_1D) + + e_l = boundary_restriction(g_1D,op.eClosure,id_l) + e_r = boundary_restriction(g_1D,op.eClosure,id_r) + e_dict = Dict(Pair(id_l,e_l),Pair(id_r,e_r)) + + d_l = normal_derivative(g_1D,op.dClosure,id_l) + d_r = normal_derivative(g_1D,op.dClosure,id_r) + d_dict = Dict(Pair(id_l,d_l),Pair(id_r,d_r)) + + H_l = inner_product(boundary_grid(g_1D,id_l),op.quadratureClosure) + H_r = inner_product(boundary_grid(g_1D,id_r),op.quadratureClosure) + Hb_dict = Dict(Pair(id_l,H_l),Pair(id_r,H_r)) + + L = Laplace(g_1D, sbp_operators_path()*"standard_diagonal.toml"; order=4) + @test cmp_fields(L,Laplace(Δ,H,Hi,e_dict,d_dict,Hb_dict)) + @test L isa TensorMapping{T,1,1} where T + @inferred Laplace(Δ,H,Hi,e_dict,d_dict,Hb_dict) + end + @testset "3D" begin + # Create all tensor mappings included in Laplace + Δ = laplace(g_3D, op.innerStencil, op.closureStencils) + H = inner_product(g_3D, op.quadratureClosure) + Hi = inverse_inner_product(g_3D, op.quadratureClosure) + + (id_l, id_r, id_s, id_n, id_b, id_t) = boundary_identifiers(g_3D) + + e_l = boundary_restriction(g_3D,op.eClosure,id_l) + e_r = boundary_restriction(g_3D,op.eClosure,id_r) + e_s = boundary_restriction(g_3D,op.eClosure,id_s) + e_n = boundary_restriction(g_3D,op.eClosure,id_n) + e_b = boundary_restriction(g_3D,op.eClosure,id_b) + e_t = boundary_restriction(g_3D,op.eClosure,id_t) + e_dict = Dict(Pair(id_l,e_l),Pair(id_r,e_r), + Pair(id_s,e_s),Pair(id_n,e_n), + Pair(id_b,e_b),Pair(id_t,e_t)) + + d_l = normal_derivative(g_3D,op.dClosure,id_l) + d_r = normal_derivative(g_3D,op.dClosure,id_r) + d_s = normal_derivative(g_3D,op.dClosure,id_s) + d_n = normal_derivative(g_3D,op.dClosure,id_n) + d_b = normal_derivative(g_3D,op.dClosure,id_b) + d_t = normal_derivative(g_3D,op.dClosure,id_t) + d_dict = Dict(Pair(id_l,d_l),Pair(id_r,d_r), + Pair(id_s,d_s),Pair(id_n,d_n), + Pair(id_b,d_b),Pair(id_t,d_t)) + + H_l = inner_product(boundary_grid(g_3D,id_l),op.quadratureClosure) + H_r = inner_product(boundary_grid(g_3D,id_r),op.quadratureClosure) + H_s = inner_product(boundary_grid(g_3D,id_s),op.quadratureClosure) + H_n = inner_product(boundary_grid(g_3D,id_n),op.quadratureClosure) + H_b = inner_product(boundary_grid(g_3D,id_b),op.quadratureClosure) + H_t = inner_product(boundary_grid(g_3D,id_t),op.quadratureClosure) + Hb_dict = Dict(Pair(id_l,H_l),Pair(id_r,H_r), + Pair(id_s,H_s),Pair(id_n,H_n), + Pair(id_b,H_b),Pair(id_t,H_t)) + + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=4) + @test cmp_fields(L,Laplace(Δ,H,Hi,e_dict,d_dict,Hb_dict)) + @test L isa TensorMapping{T,3,3} where T + @inferred Laplace(Δ,H,Hi,e_dict,d_dict,Hb_dict) + end + end + + @testset "laplace" begin @testset "1D" begin L = laplace(g_1D, op.innerStencil, op.closureStencils) @test L == second_derivative(g_1D, op.innerStencil, op.closureStencils) @@ -21,9 +104,68 @@ Dyy = second_derivative(g_3D, op.innerStencil, op.closureStencils,2) Dzz = second_derivative(g_3D, op.innerStencil, op.closureStencils,3) @test L == Dxx + Dyy + Dzz + @test L isa TensorMapping{T,3,3} where T end end + @testset "inner_product" begin + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=4) + @test inner_product(L) == inner_product(g_3D,op.quadratureClosure) + end + + @testset "inverse_inner_product" begin + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=4) + @test inverse_inner_product(L) == inverse_inner_product(g_3D,op.quadratureClosure) + end + + @testset "boundary_restriction" begin + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=4) + id_l = CartesianBoundary{1,Lower}() + id_r = CartesianBoundary{1,Upper}() + id_s = CartesianBoundary{2,Lower}() + id_n = CartesianBoundary{2,Upper}() + id_b = CartesianBoundary{3,Lower}() + id_t = CartesianBoundary{3,Upper}() + @test boundary_restriction(L,id_l) == boundary_restriction(g_3D,op.eClosure,id_l) + @test boundary_restriction(L,id_r) == boundary_restriction(g_3D,op.eClosure,id_r) + @test boundary_restriction(L,id_s) == boundary_restriction(g_3D,op.eClosure,id_s) + @test boundary_restriction(L,id_n) == boundary_restriction(g_3D,op.eClosure,id_n) + @test boundary_restriction(L,id_b) == boundary_restriction(g_3D,op.eClosure,id_b) + @test boundary_restriction(L,id_t) == boundary_restriction(g_3D,op.eClosure,id_t) + end + + @testset "normal_derivative" begin + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=4) + id_l = CartesianBoundary{1,Lower}() + id_r = CartesianBoundary{1,Upper}() + id_s = CartesianBoundary{2,Lower}() + id_n = CartesianBoundary{2,Upper}() + id_b = CartesianBoundary{3,Lower}() + id_t = CartesianBoundary{3,Upper}() + @test normal_derivative(L,id_l) == normal_derivative(g_3D,op.dClosure,id_l) + @test normal_derivative(L,id_r) == normal_derivative(g_3D,op.dClosure,id_r) + @test normal_derivative(L,id_s) == normal_derivative(g_3D,op.dClosure,id_s) + @test normal_derivative(L,id_n) == normal_derivative(g_3D,op.dClosure,id_n) + @test normal_derivative(L,id_b) == normal_derivative(g_3D,op.dClosure,id_b) + @test normal_derivative(L,id_t) == normal_derivative(g_3D,op.dClosure,id_t) + end + + @testset "boundary_quadrature" begin + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=4) + id_l = CartesianBoundary{1,Lower}() + id_r = CartesianBoundary{1,Upper}() + id_s = CartesianBoundary{2,Lower}() + id_n = CartesianBoundary{2,Upper}() + id_b = CartesianBoundary{3,Lower}() + id_t = CartesianBoundary{3,Upper}() + @test boundary_quadrature(L,id_l) == inner_product(boundary_grid(g_3D,id_l),op.quadratureClosure) + @test boundary_quadrature(L,id_r) == inner_product(boundary_grid(g_3D,id_r),op.quadratureClosure) + @test boundary_quadrature(L,id_s) == inner_product(boundary_grid(g_3D,id_s),op.quadratureClosure) + @test boundary_quadrature(L,id_n) == inner_product(boundary_grid(g_3D,id_n),op.quadratureClosure) + @test boundary_quadrature(L,id_b) == inner_product(boundary_grid(g_3D,id_b),op.quadratureClosure) + @test boundary_quadrature(L,id_t) == inner_product(boundary_grid(g_3D,id_t),op.quadratureClosure) + end + # Exact differentiation is measured point-wise. In other cases # the error is measured in the l2-norm. @testset "Accuracy" begin @@ -40,8 +182,7 @@ # 2nd order interior stencil, 1st order boundary stencil, # implies that L*v should be exact for binomials up to order 2. @testset "2nd order" begin - op = read_D2_operator(sbp_operators_path()*"standard_diagonal.toml"; order=2) - L = laplace(g_3D,op.innerStencil,op.closureStencils) + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=2) @test L*polynomials[1] ≈ zeros(Float64, size(g_3D)...) atol = 5e-9 @test L*polynomials[2] ≈ zeros(Float64, size(g_3D)...) atol = 5e-9 @test L*polynomials[3] ≈ polynomials[1] atol = 5e-9 @@ -51,8 +192,7 @@ # 4th order interior stencil, 2nd order boundary stencil, # implies that L*v should be exact for binomials up to order 3. @testset "4th order" begin - op = read_D2_operator(sbp_operators_path()*"standard_diagonal.toml"; order=4) - L = laplace(g_3D,op.innerStencil,op.closureStencils) + L = Laplace(g_3D, sbp_operators_path()*"standard_diagonal.toml"; order=4) # NOTE: high tolerances for checking the "exact" differentiation # due to accumulation of round-off errors/cancellation errors? @test L*polynomials[1] ≈ zeros(Float64, size(g_3D)...) atol = 5e-9
--- a/test/runtests.jl Fri Mar 19 16:40:30 2021 +0100 +++ b/test/runtests.jl Fri Mar 19 16:52:53 2021 +0100 @@ -1,3 +1,4 @@ +include("test_utils.jl") using Test using Glob
--- /dev/null Thu Jan 01 00:00:00 1970 +0000 +++ b/test/test_utils.jl Fri Mar 19 16:52:53 2021 +0100 @@ -0,0 +1,9 @@ +""" + cmp_fields(s1,s2) + +Compares the fields of two structs s1, s2, using the == operator. +""" +function cmp_fields(s1::T,s2::T) where T + f = fieldnames(T) + return getfield.(Ref(s1),f) == getfield.(Ref(s2),f) +end