changeset 1660:6d196fb85133 feature/sbp_operators/laplace_curvilinear

Merge feature/grids/curvilinear
author Jonatan Werpers <jonatan@werpers.com>
date Fri, 28 Jun 2024 17:04:05 +0200
parents cc9d18a5ff2d (diff) 3bbcd496e021 (current diff)
children bdb4becac704
files src/Grids/Grids.jl src/Grids/mapped_grid.jl test/Grids/mapped_grid_test.jl
diffstat 21 files changed, 731 insertions(+), 48 deletions(-) [+]
line wrap: on
line diff
--- a/Manifest.toml	Fri Jun 28 17:02:47 2024 +0200
+++ b/Manifest.toml	Fri Jun 28 17:04:05 2024 +0200
@@ -2,7 +2,7 @@
 
 julia_version = "1.10.4"
 manifest_format = "2.0"
-project_hash = "a36735c53cfa4453f39635046eeaa47a4ea1231b"
+project_hash = "32fac879810099260f177c27318d3f26de4a00cc"
 
 [[deps.Adapt]]
 deps = ["LinearAlgebra", "Requires"]
--- a/Project.toml	Fri Jun 28 17:02:47 2024 +0200
+++ b/Project.toml	Fri Jun 28 17:04:05 2024 +0200
@@ -4,9 +4,18 @@
 version = "0.1.1"
 
 [deps]
+LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
 StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"
 TOML = "fa267f1f-6049-4f14-aa54-33bafae1ed76"
 TiledIteration = "06e1c1a7-607b-532d-9fad-de7d9aa2abac"
 
+[weakdeps]
+Makie = "ee78f7c6-11fb-53f2-987a-cfe4a2b5a57a"
+Plots = "91a5bcdd-55d7-5caf-9e0b-520d859cae80"
+
+[extensions]
+SbplibMakieExt = "Makie"
+SbplibPlotsExt = "Plots"
+
 [compat]
 julia = "1.5"
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/ext/SbplibMakieExt.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -0,0 +1,91 @@
+module SbplibMakieExt
+
+using Sbplib.Grids
+using Makie
+using StaticArrays
+
+
+function verticies_and_faces_and_values(g::Grid{<:Any,2}, gf::AbstractArray{<:Any, 2})
+    ps = map(Tuple, g)[:]
+    values = gf[:]
+
+    N = length(ps)
+
+    faces = Vector{NTuple{3,Int}}()
+
+    n,m = size(g)
+    Li = LinearIndices((1:n, 1:m))
+    for i ∈ 1:n-1, j = 1:m-1
+
+        # Add point in the middle of the patch to preserve symmetries
+        push!(ps, Tuple((g[i,j] + g[i+1,j] + g[i+1,j+1] + g[i,j+1])/4))
+        push!(values, (gf[i,j] + gf[i+1,j] + gf[i+1,j+1] + gf[i,j+1])/4)
+
+        push!(faces, (Li[i,j],     Li[i+1,j],   length(ps)))
+        push!(faces, (Li[i+1,j],   Li[i+1,j+1], length(ps)))
+        push!(faces, (Li[i+1,j+1], Li[i,j+1],   length(ps)))
+        push!(faces, (Li[i,j+1],   Li[i,j],     length(ps)))
+    end
+
+    verticies = permutedims(reinterpret(reshape,eltype(eltype(ps)), ps))
+    faces = permutedims(reinterpret(reshape,Int, faces))
+
+    return verticies, faces, values
+end
+
+
+## Grids
+
+Makie.convert_arguments(::Type{<:Scatter}, g::Grid) = (reshape(map(Point,g),:),) # (map(Point,collect(g)[:]),)
+function Makie.convert_arguments(::Type{<:Lines}, g::Grid{<:Any,2})
+    M = collect(g)
+
+    function cat_with_NaN(a,b)
+        vcat(a,[@SVector[NaN,NaN]],b)
+    end
+
+    xlines = reduce(cat_with_NaN, eachrow(M))
+    ylines = reduce(cat_with_NaN, eachcol(M))
+
+    return (cat_with_NaN(xlines,ylines),)
+end
+
+Makie.plot!(plot::Plot(Grid{<:Any,2})) = lines!(plot, plot.attributes, plot[1])
+
+
+## Grid functions
+
+### 1D
+function Makie.convert_arguments(::Type{<:Lines}, g::Grid{<:Any,1}, gf::AbstractArray{<:Any, 1})
+    (collect(g), gf)
+end
+
+function Makie.convert_arguments(::Type{<:Scatter}, g::Grid{<:Any,1}, gf::AbstractArray{<:Any, 1})
+    (collect(g), gf)
+end
+
+Makie.plot!(plot::Plot(Grid{<:Any,1}, AbstractArray{<:Any,1})) = lines!(plot, plot.attributes, plot[1], plot[2])
+
+### 2D
+function Makie.convert_arguments(::Type{<:Surface}, g::Grid{<:Any,2}, gf::AbstractArray{<:Any, 2})
+    (getindex.(g,1), getindex.(g,2), gf)
+end
+
+function Makie.plot!(plot::Plot(Grid{<:Any,2},AbstractArray{<:Any, 2}))
+    r = @lift verticies_and_faces_and_values($(plot[1]), $(plot[2]))
+    v,f,c = (@lift $r[1]), (@lift $r[2]), (@lift $r[3])
+    mesh!(plot, plot.attributes, v, f;
+        color=c,
+        shading = NoShading,
+    )
+end
+# TBD: Can we define `mesh` instead of the above function and then forward plot! to that?
+
+function Makie.convert_arguments(::Type{<:Scatter}, g::Grid{<:Any,2}, gf::AbstractArray{<:Any, 2})
+    ps = map(g,gf) do (x,y), z
+        @SVector[x,y,z]
+    end
+    (reshape(ps,:),)
+end
+
+end
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/ext/SbplibPlotsExt.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -0,0 +1,72 @@
+module SbplibPlotsExt
+
+using Sbplib.Grids
+using Plots
+
+@recipe f(::Type{<:Grid}, g::Grid) = map(Tuple,g)[:]
+
+@recipe function f(c::Chart{2,<:Rectangle}, n=5, m=n; draw_border=true, bordercolor=1)
+    Ξ = parameterspace(c)
+    ξs = range(limits(Ξ,1)..., n)
+    ηs = range(limits(Ξ,2)..., m)
+
+    label := false
+    seriescolor --> 2
+    for ξ ∈ ξs
+        @series adapted_curve_grid(η->c((ξ,η)),limits(Ξ,1))
+    end
+
+    for η ∈ ηs
+        @series adapted_curve_grid(ξ->c((ξ,η)),limits(Ξ,2))
+    end
+
+    if ~draw_border
+        return
+    end
+
+    for ξ ∈ limits(Ξ,1)
+        @series begin
+            linewidth --> 3
+            seriescolor := bordercolor
+            adapted_curve_grid(η->c((ξ,η)),limits(Ξ,1))
+        end
+    end
+
+    for η ∈ limits(Ξ,2)
+        @series begin
+            linewidth --> 3
+            seriescolor := bordercolor
+            adapted_curve_grid(ξ->c((ξ,η)),limits(Ξ,2))
+        end
+    end
+end
+
+function adapted_curve_grid(g, minmax)
+    t1, _ = PlotUtils.adapted_grid(t->g(t)[1], minmax)
+    t2, _ = PlotUtils.adapted_grid(t->g(t)[2], minmax)
+
+    ts = sort(vcat(t1,t2))
+
+    x = map(ts) do t
+        g(t)[1]
+    end
+    y = map(ts) do t
+        g(t)[2]
+    end
+
+    return x, y
+end
+
+# get_axis_limits(plt, :x)
+
+
+# ReicpesPipline/src/user_recipe.jl
+# @recipe function f(f::FuncOrFuncs{F}) where {F<:Function}
+
+# @recipe function f(f::Function, xmin::Number, xmax::Number)
+
+# _scaled_adapted_grid(f, xscale, yscale, xmin, xmax)
+
+end
+
+
--- a/src/Grids/Grids.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/Grids/Grids.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -4,6 +4,31 @@
 using Sbplib.RegionIndices
 using Sbplib.LazyTensors
 using StaticArrays
