Mercurial > repos > public > sbplib_julia
comparison test/Grids/mapped_grid_test.jl @ 2057:8a2a0d678d6f feature/lazy_tensors/pretty_printing
Merge default
| author | Jonatan Werpers <jonatan@werpers.com> |
|---|---|
| date | Tue, 10 Feb 2026 22:41:19 +0100 |
| parents | 04c251bccbd4 |
| children |
comparison
equal
deleted
inserted
replaced
| 1110:c0bff9f6e0fb | 2057:8a2a0d678d6f |
|---|---|
| 1 using Diffinitive.Grids | |
| 2 using Diffinitive.RegionIndices | |
| 3 using Test | |
| 4 using StaticArrays | |
| 5 using LinearAlgebra | |
| 6 | |
| 7 | |
| 8 _skew_mapping(a,b) = (ξ̄->ξ̄[1]*a + ξ̄[2]*b, ξ̄->[a b]) | |
| 9 | |
| 10 function _partially_curved_mapping() | |
| 11 x̄((ξ, η)) = @SVector[ξ, η*(1+ξ*(ξ-1))] | |
| 12 J((ξ, η)) = @SMatrix[ | |
| 13 1 0; | |
| 14 η*(2ξ-1) 1+ξ*(ξ-1); | |
| 15 ] | |
| 16 | |
| 17 return x̄, J | |
| 18 end | |
| 19 | |
| 20 function _fully_curved_mapping() | |
| 21 x̄((ξ, η)) = @SVector[2ξ + η*(1-η), 3η+(1+η/2)*ξ^2] | |
| 22 J((ξ, η)) = @SMatrix[ | |
| 23 2 1-2η; | |
| 24 (2+η)*ξ 3+1/2*ξ^2; | |
| 25 ] | |
| 26 | |
| 27 return x̄, J | |
| 28 end | |
| 29 | |
| 30 @testset "MappedGrid" begin | |
| 31 @testset "Constructor" begin | |
| 32 lg = equidistant_grid((0,0), (1,1), 11, 21) | |
| 33 | |
| 34 x̄ = map(ξ̄ -> 2ξ̄, lg) | |
| 35 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) | |
| 36 mg = MappedGrid(lg, x̄, J) | |
| 37 | |
| 38 @test mg isa Grid{SVector{2, Float64},2} | |
| 39 @test jacobian(mg) isa Array{<:AbstractMatrix} | |
| 40 @test logical_grid(mg) isa Grid | |
| 41 | |
| 42 @test collect(mg) == x̄ | |
| 43 @test jacobian(mg) == J | |
| 44 @test logical_grid(mg) == lg | |
| 45 | |
| 46 | |
| 47 x̄ = map(ξ̄ -> @SVector[ξ̄[1],ξ̄[2], ξ̄[1] + ξ̄[2]], lg) | |
| 48 J = map(ξ̄ -> @SMatrix[1 0; 0 1; 1 1], lg) | |
| 49 mg = MappedGrid(lg, x̄, J) | |
| 50 | |
| 51 @test mg isa Grid{SVector{3, Float64},2} | |
| 52 @test jacobian(mg) isa Array{<:AbstractMatrix} | |
| 53 @test logical_grid(mg) isa Grid | |
| 54 | |
| 55 @test collect(mg) == x̄ | |
| 56 @test jacobian(mg) == J | |
| 57 @test logical_grid(mg) == lg | |
| 58 | |
| 59 sz1 = (10,11) | |
| 60 sz2 = (10,12) | |
| 61 @test_throws ArgumentError("Sizes must match") MappedGrid( | |
| 62 equidistant_grid((0,0), (1,1), sz2...), | |
| 63 rand(SVector{2},sz1...), | |
| 64 rand(SMatrix{2,2},sz1...), | |
| 65 ) | |
| 66 | |
| 67 @test_throws ArgumentError("Sizes must match") MappedGrid( | |
| 68 equidistant_grid((0,0), (1,1), sz1...), | |
| 69 rand(SVector{2},sz2...), | |
| 70 rand(SMatrix{2,2},sz1...), | |
| 71 ) | |
| 72 | |
| 73 @test_throws ArgumentError("Sizes must match") MappedGrid( | |
| 74 equidistant_grid((0,0), (1,1), sz1...), | |
| 75 rand(SVector{2},sz1...), | |
| 76 rand(SMatrix{2,2},sz2...), | |
| 77 ) | |
| 78 | |
| 79 err_str = "The size of the jacobian must match the dimensions of the grid and coordinates" | |
| 80 @test_throws ArgumentError(err_str) MappedGrid( | |
| 81 equidistant_grid((0,0), (1,1), 10, 11), | |
| 82 rand(SVector{3}, 10, 11), | |
| 83 rand(SMatrix{3,4}, 10, 11), | |
| 84 ) | |
| 85 | |
| 86 @test_throws ArgumentError(err_str) MappedGrid( | |
| 87 equidistant_grid((0,0), (1,1), 10, 11), | |
| 88 rand(SVector{3}, 10, 11), | |
| 89 rand(SMatrix{4,2}, 10, 11), | |
| 90 ) | |
| 91 end | |
| 92 | |
| 93 @testset "Indexing Interface" begin | |
| 94 lg = equidistant_grid((0,0), (1,1), 11, 21) | |
| 95 x̄ = map(ξ̄ -> 2ξ̄, lg) | |
| 96 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) | |
| 97 mg = MappedGrid(lg, x̄, J) | |
| 98 @test mg[1,1] == [0.0, 0.0] | |
| 99 @test mg[4,2] == [0.6, 0.1] | |
| 100 @test mg[6,10] == [1., 0.9] | |
| 101 | |
| 102 @test mg[begin, begin] == [0.0, 0.0] | |
| 103 @test mg[end,end] == [2.0, 2.0] | |
| 104 @test mg[begin,end] == [0., 2.] | |
| 105 | |
| 106 @test axes(mg) == (1:11, 1:21) | |
| 107 | |
| 108 @testset "cartesian indexing" begin | |
| 109 cases = [ | |
| 110 (1,1) , | |
| 111 (3,5) , | |
| 112 (10,6), | |
| 113 (1,1) , | |
| 114 (3,2) , | |
| 115 ] | |
| 116 | |
| 117 @testset "i = $is" for (lg, is) ∈ cases | |
| 118 @test mg[CartesianIndex(is...)] == mg[is...] | |
| 119 end | |
| 120 end | |
| 121 | |
| 122 @testset "eachindex" begin | |
| 123 @test eachindex(mg) == CartesianIndices((11,21)) | |
| 124 end | |
| 125 | |
| 126 @testset "firstindex" begin | |
| 127 @test firstindex(mg, 1) == 1 | |
| 128 @test firstindex(mg, 2) == 1 | |
| 129 end | |
| 130 | |
| 131 @testset "lastindex" begin | |
| 132 @test lastindex(mg, 1) == 11 | |
| 133 @test lastindex(mg, 2) == 21 | |
| 134 end | |
| 135 end | |
| 136 | |
| 137 @testset "Iterator interface" begin | |
| 138 lg = equidistant_grid((0,0), (1,1), 11, 21) | |
| 139 x̄ = map(ξ̄ -> 2ξ̄, lg) | |
| 140 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) | |
| 141 | |
| 142 mg = MappedGrid(lg, x̄, J) | |
| 143 | |
| 144 lg2 = equidistant_grid((0,0), (1,1), 15, 11) | |
| 145 sg = MappedGrid( | |
| 146 equidistant_grid((0,0), (1,1), 15, 11), | |
| 147 map(ξ̄ -> @SArray[ξ̄[1], ξ̄[2], -ξ̄[1]], lg2), rand(SMatrix{3,2,Float64},15,11) | |
| 148 ) | |
| 149 | |
| 150 @test eltype(mg) == SVector{2,Float64} | |
| 151 @test eltype(sg) == SVector{3,Float64} | |
| 152 | |
| 153 @test eltype(typeof(mg)) == SVector{2,Float64} | |
| 154 @test eltype(typeof(sg)) == SVector{3,Float64} | |
| 155 | |
| 156 @test size(mg) == (11,21) | |
| 157 @test size(sg) == (15,11) | |
| 158 | |
| 159 @test size(mg,2) == 21 | |
| 160 @test size(sg,2) == 11 | |
| 161 | |
| 162 @test length(mg) == 231 | |
| 163 @test length(sg) == 165 | |
| 164 | |
| 165 @test Base.IteratorSize(mg) == Base.HasShape{2}() | |
| 166 @test Base.IteratorSize(typeof(mg)) == Base.HasShape{2}() | |
| 167 | |
| 168 @test Base.IteratorSize(sg) == Base.HasShape{2}() | |
| 169 @test Base.IteratorSize(typeof(sg)) == Base.HasShape{2}() | |
| 170 | |
| 171 element, state = iterate(mg) | |
| 172 @test element == lg[1,1].*2 | |
| 173 element, _ = iterate(mg, state) | |
| 174 @test element == lg[2,1].*2 | |
| 175 | |
| 176 element, state = iterate(sg) | |
| 177 @test element == sg.physicalcoordinates[1,1] | |
| 178 element, _ = iterate(sg, state) | |
| 179 @test element == sg.physicalcoordinates[2,1] | |
| 180 | |
| 181 @test collect(mg) == 2 .* lg | |
| 182 end | |
| 183 | |
| 184 @testset "Base" begin | |
| 185 lg = equidistant_grid((0,0), (1,1), 11, 21) | |
| 186 x̄ = map(ξ̄ -> 2ξ̄, lg) | |
| 187 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) | |
| 188 mg = MappedGrid(lg, x̄, J) | |
| 189 | |
| 190 @test ndims(mg) == 2 | |
| 191 end | |
| 192 | |
| 193 @testset "==" begin | |
| 194 sz = (15,11) | |
| 195 lg = equidistant_grid((0,0), (1,1), sz...) | |
| 196 x = rand(SVector{3,Float64}, sz...) | |
| 197 J = rand(SMatrix{3,2,Float64}, sz...) | |
| 198 | |
| 199 sg = MappedGrid(lg, x, J) | |
| 200 | |
| 201 sg1 = MappedGrid(equidistant_grid((0,0), (1,1), sz...), copy(x), copy(J)) | |
| 202 | |
| 203 sz2 = (15,12) | |
| 204 lg2 = equidistant_grid((0,0), (1,1), sz2...) | |
| 205 x2 = rand(SVector{3,Float64}, sz2...) | |
| 206 J2 = rand(SMatrix{3,2,Float64}, sz2...) | |
| 207 sg2 = MappedGrid(lg2, x2, J2) | |
| 208 | |
| 209 sg3 = MappedGrid(lg, rand(SVector{3,Float64}, sz...), J) | |
| 210 sg4 = MappedGrid(lg, x, rand(SMatrix{3,2,Float64}, sz...)) | |
| 211 | |
| 212 @test sg == sg1 | |
| 213 @test sg != sg2 # Different size | |
| 214 @test sg != sg3 # Different coordinates | |
| 215 @test sg != sg4 # Different jacobian | |
| 216 end | |
| 217 | |
| 218 @testset "boundary_identifiers" begin | |
| 219 lg = equidistant_grid((0,0), (1,1), 11, 15) | |
| 220 x̄ = map(ξ̄ -> 2ξ̄, lg) | |
| 221 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) | |
| 222 mg = MappedGrid(lg, x̄, J) | |
| 223 @test boundary_identifiers(mg) == boundary_identifiers(lg) | |
| 224 end | |
| 225 | |
| 226 @testset "boundary_indices" begin | |
| 227 lg = equidistant_grid((0,0), (1,1), 11, 15) | |
| 228 x̄ = map(ξ̄ -> 2ξ̄, lg) | |
| 229 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) | |
| 230 mg = MappedGrid(lg, x̄, J) | |
| 231 | |
| 232 @test boundary_indices(mg, CartesianBoundary{1,LowerBoundary}()) == boundary_indices(lg,CartesianBoundary{1,LowerBoundary}()) | |
| 233 @test boundary_indices(mg, CartesianBoundary{2,LowerBoundary}()) == boundary_indices(lg,CartesianBoundary{2,LowerBoundary}()) | |
| 234 @test boundary_indices(mg, CartesianBoundary{1,UpperBoundary}()) == boundary_indices(lg,CartesianBoundary{1,UpperBoundary}()) | |
| 235 end | |
| 236 | |
| 237 @testset "boundary_grid" begin | |
| 238 x̄, J = _partially_curved_mapping() | |
| 239 mg = mapped_grid(x̄, J, 10, 11) | |
| 240 J1((ξ, η)) = @SMatrix[ | |
| 241 1 ; | |
| 242 η*(2ξ-1); | |
| 243 ] | |
| 244 J2((ξ, η)) = @SMatrix[ | |
| 245 0; | |
| 246 1+ξ*(ξ-1); | |
| 247 ] | |
| 248 | |
| 249 function expected_bg(mg, bId, Jb) | |
| 250 lg = logical_grid(mg) | |
| 251 return MappedGrid( | |
| 252 boundary_grid(lg, bId), | |
| 253 map(x̄, boundary_grid(lg, bId)), | |
| 254 map(Jb, boundary_grid(lg, bId)), | |
| 255 ) | |
| 256 end | |
| 257 | |
| 258 let bid = TensorGridBoundary{1, LowerBoundary}() | |
| 259 @test boundary_grid(mg, bid) == expected_bg(mg, bid, J2) | |
| 260 end | |
| 261 | |
| 262 let bid = TensorGridBoundary{1, UpperBoundary}() | |
| 263 @test boundary_grid(mg, bid) == expected_bg(mg, bid, J2) | |
| 264 end | |
| 265 | |
| 266 let bid = TensorGridBoundary{2, LowerBoundary}() | |
| 267 @test boundary_grid(mg, bid) == expected_bg(mg, bid, J1) | |
| 268 end | |
| 269 | |
| 270 let bid = TensorGridBoundary{2, UpperBoundary}() | |
| 271 @test boundary_grid(mg, bid) == expected_bg(mg, bid, J1) | |
| 272 end | |
| 273 end | |
| 274 end | |
| 275 | |
| 276 @testset "mapped_grid" begin | |
| 277 x̄, J = _partially_curved_mapping() | |
| 278 mg = mapped_grid(x̄, J, 10, 11) | |
| 279 @test mg isa MappedGrid{SVector{2,Float64}, 2} | |
| 280 | |
| 281 lg = equidistant_grid(unitsquare(), 10, 11) | |
| 282 @test logical_grid(mg) == lg | |
| 283 @test collect(mg) == map(x̄, lg) | |
| 284 | |
| 285 @test mapped_grid(x̄, J, lg) == mg | |
| 286 | |
| 287 @test mapped_grid(x̄, J, unitsquare(), 10, 11) == mg | |
| 288 end | |
| 289 | |
| 290 @testset "metric_tensor" begin | |
| 291 x̄((ξ, η)) = @SVector[ξ*η, ξ + η^2] | |
| 292 J((ξ, η)) = @SMatrix[ | |
| 293 η ξ; | |
| 294 1 2η; | |
| 295 ] | |
| 296 | |
| 297 g = mapped_grid(x̄, J, 10, 11) | |
| 298 G = map(logical_grid(g)) do (ξ,η) | |
| 299 @SMatrix[ | |
| 300 1+η^2 ξ*η+2η; | |
| 301 ξ*η+2η ξ^2 + 4η^2; | |
| 302 ] | |
| 303 end | |
| 304 @test metric_tensor(g) ≈ G | |
| 305 end | |
| 306 | |
| 307 @testset "min_spacing" begin | |
| 308 let g = mapped_grid(identity, x->@SMatrix[1], 11) | |
| 309 @test min_spacing(g) ≈ 0.1 | |
| 310 end | |
| 311 | |
| 312 let g = mapped_grid(x->x+x.^2/2, x->@SMatrix[1 .+ x], 11) | |
| 313 @test min_spacing(g) ≈ 0.105 | |
| 314 end | |
| 315 | |
| 316 let g = mapped_grid(x->x + x.*(1 .- x)/2, x->@SMatrix[1.5 .- x], 11) | |
| 317 @test min_spacing(g) ≈ 0.055 | |
| 318 end | |
| 319 | |
| 320 let g = mapped_grid(identity, x->@SMatrix[1 0; 0 1], 11,11) | |
| 321 @test min_spacing(g) ≈ 0.1 | |
| 322 end | |
| 323 | |
| 324 let g = mapped_grid(identity, x->@SMatrix[1 0; 0 1], 11,21) | |
| 325 @test min_spacing(g) ≈ 0.05 | |
| 326 end | |
| 327 | |
| 328 | |
| 329 @testset let a = @SVector[1,0], b = @SVector[1,1]/√2 | |
| 330 g = mapped_grid(_skew_mapping(a,b)...,11,11) | |
| 331 | |
| 332 @test min_spacing(g) ≈ 0.1*norm(b-a) | |
| 333 end | |
| 334 | |
| 335 @testset let a = @SVector[1,0], b = @SVector[-1,1]/√2 | |
| 336 g = mapped_grid(_skew_mapping(a,b)...,11,11) | |
| 337 | |
| 338 @test min_spacing(g) ≈ 0.1*norm(a+b) | |
| 339 end | |
| 340 end | |
| 341 | |
| 342 @testset "normal" begin | |
| 343 g = mapped_grid(_partially_curved_mapping()...,10, 11) | |
| 344 | |
| 345 @test normal(g, CartesianBoundary{1,LowerBoundary}()) == fill(@SVector[-1,0], 11) | |
| 346 @test normal(g, CartesianBoundary{1,UpperBoundary}()) == fill(@SVector[1,0], 11) | |
| 347 @test normal(g, CartesianBoundary{2,LowerBoundary}()) == fill(@SVector[0,-1], 10) | |
| 348 @test normal(g, CartesianBoundary{2,UpperBoundary}()) ≈ map(boundary_grid(g,CartesianBoundary{2,UpperBoundary}())|>logical_grid) do ξ̄ | |
| 349 α = 1-2ξ̄[1] | |
| 350 @SVector[α,1]/√(α^2 + 1) | |
| 351 end | |
| 352 | |
| 353 g = mapped_grid(_fully_curved_mapping()...,5,4) | |
| 354 | |
| 355 unit(v) = v/norm(v) | |
| 356 @testset let bId = CartesianBoundary{1,LowerBoundary}() | |
| 357 lbg = boundary_grid(logical_grid(g), bId) | |
| 358 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
| 359 -unit(@SVector[1/2, η/3-1/6]) | |
| 360 end | |
| 361 end | |
| 362 | |
| 363 @testset let bId = CartesianBoundary{1,UpperBoundary}() | |
| 364 lbg = boundary_grid(logical_grid(g), bId) | |
| 365 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
| 366 unit(@SVector[7/2, 2η-1]/(5 + 3η + 2η^2)) | |
| 367 end | |
| 368 end | |
| 369 | |
| 370 @testset let bId = CartesianBoundary{2,LowerBoundary}() | |
| 371 lbg = boundary_grid(logical_grid(g), bId) | |
| 372 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
| 373 -unit(@SVector[-2ξ, 2]/(6 + ξ^2 - 2ξ)) | |
| 374 end | |
| 375 end | |
| 376 | |
| 377 @testset let bId = CartesianBoundary{2,UpperBoundary}() | |
| 378 lbg = boundary_grid(logical_grid(g), bId) | |
| 379 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
| 380 unit(@SVector[-3ξ, 2]/(6 + ξ^2 + 3ξ)) | |
| 381 end | |
| 382 end | |
| 383 end |
