Mercurial > repos > public > sbplib_julia
changeset 1694:7a1eaaa0cb5a feature/grids/curvilinear
Merge default
author | Jonatan Werpers <jonatan@werpers.com> |
---|---|
date | Wed, 28 Aug 2024 10:20:28 +0200 |
parents | 5bf4a35a78c5 (diff) 5a1f51b4e3d9 (current diff) |
children | a4c52ae93b11 ff7aeda1b292 |
files | Manifest.toml |
diffstat | 8 files changed, 438 insertions(+), 2 deletions(-) [+] |
line wrap: on
line diff
--- a/Manifest.toml Wed Aug 28 10:12:37 2024 +0200 +++ b/Manifest.toml Wed Aug 28 10:20:28 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 Wed Aug 28 10:12:37 2024 +0200 +++ b/Project.toml Wed Aug 28 10:20:28 2024 +0200 @@ -4,6 +4,7 @@ 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"
--- a/src/Grids/Grids.jl Wed Aug 28 10:12:37 2024 +0200 +++ b/src/Grids/Grids.jl Wed Aug 28 10:20:28 2024 +0200 @@ -1,8 +1,10 @@ +# TODO: Double check that the interfaces for indexing and iterating are fully implemented and tested for all grids. module Grids using Sbplib.RegionIndices using Sbplib.LazyTensors using StaticArrays +using LinearAlgebra # Grid export Grid @@ -19,6 +21,7 @@ export eval_on export componentview export ArrayComponentView +export normal export BoundaryIdentifier export TensorGridBoundary @@ -33,11 +36,21 @@ export equidistant_grid +# MappedGrid +export MappedGrid +export jacobian +export logicalgrid +export mapped_grid +export jacobian_determinant +export metric_tensor +export metric_tensor_inverse + abstract type BoundaryIdentifier end include("grid.jl") include("tensor_grid.jl") include("equidistant_grid.jl") include("zero_dim_grid.jl") +include("mapped_grid.jl") end # module
--- a/src/Grids/equidistant_grid.jl Wed Aug 28 10:12:37 2024 +0200 +++ b/src/Grids/equidistant_grid.jl Wed Aug 28 10:20:28 2024 +0200 @@ -121,7 +121,7 @@ Constructs a 1D equidistant grid. """ -function equidistant_grid(limit_lower::T, limit_upper::T, size::Int) where T +function equidistant_grid(limit_lower::Number, limit_upper::Number, size::Int) if any(size .<= 0) throw(DomainError("size must be postive")) end
--- /dev/null Thu Jan 01 00:00:00 1970 +0000 +++ b/src/Grids/mapped_grid.jl Wed Aug 28 10:20:28 2024 +0200 @@ -0,0 +1,141 @@ +struct MappedGrid{T,D, GT<:Grid{<:Any,D}, CT<:AbstractArray{T,D}, JT<:AbstractArray{<:AbstractArray{<:Any, 2}, D}} <: Grid{T,D} + logicalgrid::GT + physicalcoordinates::CT + jacobian::JT +end + +jacobian(g::MappedGrid) = g.jacobian +logicalgrid(g::MappedGrid) = g.logicalgrid + + +# Indexing interface +Base.getindex(g::MappedGrid, I::Vararg{Int}) = g.physicalcoordinates[I...] +Base.eachindex(g::MappedGrid) = eachindex(g.logicalgrid) + +Base.firstindex(g::MappedGrid, d) = firstindex(g.logicalgrid, d) +Base.lastindex(g::MappedGrid, d) = lastindex(g.logicalgrid, d) +# TODO: axes(...)? + +# Iteration interface + +Base.iterate(g::MappedGrid) = iterate(g.physicalcoordinates) +Base.iterate(g::MappedGrid, state) = iterate(g.physicalcoordinates, state) + +Base.IteratorSize(::Type{<:MappedGrid{<:Any, D}}) where D = Base.HasShape{D}() +Base.length(g::MappedGrid) = length(g.logicalgrid) +Base.size(g::MappedGrid) = size(g.logicalgrid) +Base.size(g::MappedGrid, d) = size(g.logicalgrid, d) + +boundary_identifiers(g::MappedGrid) = boundary_identifiers(g.logicalgrid) +boundary_indices(g::MappedGrid, id::TensorGridBoundary) = boundary_indices(g.logicalgrid, id) + +function boundary_grid(g::MappedGrid, id::TensorGridBoundary) + b_indices = boundary_indices(g.logicalgrid, id) + + # Calculate indices of needed jacobian components + D = ndims(g) + all_indices = SVector{D}(1:D) + free_variable_indices = deleteat(all_indices, grid_id(id)) + jacobian_components = (:, free_variable_indices) + + # Create grid function for boundary grid jacobian + boundary_jacobian = componentview((@view g.jacobian[b_indices...]) , jacobian_components...) + boundary_physicalcoordinates = @view g.physicalcoordinates[b_indices...] + + return MappedGrid( + boundary_grid(g.logicalgrid, id), + boundary_physicalcoordinates, + boundary_jacobian, + ) +end + +# TBD: refine and coarsen could be implemented once we have a simple manifold implementation. +# Before we do, we should consider the overhead of including such a field in the mapped grid struct. + +function mapped_grid(x, J, size...) + D = length(size) + lg = equidistant_grid(ntuple(i->0., D), ntuple(i->1., D), size...) + return MappedGrid( + lg, + map(x,lg), + map(J,lg), + ) +end +# TODO: Delete this function + +function jacobian_determinant(g::MappedGrid) + return map(jacobian(g)) do ∂x∂ξ + det(∂x∂ξ) + end +end + +function metric_tensor(g::MappedGrid) + return map(jacobian(g)) do ∂x∂ξ + ∂x∂ξ'*∂x∂ξ + end +end + +function metric_tensor_inverse(g::MappedGrid) + return map(jacobian(g)) do ∂x∂ξ + inv(∂x∂ξ'*∂x∂ξ) + end +end + +function min_spacing(g::MappedGrid{T,1} where T) + n, = size(g) + + ms = Inf + for i ∈ 1:n-1 + ms = min(ms, norm(g[i+1]-g[i])) + end + + return ms +end + +function min_spacing(g::MappedGrid{T,2} where T) + n, m = size(g) + + ms = Inf + for i ∈ 1:n-1, j ∈ 1:m-1 # loop over each cell of the grid + + ms = min( + ms, + norm(g[i+1,j]-g[i,j]), + norm(g[i,j+1]-g[i,j]), + + norm(g[i+1,j]-g[i+1,j+1]), + norm(g[i,j+1]-g[i+1,j+1]), + + norm(g[i+1,j+1]-g[i,j]), + norm(g[i+1,j]-g[i,j+1]), + ) + # NOTE: This could be optimized to avoid checking all interior edges twice. + end + + return ms +end + +""" + normal(g::MappedGrid, boundary) + +The outward pointing normal as a grid function on the boundary +""" +function normal(g::MappedGrid, boundary) + b_indices = boundary_indices(g, boundary) + σ =_boundary_sign(component_type(g), boundary) + return map(jacobian(g)[b_indices...]) do ∂x∂ξ + ∂ξ∂x = inv(∂x∂ξ) + k = grid_id(boundary) + σ*∂ξ∂x[k,:]/norm(∂ξ∂x[k,:]) + end +end + +function _boundary_sign(T, boundary) + if boundary_id(boundary) == Upper() + return one(T) + elseif boundary_id(boundary) == Lower() + return -one(T) + else + throw(ArgumentError("The boundary identifier must be either `Lower()` or `Upper()`")) + end +end
--- a/src/Grids/tensor_grid.jl Wed Aug 28 10:12:37 2024 +0200 +++ b/src/Grids/tensor_grid.jl Wed Aug 28 10:20:28 2024 +0200 @@ -1,3 +1,5 @@ +# TODO: Check this file and other grids for duplicate implementation of general methods implemented for Grid + """ TensorGrid{T,D} <: Grid{T,D}
--- a/test/Grids/equidistant_grid_test.jl Wed Aug 28 10:12:37 2024 +0200 +++ b/test/Grids/equidistant_grid_test.jl Wed Aug 28 10:20:28 2024 +0200 @@ -112,6 +112,7 @@ @testset "equidistant_grid" begin @test equidistant_grid(0.0,1.0, 4) isa EquidistantGrid @test equidistant_grid((0.0,0.0),(8.0,5.0), 4, 3) isa TensorGrid + @test equidistant_grid((0.0,),(8.0,), 4) isa TensorGrid # constuctor @test_throws DomainError equidistant_grid(0.0, 1.0, 0)
--- /dev/null Thu Jan 01 00:00:00 1970 +0000 +++ b/test/Grids/mapped_grid_test.jl Wed Aug 28 10:20:28 2024 +0200 @@ -0,0 +1,278 @@ +using Sbplib.Grids +using Sbplib.RegionIndices +using Test +using StaticArrays +using LinearAlgebra + + +_skew_mapping(a,b) = (ξ̄->ξ̄[1]*a + ξ̄[2]*b, ξ̄->[a b]) + +function _partially_curved_mapping() + x̄((ξ, η)) = @SVector[ξ, η*(1+ξ*(ξ-1))] + J((ξ, η)) = @SMatrix[ + 1 0; + η*(2ξ-1) 1+ξ*(ξ-1); + ] + + return x̄, J +end + +function _fully_curved_mapping() + x̄((ξ, η)) = @SVector[2ξ + η*(1-η), 3η+(1+η/2)*ξ^2] + J((ξ, η)) = @SMatrix[ + 2 1-2η; + (2+η)*ξ 3+1/2*ξ^2; + ] + + return x̄, J +end + +@testset "MappedGrid" begin + lg = equidistant_grid((0,0), (1,1), 11, 11) # TODO: Change dims of the grid to be different + x̄ = map(ξ̄ -> 2ξ̄, lg) + J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) + mg = MappedGrid(lg, x̄, J) + + # TODO: Test constructor for different dims of range and domain for the coordinates + # TODO: Test constructor with different type than TensorGrid. a dummy type? + + @test_broken false # @test_throws ArgumentError("Sizes must match") MappedGrid(lg, map(ξ̄ -> @SArray[ξ̄[1], ξ̄[2], -ξ̄[1]], lg), rand(SMatrix{2,3,Float64},15,11)) + + + @test mg isa Grid{SVector{2, Float64},2} + + @test jacobian(mg) isa Array{<:AbstractMatrix} + @test logicalgrid(mg) isa Grid + + @testset "Indexing Interface" begin + mg = MappedGrid(lg, x̄, J) + @test mg[1,1] == [0.0, 0.0] + @test mg[4,2] == [0.6, 0.2] + @test mg[6,10] == [1., 1.8] + + @test mg[begin, begin] == [0.0, 0.0] + @test mg[end,end] == [2.0, 2.0] + @test mg[begin,end] == [0., 2.] + + @test eachindex(mg) == CartesianIndices((11,11)) + + @testset "cartesian indexing" begin + cases = [ + (1,1) , + (3,5) , + (10,6), + (1,1) , + (3,2) , + ] + + @testset "i = $is" for (lg, is) ∈ cases + @test mg[CartesianIndex(is...)] == mg[is...] + end + end + + @testset "eachindex" begin + @test eachindex(mg) == CartesianIndices((11,11)) + end + + @testset "firstindex" begin + @test firstindex(mg, 1) == 1 + @test firstindex(mg, 2) == 1 + end + + @testset "lastindex" begin + @test lastindex(mg, 1) == 11 + @test lastindex(mg, 2) == 11 + end + end + # TODO: Test with different types of logical grids + + @testset "Iterator interface" begin + sg = MappedGrid( + equidistant_grid((0,0), (1,1), 15, 11), + map(ξ̄ -> @SArray[ξ̄[1], ξ̄[2], -ξ̄[1]], lg), rand(SMatrix{2,3,Float64},15,11) + ) + + @test eltype(mg) == SVector{2,Float64} + @test eltype(sg) == SVector{3,Float64} + + @test eltype(typeof(mg)) == SVector{2,Float64} + @test eltype(typeof(sg)) == SVector{3,Float64} + + @test size(mg) == (11,11) + @test size(sg) == (15,11) + + @test size(mg,2) == 11 + @test size(sg,2) == 11 + + @test length(mg) == 121 + @test length(sg) == 165 + + @test Base.IteratorSize(mg) == Base.HasShape{2}() + @test Base.IteratorSize(typeof(mg)) == Base.HasShape{2}() + + @test Base.IteratorSize(sg) == Base.HasShape{2}() + @test Base.IteratorSize(typeof(sg)) == Base.HasShape{2}() + + element, state = iterate(mg) + @test element == lg[1,1].*2 + element, _ = iterate(mg, state) + @test element == lg[2,1].*2 + + element, state = iterate(sg) + @test element == sg.physicalcoordinates[1,1] + element, _ = iterate(sg, state) + @test element == sg.physicalcoordinates[2,1] + + @test collect(mg) == 2 .* lg + end + + @testset "Base" begin + @test ndims(mg) == 2 + end + + @testset "boundary_identifiers" begin + @test boundary_identifiers(mg) == boundary_identifiers(lg) + end + + @testset "boundary_indices" begin + @test boundary_indices(mg, CartesianBoundary{1,Lower}()) == boundary_indices(lg,CartesianBoundary{1,Lower}()) + @test boundary_indices(mg, CartesianBoundary{2,Lower}()) == boundary_indices(lg,CartesianBoundary{2,Lower}()) + @test boundary_indices(mg, CartesianBoundary{1,Upper}()) == boundary_indices(lg,CartesianBoundary{1,Upper}()) + end + + @testset "boundary_grid" begin + x̄, J = _partially_curved_mapping() + mg = mapped_grid(x̄, J, 10, 11) + J1((ξ, η)) = @SMatrix[ + 1 ; + η*(2ξ-1); + ] + J2((ξ, η)) = @SMatrix[ + 0; + 1+ξ*(ξ-1); + ] + + function test_boundary_grid(mg, bId, Jb) + bg = boundary_grid(mg, bId) + + lg = logicalgrid(mg) + expected_bg = MappedGrid( + boundary_grid(lg, bId), + map(x̄, boundary_grid(lg, bId)), + map(Jb, boundary_grid(lg, bId)), + ) + + @testset let bId=bId, bg=bg, expected_bg=expected_bg + @test collect(bg) == collect(expected_bg) + @test logicalgrid(bg) == logicalgrid(expected_bg) + @test jacobian(bg) == jacobian(expected_bg) + # TODO: Implement equality of a curvilinear grid and simlify the above + end + end + + @testset test_boundary_grid(mg, TensorGridBoundary{1, Lower}(), J2) + @testset test_boundary_grid(mg, TensorGridBoundary{1, Upper}(), J2) + @testset test_boundary_grid(mg, TensorGridBoundary{2, Lower}(), J1) + @testset test_boundary_grid(mg, TensorGridBoundary{2, Upper}(), J1) + end +end + +@testset "mapped_grid" begin + x̄, J = _partially_curved_mapping() + mg = mapped_grid(x̄, J, 10, 11) + @test mg isa MappedGrid{SVector{2,Float64}, 2} + + lg = equidistant_grid((0,0), (1,1), 10, 11) + @test logicalgrid(mg) == lg + @test collect(mg) == map(x̄, lg) +end + +@testset "jacobian_determinant" begin + @test_broken false +end + +@testset "metric_tensor" begin + @test_broken false +end + +@testset "metric_tensor_inverse" begin + @test_broken false +end + +@testset "min_spacing" begin + let g = mapped_grid(identity, x->@SMatrix[1], 11) + @test min_spacing(g) ≈ 0.1 + end + + let g = mapped_grid(x->x+x.^2/2, x->@SMatrix[1 .+ x], 11) + @test min_spacing(g) ≈ 0.105 + end + + let g = mapped_grid(x->x + x.*(1 .- x)/2, x->@SMatrix[1.5 .- x], 11) + @test min_spacing(g) ≈ 0.055 + end + + let g = mapped_grid(identity, x->@SMatrix[1 0; 0 1], 11,11) + @test min_spacing(g) ≈ 0.1 + end + + let g = mapped_grid(identity, x->@SMatrix[1 0; 0 1], 11,21) + @test min_spacing(g) ≈ 0.05 + end + + + @testset let a = @SVector[1,0], b = @SVector[1,1]/√2 + g = mapped_grid(_skew_mapping(a,b)...,11,11) + + @test min_spacing(g) ≈ 0.1*norm(b-a) + end + + @testset let a = @SVector[1,0], b = @SVector[-1,1]/√2 + g = mapped_grid(_skew_mapping(a,b)...,11,11) + + @test min_spacing(g) ≈ 0.1*norm(a+b) + end +end + +@testset "normal" begin + g = mapped_grid(_partially_curved_mapping()...,10, 11) + + @test normal(g, CartesianBoundary{1,Lower}()) == fill(@SVector[-1,0], 11) + @test normal(g, CartesianBoundary{1,Upper}()) == fill(@SVector[1,0], 11) + @test normal(g, CartesianBoundary{2,Lower}()) == fill(@SVector[0,-1], 10) + @test normal(g, CartesianBoundary{2,Upper}()) ≈ map(boundary_grid(g,CartesianBoundary{2,Upper}())|>logicalgrid) do ξ̄ + α = 1-2ξ̄[1] + @SVector[α,1]/√(α^2 + 1) + end + + g = mapped_grid(_fully_curved_mapping()...,5,4) + + unit(v) = v/norm(v) + @testset let bId = CartesianBoundary{1,Lower}() + lbg = boundary_grid(logicalgrid(g), bId) + @test normal(g, bId) ≈ map(lbg) do (ξ, η) + -unit(@SVector[1/2, η/3-1/6]) + end + end + + @testset let bId = CartesianBoundary{1,Upper}() + lbg = boundary_grid(logicalgrid(g), bId) + @test normal(g, bId) ≈ map(lbg) do (ξ, η) + unit(@SVector[7/2, 2η-1]/(5 + 3η + 2η^2)) + end + end + + @testset let bId = CartesianBoundary{2,Lower}() + lbg = boundary_grid(logicalgrid(g), bId) + @test normal(g, bId) ≈ map(lbg) do (ξ, η) + -unit(@SVector[-2ξ, 2]/(6 + ξ^2 - 2ξ)) + end + end + + @testset let bId = CartesianBoundary{2,Upper}() + lbg = boundary_grid(logicalgrid(g), bId) + @test normal(g, bId) ≈ map(lbg) do (ξ, η) + unit(@SVector[-3ξ, 2]/(6 + ξ^2 + 3ξ)) + end + end +end