+using LinearAlgebra
+
+export ParameterSpace
+export HyperBox
+export Simplex
+export Interval
+export Rectangle
+export Box
+export Triangle
+export Tetrahedron
+
+export limits
+export unitinterval
+export unitsquare
+export unitcube
+export unithyperbox
+
+export verticies
+export unittriangle
+export unittetrahedron
+export unitsimplex
+
+export Chart
+export ConcreteChart
+export parameterspace
 
 # Grid
 export Grid
@@ -45,6 +70,7 @@
 
 abstract type BoundaryIdentifier end
 
+include("manifolds.jl")
 include("grid.jl")
 include("tensor_grid.jl")
 include("equidistant_grid.jl")
--- a/src/Grids/equidistant_grid.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/Grids/equidistant_grid.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -130,6 +130,20 @@
 	return EquidistantGrid(range(limit_lower, limit_upper, length=size)) # TBD: Should it use LinRange instead?
 end
 
+
+equidistant_grid(hb::HyperBox, dims::Vararg{Int}) = equidistant_grid(hb.a, hb.b, dims...)
+# TODO: One dimensional grids shouldn't have vector eltype right?, Change here or in HyperBox?
+
+function equidistant_grid(c::Chart, dims::Vararg{Int})
+    lg = equidistant_grid(parameterspace(c), dims...)
+    return MappedGrid(
+        lg,
+        map(c,lg),
+        map(ξ->jacobian(c, ξ), lg),
+    )
+end
+
+
 CartesianBoundary{D,BID} = TensorGridBoundary{D,BID} # TBD: What should we do about the naming of this boundary?
 
 
