Mercurial > repos > public > sbplib_julia
comparison test/Grids/mapped_grid_test.jl @ 1695:a4c52ae93b11 feature/grids/manifolds
Merge feature/grids/curvilinear
author | Jonatan Werpers <jonatan@werpers.com> |
---|---|
date | Wed, 28 Aug 2024 10:35:08 +0200 |
parents | 13a7a4ff49e3 5bf4a35a78c5 |
children | 03894fd7b132 |
comparison
equal
deleted
inserted
replaced
1686:b996c35dd647 | 1695:a4c52ae93b11 |
---|---|
1 using Sbplib.Grids | 1 using Sbplib.Grids |
2 using Sbplib.RegionIndices | 2 using Sbplib.RegionIndices |
3 using Test | 3 using Test |
4 using StaticArrays | 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 | |
5 | 29 |
6 @testset "MappedGrid" begin | 30 @testset "MappedGrid" begin |
7 lg = equidistant_grid((0,0), (1,1), 11, 11) # TODO: Change dims of the grid to be different | 31 lg = equidistant_grid((0,0), (1,1), 11, 11) # TODO: Change dims of the grid to be different |
8 x̄ = map(ξ̄ -> 2ξ̄, lg) | 32 x̄ = map(ξ̄ -> 2ξ̄, lg) |
9 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) | 33 J = map(ξ̄ -> @SArray(fill(2., 2, 2)), lg) |
115 @test boundary_indices(mg, CartesianBoundary{2,Lower}()) == boundary_indices(lg,CartesianBoundary{2,Lower}()) | 139 @test boundary_indices(mg, CartesianBoundary{2,Lower}()) == boundary_indices(lg,CartesianBoundary{2,Lower}()) |
116 @test boundary_indices(mg, CartesianBoundary{1,Upper}()) == boundary_indices(lg,CartesianBoundary{1,Upper}()) | 140 @test boundary_indices(mg, CartesianBoundary{1,Upper}()) == boundary_indices(lg,CartesianBoundary{1,Upper}()) |
117 end | 141 end |
118 | 142 |
119 @testset "boundary_grid" begin | 143 @testset "boundary_grid" begin |
120 x̄((ξ, η)) = @SVector[ξ, η*(1+ξ*(ξ-1))] | 144 x̄, J = _partially_curved_mapping() |
121 J((ξ, η)) = @SMatrix[ | |
122 1 0; | |
123 η*(2ξ-1) 1+ξ*(ξ-1); | |
124 ] | |
125 | |
126 mg = mapped_grid(x̄, J, 10, 11) | 145 mg = mapped_grid(x̄, J, 10, 11) |
127 J1((ξ, η)) = @SMatrix[ | 146 J1((ξ, η)) = @SMatrix[ |
128 1 ; | 147 1 ; |
129 η*(2ξ-1); | 148 η*(2ξ-1); |
130 ] | 149 ] |
154 @testset test_boundary_grid(mg, TensorGridBoundary{1, Lower}(), J2) | 173 @testset test_boundary_grid(mg, TensorGridBoundary{1, Lower}(), J2) |
155 @testset test_boundary_grid(mg, TensorGridBoundary{1, Upper}(), J2) | 174 @testset test_boundary_grid(mg, TensorGridBoundary{1, Upper}(), J2) |
156 @testset test_boundary_grid(mg, TensorGridBoundary{2, Lower}(), J1) | 175 @testset test_boundary_grid(mg, TensorGridBoundary{2, Lower}(), J1) |
157 @testset test_boundary_grid(mg, TensorGridBoundary{2, Upper}(), J1) | 176 @testset test_boundary_grid(mg, TensorGridBoundary{2, Upper}(), J1) |
158 end | 177 end |
159 | |
160 @testset "jacobian_determinant" begin | |
161 @test_broken false | |
162 end | |
163 | |
164 @testset "geometric_tensor" begin | |
165 @test_broken false | |
166 end | |
167 | |
168 @testset "geometric_tensor_inverse" begin | |
169 @test_broken false | |
170 end | |
171 | |
172 end | 178 end |
173 | 179 |
174 @testset "mapped_grid" begin | 180 @testset "mapped_grid" begin |
175 x̄((ξ, η)) = @SVector[ξ, η*(1+ξ*(ξ-1))] | 181 x̄, J = _partially_curved_mapping() |
176 J((ξ, η)) = @SMatrix[ | |
177 1 0; | |
178 η*(2ξ-1) 1+ξ*(ξ-1); | |
179 ] | |
180 mg = mapped_grid(x̄, J, 10, 11) | 182 mg = mapped_grid(x̄, J, 10, 11) |
181 @test mg isa MappedGrid{SVector{2,Float64}, 2} | 183 @test mg isa MappedGrid{SVector{2,Float64}, 2} |
182 | 184 |
183 lg = equidistant_grid((0,0), (1,1), 10, 11) | 185 lg = equidistant_grid((0,0), (1,1), 10, 11) |
184 @test logicalgrid(mg) == lg | 186 @test logicalgrid(mg) == lg |
185 @test collect(mg) == map(x̄, lg) | 187 @test collect(mg) == map(x̄, lg) |
186 | 188 end |
187 | 189 |
188 @testset "normal" begin | 190 @testset "jacobian_determinant" begin |
189 @test normal(mg, CartesianBoundary{1,Lower}()) == fill(@SVector[-1,0], 11) | 191 @test_broken false |
190 @test normal(mg, CartesianBoundary{1,Upper}()) == fill(@SVector[1,0], 11) | 192 end |
191 @test normal(mg, CartesianBoundary{2,Lower}()) == fill(@SVector[0,-1], 10) | 193 |
192 @test normal(mg, CartesianBoundary{2,Upper}()) ≈ map(boundary_grid(mg,CartesianBoundary{2,Upper}())|>logicalgrid) do ξ̄ | 194 @testset "metric_tensor" begin |
193 α = 1-2ξ̄[1] | 195 @test_broken false |
194 @SVector[α,1]/√(α^2 + 1) | 196 end |
195 end | 197 |
196 end | 198 @testset "metric_tensor_inverse" begin |
197 end | 199 @test_broken false |
200 end | |
201 | |
202 @testset "min_spacing" begin | |
203 let g = mapped_grid(identity, x->@SMatrix[1], 11) | |
204 @test min_spacing(g) ≈ 0.1 | |
205 end | |
206 | |
207 let g = mapped_grid(x->x+x.^2/2, x->@SMatrix[1 .+ x], 11) | |
208 @test min_spacing(g) ≈ 0.105 | |
209 end | |
210 | |
211 let g = mapped_grid(x->x + x.*(1 .- x)/2, x->@SMatrix[1.5 .- x], 11) | |
212 @test min_spacing(g) ≈ 0.055 | |
213 end | |
214 | |
215 let g = mapped_grid(identity, x->@SMatrix[1 0; 0 1], 11,11) | |
216 @test min_spacing(g) ≈ 0.1 | |
217 end | |
218 | |
219 let g = mapped_grid(identity, x->@SMatrix[1 0; 0 1], 11,21) | |
220 @test min_spacing(g) ≈ 0.05 | |
221 end | |
222 | |
223 | |
224 @testset let a = @SVector[1,0], b = @SVector[1,1]/√2 | |
225 g = mapped_grid(_skew_mapping(a,b)...,11,11) | |
226 | |
227 @test min_spacing(g) ≈ 0.1*norm(b-a) | |
228 end | |
229 | |
230 @testset let a = @SVector[1,0], b = @SVector[-1,1]/√2 | |
231 g = mapped_grid(_skew_mapping(a,b)...,11,11) | |
232 | |
233 @test min_spacing(g) ≈ 0.1*norm(a+b) | |
234 end | |
235 end | |
236 | |
237 @testset "normal" begin | |
238 g = mapped_grid(_partially_curved_mapping()...,10, 11) | |
239 | |
240 @test normal(g, CartesianBoundary{1,Lower}()) == fill(@SVector[-1,0], 11) | |
241 @test normal(g, CartesianBoundary{1,Upper}()) == fill(@SVector[1,0], 11) | |
242 @test normal(g, CartesianBoundary{2,Lower}()) == fill(@SVector[0,-1], 10) | |
243 @test normal(g, CartesianBoundary{2,Upper}()) ≈ map(boundary_grid(g,CartesianBoundary{2,Upper}())|>logicalgrid) do ξ̄ | |
244 α = 1-2ξ̄[1] | |
245 @SVector[α,1]/√(α^2 + 1) | |
246 end | |
247 | |
248 g = mapped_grid(_fully_curved_mapping()...,5,4) | |
249 | |
250 unit(v) = v/norm(v) | |
251 @testset let bId = CartesianBoundary{1,Lower}() | |
252 lbg = boundary_grid(logicalgrid(g), bId) | |
253 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
254 -unit(@SVector[1/2, η/3-1/6]) | |
255 end | |
256 end | |
257 | |
258 @testset let bId = CartesianBoundary{1,Upper}() | |
259 lbg = boundary_grid(logicalgrid(g), bId) | |
260 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
261 unit(@SVector[7/2, 2η-1]/(5 + 3η + 2η^2)) | |
262 end | |
263 end | |
264 | |
265 @testset let bId = CartesianBoundary{2,Lower}() | |
266 lbg = boundary_grid(logicalgrid(g), bId) | |
267 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
268 -unit(@SVector[-2ξ, 2]/(6 + ξ^2 - 2ξ)) | |
269 end | |
270 end | |
271 | |
272 @testset let bId = CartesianBoundary{2,Upper}() | |
273 lbg = boundary_grid(logicalgrid(g), bId) | |
274 @test normal(g, bId) ≈ map(lbg) do (ξ, η) | |
275 unit(@SVector[-3ξ, 2]/(6 + ξ^2 + 3ξ)) | |
276 end | |
277 end | |
278 end |