--- a/src/Grids/grid.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/Grids/grid.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -107,6 +107,8 @@
 """
 function boundary_identifiers end
 
+# TBD: Boundary identifiers for charts and atlases?
+
 """
     boundary_grid(g::Grid, id::BoundaryIdentifier)
 
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/src/Grids/manifolds.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -0,0 +1,207 @@
+"""
+    ParameterSpace{D}
+
+A space of parameters of dimension `D`. Used with `Chart` to indicate which
+parameters are valid for that chart.
+
+Common parameter spaces are created using the functions unit sized spaces
+* `unitinterval`
+* `unitrectangle`
+* `unitbox`
+* `unittriangle`
+* `unittetrahedron`
+* `unithyperbox`
+* `unitsimplex`
+
+See also: [`Interval`](@ref), [`Rectangle`](@ref), [`Box`](@ref),
+[`Triangle`](@ref), [`Tetrahedron`](@ref), [`HyperBox`](@ref),
+[`Simplex`](@ref),
+"""
+abstract type ParameterSpace{D} end
+Base.ndims(::ParameterSpace{D}) where D = D
+# TBD:  Should implement domain_dim?
+
+struct HyperBox{T,D} <: ParameterSpace{D}
+    a::SVector{D,T}
+    b::SVector{D,T}
+end
+
+function HyperBox(a,b)
+    T = SVector{length(a)}
+    HyperBox(convert(T,a), convert(T,b))
+end
+
+Interval{T} = HyperBox{T,1}
+Rectangle{T} = HyperBox{T,2}
+Box{T} = HyperBox{T,3}
+
+limits(box::HyperBox, d) = (box.a[d], box.b[d])
+limits(box::HyperBox) = (box.a, box.b)
+
+unitinterval(T=Float64) = unithyperbox(T,1)
+unitsquare(T=Float64) = unithyperbox(T,2)
+unitcube(T=Float64) = unithyperbox(T,3)
+unithyperbox(T, D) = HyperBox((@SVector zeros(T,D)), (@SVector ones(T,D)))
+unithyperbox(D) = unithyperbox(Float64,D)
+
+
+struct Simplex{T,D,NV} <: ParameterSpace{D}
+    verticies::NTuple{NV,SVector{D,T}}
+end
+
+Simplex(verticies::Vararg{AbstractArray}) = Simplex(Tuple(SVector(v...) for v ∈ verticies))
+
+verticies(s::Simplex) = s.verticies
+
+Triangle{T} = Simplex{T,2}
+Tetrahedron{T} = Simplex{T,3}
+
+unittriangle(T=Float64) = unitsimplex(T,2)
+unittetrahedron(T=Float64) = unitsimplex(T,3)
+function unitsimplex(T,D)
+    z = @SVector zeros(T,D)
+    unitelement = one(eltype(z))
+    verticies = ntuple(i->setindex(z, unitelement, i), D)
+    return Simplex((z,verticies...))
+end
+unitsimplex(D) = unitsimplex(Float64, D)
+
+"""
+    Chart{D}
+
+A parametrized description of a manifold or part of a manifold.
+"""
+struct Chart{D, PST<:ParameterSpace{D}, MT}
+    mapping::MT
+    parameterspace::PST
+end
+
+domain_dim(::Chart{D}) where D = D
+(c::Chart)(ξ) = c.mapping(ξ)
+parameterspace(c::Chart) = c.parameterspace
+
+"""
+    jacobian(c::Chart, ξ)
+
+The jacobian of the mapping evaluated at `ξ`. This defers to the
+implementation of `jacobian` for the mapping itself. If no implementation is
+available one can easily be specified for either the mapping function or the
+chart itself.
+```julia
+c = Chart(f, ps)
+jacobian(f::typeof(f), ξ) = f′(ξ)
+```
+or
+```julia
+c = Chart(f, ps)
+jacobian(c::typeof(c),ξ) = f′(ξ)
+```
+which will both allow calling `jacobian(c,ξ)`.
+"""
+jacobian(c::Chart, ξ) = jacobian(c.mapping, ξ)
+# TBD: Can we register a error hint for when jacobian is called with a function that doesn't have a registered jacobian?
+
+
+# TBD: Should Charts, parameterspaces have boundary names?
+
+"""
+    Atlas
+
+A collection of charts and their connections.
+Should implement methods for `charts` and
+"""
+abstract type Atlas end
+
+"""
+    charts(::Atlas)
+
+The colloction of charts in the atlas.
+"""
+function charts end
+
+"""
+    connections
+
+TBD: What exactly should this return?
+
+"""
+
+struct CartesianAtlas <: Atlas
+    charts::Matrix{Chart}
+end
+
+charts(a::CartesianAtlas) = a.charts
+
+struct UnstructuredAtlas <: Atlas
+    charts::Vector{Chart}
+    connections
+end
+
+charts(a::UnstructuredAtlas) = a.charts
+
+
+###
+# Geometry
+###
+
+abstract type Curve end
+abstract type Surface end
+
+
+struct Line{PT} <: Curve
+    p::PT
+    tangent::PT
+end
+
+(c::Line)(s) = c.p + s*c.tangent
+
+
+struct LineSegment{PT} <: Curve
+    a::PT
+    b::PT
+end
+
+(c::LineSegment)(s) = (1-s)*c.a + s*c.b
+
+
+function linesegments(ps...)
+    return [LineSegment(ps[i], ps[i+1]) for i ∈ 1:length(ps)-1]
+end
+
+
+function polygon_edges(ps...)
+    n = length(ps)
+    return [LineSegment(ps[i], ps[mod1(i+1,n)]) for i ∈ eachindex(Ps)]
+end
+
+struct Circle{T,PT} <: Curve
+    c::PT
+    r::T
+end
+
+(c::Circle)(θ) = c.c + r*@SVector[cos(Θ), sin(Θ)]
+
+struct TransfiniteInterpolationSurface{T1,T2,T3,T4} <: Surface
+    c₁::T1
+    c₂::T2
+    c₃::T3
+    c₄::T4
+end
+
+function (s::TransfiniteInterpolationSurface)(u,v)
+    c₁, c₂, c₃, c₄ = s.c₁, s.c₂, s.c₃, s.c₄
+    P₀₀ = c₁(0)
+    P₁₀ = c₂(0)
+    P₁₁ = c₃(0)
+    P₀₁ = c₄(0)
+    return (1-v)*c₁(u) + u*c₂(v) + v*c₃(1-u) + (1-u)*c₄(1-v) - (
+        (1-u)*(1-v)*P₀₀ + u*(1-v)*P₁₀ + u*v*P₁₁ + (1-u)*v*P₀₁
+    )
+end
+
+function (s::TransfiniteInterpolationSurface)(ξ̄::AbstractArray)
+    s(ξ̄...)
+end
+
+# TODO: Implement jacobian() for the different mapping helpers
+
--- a/src/Grids/mapped_grid.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/Grids/mapped_grid.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -62,6 +62,7 @@
     )
 end
 
+
 function jacobian_determinant(g::MappedGrid)
     return map(jacobian(g)) do ∂x∂ξ
         det(∂x∂ξ)
--- a/src/SbpOperators/boundaryops/boundary_restriction.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/SbpOperators/boundaryops/boundary_restriction.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -25,3 +25,7 @@
     converted_stencil = convert(Stencil{eltype(g)}, closure_stencil)
     return BoundaryOperator(g, converted_stencil, boundary)
 end
+
+function boundary_restriction(g::MappedGrid, stencil_set::StencilSet, boundary)
+    return boundary_restriction(logicalgrid(g), stencil_set, boundary)
+end
--- a/src/SbpOperators/boundaryops/normal_derivative.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/SbpOperators/boundaryops/normal_derivative.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -28,3 +28,30 @@
     scaled_stencil = scale(closure_stencil,h_inv)
     return BoundaryOperator(g, scaled_stencil, boundary)
 end
+
+function normal_derivative(g::MappedGrid, stencil_set::StencilSet, boundary)
+    k = grid_id(boundary)
+    b_indices = boundary_indices(g, boundary)
+
+    # Compute the weights for the logical derivatives
+    g⁻¹ = geometric_tensor_inverse(g)
+    α = map(CartesianIndices(g⁻¹)[b_indices...]) do I # TODO: Fix iterator here
+        gᵏⁱ = g⁻¹[I][k,:]
+        gᵏᵏ = g⁻¹[I][k,k]
+
+        gᵏⁱ./sqrt(gᵏᵏ)
+    end
+
+    # Assemble difference operator
+    mapreduce(+,1:ndims(g)) do i
+        if i == k
+            ∂_ξᵢ = normal_derivative(logicalgrid(g), stencil_set, boundary)
+        else
+            e = boundary_restriction(logicalgrid(g), stencil_set, boundary)
+            ∂_ξᵢ = e ∘ first_derivative(logicalgrid(g), stencil_set, i)
+        end
+
+        αᵢ = componentview(α,i)
+        DiagonalTensor(αᵢ) ∘ ∂_ξᵢ
+    end
+end
--- a/src/SbpOperators/volumeops/inner_products/inner_product.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/SbpOperators/volumeops/inner_products/inner_product.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -50,3 +50,8 @@
 """
 inner_product(g::ZeroDimGrid, stencil_set::StencilSet) = IdentityTensor{component_type(g)}()
 
+
+function inner_product(g::MappedGrid, stencil_set)
+    J = jacobian_determinant(g)
+    DiagonalTensor(J)∘inner_product(logicalgrid(g), stencil_set)
+end
--- a/src/SbpOperators/volumeops/inner_products/inverse_inner_product.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/SbpOperators/volumeops/inner_products/inverse_inner_product.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -49,3 +49,8 @@
 Implemented to simplify 1D code for SBP operators.
 """
 inverse_inner_product(g::ZeroDimGrid, stencil_set::StencilSet) = IdentityTensor{component_type(g)}()
+
+function inverse_inner_product(g::MappedGrid, stencil_set)
+    J⁻¹ = map(inv, jacobian_determinant(g))
+    DiagonalTensor(J⁻¹)∘inner_product(logicalgrid(g), stencil_set)
+end
--- a/src/SbpOperators/volumeops/laplace/laplace.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/src/SbpOperators/volumeops/laplace/laplace.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -51,8 +51,31 @@
     end
     return Δ
 end
+
 laplace(g::EquidistantGrid, stencil_set) = second_derivative(g, stencil_set)
 
+function laplace(grid::MappedGrid, stencil_set)
+    J = jacobian_determinant(grid)
+    J⁻¹ = DiagonalTensor(map(inv, J))
+
+    Jg = map(*, J, geometric_tensor_inverse(grid))
+    lg = logicalgrid(grid)
+
+    return mapreduce(+, CartesianIndices(first(Jg))) do I
+        i, j = I[1], I[2]
+        Jgⁱʲ = componentview(Jg, i, j)
+
+        if i == j
+            J⁻¹∘second_derivative_variable(lg, Jgⁱʲ, stencil_set, i)
+        else
+            Dᵢ = first_derivative(lg, stencil_set, i)
+            Dⱼ = first_derivative(lg, stencil_set, j)
+            J⁻¹∘Dᵢ∘DiagonalTensor(Jgⁱʲ)∘Dⱼ
+        end
+    end
+end
+
+
 """
     sat_tensors(Δ::Laplace, g::Grid, bc::DirichletCondition; H_tuning, R_tuning)
 
--- a/test/Grids/equidistant_grid_test.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/test/Grids/equidistant_grid_test.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -2,6 +2,7 @@
 using Test
 using Sbplib.RegionIndices
 using Sbplib.LazyTensors
+using StaticArrays
 
 
 @testset "EquidistantGrid" begin
@@ -145,6 +146,23 @@
             @test [gp[i]...] ≈ [p[i]...] atol=5e-13
         end
     end
+
+
+    @testset "equidistant_grid(::ParameterSpace)" begin
+        ps = HyperBox((0,0),(2,1))
+
+        @test equidistant_grid(ps, 3,4) == equidistant_grid((0,0), (2,1), 3,4)
+    end
+
+
+    @testset "equidistant_grid(::Chart)" begin
+        c = Chart(unitsquare()) do (ξ,η)
+            @SVector[2ξ, 3η]
+        end
+        Grids.jacobian(c::typeof(c), ξ̄) = @SMatrix[2 0; 0 3]
+
+        @test equidistant_grid(c, 5, 4) isa Grid
+    end
 end
 
 
--- /dev/null	Thu Jan 01 00:00:00 1970 +0000
+++ b/test/Grids/manifolds_test.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -0,0 +1,63 @@
+using Test
+
+using Sbplib.Grids
+using Sbplib.RegionIndices
+using Sbplib.LazyTensors
+
+# using StaticArrays
+
+@testset "ParameterSpace" begin
+    @test ndims(HyperBox([1,1], [2,2])) == 2
+    @test ndims(unittetrahedron()) == 3
+end
+
+@testset "HyperBox" begin
+    @test HyperBox{<:Any, 2} <: ParameterSpace{2}
+    @test HyperBox([1,1], [2,2]) isa HyperBox{Int, 2}
+
+    @test limits(HyperBox([1,2], [3,4])) == ([1,2], [3,4])
+    @test limits(HyperBox([1,2], [3,4]), 1) == (1,3)
+    @test limits(HyperBox([1,2], [3,4]), 2) == (2,4)
+
+    @test unitinterval() isa HyperBox{Float64,1}
+    @test limits(unitinterval()) == ([0], [1])
+
+    @test unitinterval(Int) isa HyperBox{Int,1}
+    @test limits(unitinterval(Int)) == ([0], [1])
+
+    @test unitsquare() isa HyperBox{Float64,2}
+    @test limits(unitsquare()) == ([0,0],[1,1])
+
+    @test unitcube() isa HyperBox{Float64,3}
+    @test limits(unitcube()) == ([0,0,0],[1,1,1])
+
+    @test unithyperbox(4) isa HyperBox{Float64,4}
+    @test limits(unithyperbox(4)) == ([0,0,0,0],[1,1,1,1])
+end
+
+@testset "Simplex" begin
+    @test Simplex{<:Any, 3} <: ParameterSpace{3}
+    @test Simplex([1,2], [3,4]) isa Simplex{Int, 2}
+    @test Simplex([1,2,3], [4,5,6],[1,1,1]) isa Simplex{Int, 3}
+
+    @test verticies(Simplex([1,2], [3,4])) == ([1,2], [3,4])
+
+    @test unittriangle() isa Simplex{Float64,2}
+    @test verticies(unittriangle()) == ([0,0], [1,0], [0,1])
+
+    @test unittetrahedron() isa  Simplex{Float64,3}
+    @test verticies(unittetrahedron()) == ([0,0,0], [1,0,0], [0,1,0],[0,0,1])
+
+    @test unitsimplex(4) isa Simplex{Float64,4}
+end
+
+@testset "Chart" begin
+    c = Chart(x->2x, unitsquare())
+    @test c isa Chart{2}
+    @test c([3,2]) == [6,4]
+    @test parameterspace(c) == unitsquare()
+end
+
+@testset "Atlas" begin
+
+end
--- a/test/SbpOperators/boundaryops/boundary_restriction_test.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/test/SbpOperators/boundaryops/boundary_restriction_test.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -6,6 +6,8 @@
 using Sbplib.RegionIndices
 using Sbplib.SbpOperators: BoundaryOperator, Stencil
 
+using StaticArrays
+
 @testset "boundary_restriction" begin
 	stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order = 4)
 	e_closure = parse_stencil(stencil_set["e"]["closure"])
@@ -33,7 +35,7 @@
     end
 
     @testset "Application" begin
-        @testset "1D" begin
+        @testset "EquidistantGrid" begin
             e_l, e_r = boundary_restriction.(Ref(g_1D), Ref(stencil_set), boundary_identifiers(g_1D))
             v = eval_on(g_1D,x->1+x^2)
             u = fill(3.124)
@@ -43,7 +45,7 @@
             @test (e_r*v)[1] == v[end]
         end
 
-        @testset "2D" begin
+        @testset "TensorGrid" begin
             e_w, e_e, e_s, e_n = boundary_restriction.(Ref(g_2D), Ref(stencil_set), boundary_identifiers(g_2D))
             v = rand(11, 15)
             u = fill(3.124)
@@ -53,5 +55,22 @@
             @test e_s*v == v[:,1]
             @test e_n*v == v[:,end]
        end
+
+       @testset "MappedGrid" begin
+            c = Chart(unitsquare()) do (ξ,η)
+                @SVector[2ξ + η*(1-η), 3η+(1+η/2)*ξ^2]
+            end
+            Grids.jacobian(c::typeof(c), (ξ,η)) = @SMatrix[2 1-2η; (2+η)*ξ 3+ξ^2/2]
+
+            mg = equidistant_grid(c, 10,13)
+
+            e_w, e_e, e_s, e_n = boundary_restriction.(Ref(mg), Ref(stencil_set), boundary_identifiers(mg))
+            v = rand(10, 13)
+
+            @test e_w*v == v[1,:]
+            @test e_e*v == v[end,:]
+            @test e_s*v == v[:,1]
+            @test e_n*v == v[:,end]
+       end
     end
 end
--- a/test/SbpOperators/boundaryops/normal_derivative_test.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/test/SbpOperators/boundaryops/normal_derivative_test.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -6,54 +6,99 @@
 using Sbplib.RegionIndices
 import Sbplib.SbpOperators.BoundaryOperator
 
+using StaticArrays
+
 @testset "normal_derivative" begin
-    g_1D = equidistant_grid(0.0, 1.0, 11)
-    g_2D = equidistant_grid((0.0, 0.0), (1.0,1.0), 11, 12)
-    @testset "normal_derivative" begin
-    	stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=4)
-        @testset "1D" begin
-            d_l = normal_derivative(g_1D, stencil_set, Lower())
-            @test d_l == normal_derivative(g_1D, stencil_set, Lower())
-            @test d_l isa BoundaryOperator{T,Lower} where T
-            @test d_l isa LazyTensor{T,0,1} where T
-        end
-        @testset "2D" begin
-            d_w = normal_derivative(g_2D, stencil_set, CartesianBoundary{1,Lower}())
-            d_n = normal_derivative(g_2D, stencil_set, CartesianBoundary{2,Upper}())
-            Ix = IdentityTensor{Float64}((size(g_2D)[1],))
-            Iy = IdentityTensor{Float64}((size(g_2D)[2],))
-            d_l = normal_derivative(g_2D.grids[1], stencil_set, Lower())
-            d_r = normal_derivative(g_2D.grids[2], stencil_set, Upper())
-            @test d_w == normal_derivative(g_2D, stencil_set, CartesianBoundary{1,Lower}())
-            @test d_w ==  d_l⊗Iy
-            @test d_n ==  Ix⊗d_r
-            @test d_w isa LazyTensor{T,1,2} where T
-            @test d_n isa LazyTensor{T,1,2} where T
+	stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=4)
+
+    @testset "EquidistantGrid" begin
+        g_1D = equidistant_grid(0.0, 1.0, 11)
+
+        d_l = normal_derivative(g_1D, stencil_set, Lower())
+        @test d_l == normal_derivative(g_1D, stencil_set, Lower())
+        @test d_l isa BoundaryOperator{T,Lower} where T
+        @test d_l isa LazyTensor{T,0,1} where T
+    end
+
+    @testset "TensorGrid" begin
+        g_2D = equidistant_grid((0.0, 0.0), (1.0,1.0), 11, 12)
+        d_w = normal_derivative(g_2D, stencil_set, CartesianBoundary{1,Lower}())
+        d_n = normal_derivative(g_2D, stencil_set, CartesianBoundary{2,Upper}())
+        Ix = IdentityTensor{Float64}((size(g_2D)[1],))
+        Iy = IdentityTensor{Float64}((size(g_2D)[2],))
+        d_l = normal_derivative(g_2D.grids[1], stencil_set, Lower())
+        d_r = normal_derivative(g_2D.grids[2], stencil_set, Upper())
+        @test d_w == normal_derivative(g_2D, stencil_set, CartesianBoundary{1,Lower}())
+        @test d_w ==  d_l⊗Iy
+        @test d_n ==  Ix⊗d_r
+        @test d_w isa LazyTensor{T,1,2} where T
+        @test d_n isa LazyTensor{T,1,2} where T
+
+        @testset "Accuracy" begin
+            v = eval_on(g_2D, (x,y)-> x^2 + (y-1)^2 + x*y)
+            v∂x = eval_on(g_2D, (x,y)-> 2*x + y)
+            v∂y = eval_on(g_2D, (x,y)-> 2*(y-1) + x)
+            # TODO: Test for higher order polynomials?
+            @testset "2nd order" begin
+            	stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=2)
+                d_w, d_e, d_s, d_n = normal_derivative.(Ref(g_2D), Ref(stencil_set), boundary_identifiers(g_2D))
+
+                @test d_w*v ≈ -v∂x[1,:] atol = 1e-13
+                @test d_e*v ≈ v∂x[end,:] atol = 1e-13
+                @test d_s*v ≈ -v∂y[:,1] atol = 1e-13
+                @test d_n*v ≈ v∂y[:,end] atol = 1e-13
+            end
+
+            @testset "4th order" begin
+                stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=4)
+                d_w, d_e, d_s, d_n = normal_derivative.(Ref(g_2D), Ref(stencil_set), boundary_identifiers(g_2D))
+
+                @test d_w*v ≈ -v∂x[1,:] atol = 1e-13
+                @test d_e*v ≈ v∂x[end,:] atol = 1e-13
+                @test d_s*v ≈ -v∂y[:,1] atol = 1e-13
+                @test d_n*v ≈ v∂y[:,end] atol = 1e-13
+            end
         end
     end
-    @testset "Accuracy" begin
-        v = eval_on(g_2D, (x,y)-> x^2 + (y-1)^2 + x*y)
-        v∂x = eval_on(g_2D, (x,y)-> 2*x + y)
-        v∂y = eval_on(g_2D, (x,y)-> 2*(y-1) + x)
-        # TODO: Test for higher order polynomials?
-        @testset "2nd order" begin
-        	stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=2)
-            d_w, d_e, d_s, d_n = normal_derivative.(Ref(g_2D), Ref(stencil_set), boundary_identifiers(g_2D))
+
+    @testset "MappedGrid" begin
+        c = Chart(unitsquare()) do (ξ,η)
+            @SVector[2ξ + η*(1-η), 3η+(1+η/2)*ξ^2]
+        end
+        Grids.jacobian(c::typeof(c), (ξ,η)) = @SMatrix[2 1-2η; (2+η)*ξ 3+ξ^2/2]
+
+        mg = equidistant_grid(c, 10,13)
+
+        b_w, b_e, b_s, b_n = boundary_identifiers(mg)
+
+        @test_broken normal_derivative(mg, stencil_set, b_w) isa LazyTensor{<:Any, 1, 2}
+        @test_broken normal_derivative(mg, stencil_set, b_e) isa LazyTensor{<:Any, 1, 2}
+        @test_broken normal_derivative(mg, stencil_set, b_s) isa LazyTensor{<:Any, 1, 2}
+        @test_broken normal_derivative(mg, stencil_set, b_n) isa LazyTensor{<:Any, 1, 2}
+
+        @testset "Accuracy" begin
+            v = map(x̄ -> NaN, mg)
+            dₙv = map(x̄ -> NaN, mg)
 
-            @test d_w*v ≈ -v∂x[1,:] atol = 1e-13
-            @test d_e*v ≈ v∂x[end,:] atol = 1e-13
-            @test d_s*v ≈ -v∂y[:,1] atol = 1e-13
-            @test d_n*v ≈ v∂y[:,end] atol = 1e-13
-        end
+            @testset "2nd order" begin
+                stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=2)
+                d_w, d_e, d_s, d_n = normal_derivative.(Ref(mg), Ref(stencil_set), boundary_identifiers(mg))
+
+                @test_broken d_w*v ≈ dₙv atol = 1e-13
+                @test_broken d_e*v ≈ dₙv atol = 1e-13
+                @test_broken d_s*v ≈ dₙv atol = 1e-13
+                @test_broken d_n*v ≈ dₙv atol = 1e-13
+            end
 
-        @testset "4th order" begin
-            stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=4)
-            d_w, d_e, d_s, d_n = normal_derivative.(Ref(g_2D), Ref(stencil_set), boundary_identifiers(g_2D))
-            
-            @test d_w*v ≈ -v∂x[1,:] atol = 1e-13
-            @test d_e*v ≈ v∂x[end,:] atol = 1e-13
-            @test d_s*v ≈ -v∂y[:,1] atol = 1e-13
-            @test d_n*v ≈ v∂y[:,end] atol = 1e-13
+            @testset "4th order" begin
+                stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=4)
+                d_w, d_e, d_s, d_n = normal_derivative.(Ref(mg), Ref(stencil_set), boundary_identifiers(mg))
+
+                @test_broken d_w*v ≈ dₙv atol = 1e-13
+                @test_broken d_e*v ≈ dₙv atol = 1e-13
+                @test_broken d_s*v ≈ dₙv atol = 1e-13
+                @test_broken d_n*v ≈ dₙv atol = 1e-13
+            end
         end
     end
 end
--- a/test/SbpOperators/volumeops/inner_products/inner_product_test.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/test/SbpOperators/volumeops/inner_products/inner_product_test.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -6,6 +6,8 @@
 
 import Sbplib.SbpOperators.ConstantInteriorScalingOperator
 
+using StaticArrays
+
 @testset "Diagonal-stencil inner_product" begin
     Lx = π/2.
     Ly = Float64(π)
@@ -94,4 +96,17 @@
             end
         end
     end
+
+    @testset "MappedGrid" begin
+        stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=4)
+        c = Chart(unitsquare()) do (ξ,η)
+            @SVector[2ξ + η*(1-η), 3η+(1+η/2)*ξ^2]
+        end
+        Grids.jacobian(c::typeof(c), (ξ,η)) = @SMatrix[2 1-2η; (2+η)*ξ 3+ξ^2/2]
+
+        mg = equidistant_grid(c, 10,13)
+
+        @test inner_product(mg, stencil_set) isa LazyTensor{<:Any, 2,2}
+        @test_broken false # Test that it calculates the right thing
+    end
 end
--- a/test/SbpOperators/volumeops/inner_products/inverse_inner_product_test.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/test/SbpOperators/volumeops/inner_products/inverse_inner_product_test.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -6,6 +6,8 @@
 
 import Sbplib.SbpOperators.ConstantInteriorScalingOperator
 
+using StaticArrays
+
 @testset "Diagonal-stencil inverse_inner_product" begin
     Lx = π/2.
     Ly = Float64(π)
@@ -82,4 +84,17 @@
             end
         end
     end
+
+    @testset "MappedGrid" begin
+        stencil_set = read_stencil_set(sbp_operators_path()*"standard_diagonal.toml"; order=4)
+        c = Chart(unitsquare()) do (ξ,η)
+            @SVector[2ξ + η*(1-η), 3η+(1+η/2)*ξ^2]
+        end
+        Grids.jacobian(c::typeof(c), (ξ,η)) = @SMatrix[2 1-2η; (2+η)*ξ 3+ξ^2/2]
+
+        mg = equidistant_grid(c, 10,13)
+
+        @test inverse_inner_product(mg, stencil_set) isa LazyTensor{<:Any, 2,2}
+        @test_broken false # Test that it calculates the right thing
+    end
 end
--- a/test/SbpOperators/volumeops/laplace/laplace_test.jl	Fri Jun 28 17:02:47 2024 +0200
+++ b/test/SbpOperators/volumeops/laplace/laplace_test.jl	Fri Jun 28 17:04:05 2024 +0200
@@ -4,6 +4,8 @@
 using Sbplib.Grids
 using Sbplib.LazyTensors
 
+using StaticArrays
+
 @testset "Laplace" begin
     # Default stencils (4th order)
     operator_path = sbp_operators_path()*"standard_diagonal.toml"
@@ -72,12 +74,12 @@
     g_1D = equidistant_grid(0.0, 1., 101)
     g_3D = equidistant_grid((0.0, -1.0, 0.0), (1., 1., 1.), 51, 101, 52)
 
-    @testset "1D" begin
+    @testset "EquidistantGrid" begin
         Δ = laplace(g_1D, stencil_set)
         @test Δ == second_derivative(g_1D, stencil_set)
         @test Δ isa LazyTensor{Float64,1,1}
     end
-    @testset "3D" begin
+    @testset "TensorGrid" begin
         Δ = laplace(g_3D, stencil_set)
         @test Δ isa LazyTensor{Float64,3,3}
         Dxx = second_derivative(g_3D, stencil_set, 1)
@@ -86,6 +88,26 @@
         @test Δ == Dxx + Dyy + Dzz
         @test Δ isa LazyTensor{Float64,3,3}
     end
+
+    @testset "MappedGrid" begin
+        c = Chart(unitsquare()) do (ξ,η)
+            @SVector[2ξ + η*(1-η), 3η+(1+η/2)*ξ^2]
+        end
+        Grids.jacobian(c::typeof(c), (ξ,η)) = @SMatrix[2 1-2η; (2+η)*ξ 3+ξ^2/2]
+
+        g = equidistant_grid(c, 30,30)
+
+        @test laplace(g, stencil_set) isa LazyTensor{<:Any,2,2}
+
+        f((x,y)) = sin(4(x + y))
+        Δf((x,y)) = -16sin(4(x + y))
+        gf = map(f,g)
+
+        Δ = laplace(g, stencil_set)
+
+        @test collect(Δ*gf) isa Array{<:Any,2}
+        @test Δ*gf ≈ map(Δf, g)
+    end
 end
 
 @testset "sat_tensors" begin