Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/test/test_special_elixirs.jl
2055 views
1
module TestElixirs
2
3
using LinearAlgebra
4
using Test
5
using Trixi
6
using OrdinaryDiffEqSSPRK: SSPRK43
7
8
import ForwardDiff
9
10
include("test_trixi.jl")
11
12
# Start with a clean environment: remove Trixi.jl output directory if it exists
13
outdir = "out"
14
isdir(outdir) && rm(outdir, recursive = true)
15
16
EXAMPLES_DIR = examples_dir()
17
18
@testset "Special elixirs" begin
19
#! format: noindent
20
21
@testset "Convergence test" begin
22
@timed_testset "tree_2d_dgsem" begin
23
mean_convergence = convergence_test(@__MODULE__,
24
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
25
"elixir_advection_extended.jl"),
26
3, initial_refinement_level = 2)
27
@test isapprox(mean_convergence[:l2], [4.0], rtol = 0.05)
28
end
29
30
@timed_testset "structured_2d_dgsem" begin
31
mean_convergence = convergence_test(@__MODULE__,
32
joinpath(EXAMPLES_DIR,
33
"structured_2d_dgsem",
34
"elixir_advection_extended.jl"),
35
3, cells_per_dimension = (5, 9))
36
@test isapprox(mean_convergence[:l2], [4.0], rtol = 0.05)
37
end
38
39
@timed_testset "structured_2d_dgsem coupled" begin
40
mean_convergence = convergence_test(@__MODULE__,
41
joinpath(EXAMPLES_DIR,
42
"structured_2d_dgsem",
43
"elixir_advection_coupled.jl"),
44
3)
45
@test isapprox(mean_convergence[1][:l2], [4.0], rtol = 0.05)
46
@test isapprox(mean_convergence[2][:l2], [4.0], rtol = 0.05)
47
end
48
49
@timed_testset "p4est_2d_dgsem" begin
50
# Run convergence test on unrefined mesh
51
no_refine = @cfunction((p4est, which_tree, quadrant)->Cint(0), Cint,
52
(Ptr{Trixi.p4est_t}, Ptr{Trixi.p4est_topidx_t},
53
Ptr{Trixi.p4est_quadrant_t}))
54
mean_convergence = convergence_test(@__MODULE__,
55
joinpath(EXAMPLES_DIR, "p4est_2d_dgsem",
56
"elixir_euler_source_terms_nonconforming_unstructured_flag.jl"),
57
2, refine_fn_c = no_refine)
58
@test isapprox(mean_convergence[:linf], [3.2, 3.2, 4.0, 3.7], rtol = 0.05)
59
end
60
61
@timed_testset "structured_3d_dgsem" begin
62
mean_convergence = convergence_test(@__MODULE__,
63
joinpath(EXAMPLES_DIR,
64
"structured_3d_dgsem",
65
"elixir_advection_basic.jl"),
66
2, cells_per_dimension = (7, 4, 5))
67
@test isapprox(mean_convergence[:l2], [4.0], rtol = 0.05)
68
end
69
70
@timed_testset "p4est_3d_dgsem" begin
71
mean_convergence = convergence_test(@__MODULE__,
72
joinpath(EXAMPLES_DIR, "p4est_3d_dgsem",
73
"elixir_advection_unstructured_curved.jl"),
74
2, initial_refinement_level = 0)
75
@test isapprox(mean_convergence[:l2], [2.7], rtol = 0.05)
76
end
77
78
@timed_testset "paper_self_gravitating_gas_dynamics" begin
79
mean_convergence = convergence_test(@__MODULE__,
80
joinpath(EXAMPLES_DIR,
81
"paper_self_gravitating_gas_dynamics",
82
"elixir_eulergravity_convergence.jl"),
83
2, tspan = (0.0, 0.25),
84
initial_refinement_level = 1)
85
@test isapprox(mean_convergence[:l2], 4 * ones(4), atol = 0.4)
86
end
87
end
88
89
@timed_testset "Test linear structure (2D)" begin
90
trixi_include(@__MODULE__,
91
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
92
"elixir_advection_extended.jl"),
93
tspan = (0.0, 0.0), initial_refinement_level = 2)
94
A, b = linear_structure(semi)
95
λ = eigvals(Matrix(A))
96
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
97
98
trixi_include(@__MODULE__,
99
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
100
"elixir_hypdiff_lax_friedrichs.jl"),
101
tspan = (0.0, 0.0), initial_refinement_level = 2)
102
A, b = linear_structure(semi)
103
λ = eigvals(Matrix(A))
104
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
105
106
# check whether the user can modify `b` without changing `A`
107
x = vec(ode.u0)
108
Ax = A * x
109
@. b = 2 * b + x
110
@test A * x Ax
111
end
112
113
@testset "Test Jacobian of DG (1D)" begin
114
@timed_testset "TreeMesh: Linear advection" begin
115
trixi_include(@__MODULE__,
116
joinpath(EXAMPLES_DIR, "tree_1d_fdsbp",
117
"elixir_advection_upwind.jl"),
118
tspan = (0.0, 0.0))
119
120
A, _ = linear_structure(semi)
121
122
J = jacobian_ad_forward(semi)
123
@test Matrix(A) J
124
λ = eigvals(J)
125
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
126
127
J = jacobian_fd(semi)
128
@test Matrix(A) J
129
λ = eigvals(J)
130
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
131
132
# See https://github.com/trixi-framework/Trixi.jl/pull/2514
133
@test count(real.(λ) .>= -10) > 5
134
# See https://github.com/trixi-framework/Trixi.jl/pull/2522
135
t0 = zero(real(semi))
136
u0_ode = 1e9 * compute_coefficients(t0, semi)
137
J = jacobian_fd(semi; t0, u0_ode)
138
λ = eigvals(J)
139
@test count((-200 .<= real.(λ) .<= -10) .&& (-100 .<= imag.(λ) .<= 100)) == 0
140
@test count(isapprox.(imag.(λ), 0.0, atol = 10 * sqrt(eps(real(semi))))) == 2
141
end
142
143
@timed_testset "StructuredMesh: Compressible Euler equations" begin
144
trixi_include(@__MODULE__,
145
joinpath(EXAMPLES_DIR, "structured_1d_dgsem",
146
"elixir_euler_source_terms.jl"),
147
tspan = (0.0, 0.0))
148
149
J = jacobian_ad_forward(semi)
150
λ = eigvals(J)
151
@test maximum(real, λ) < 1e-13
152
153
J = jacobian_fd(semi)
154
λ = eigvals(J)
155
@test maximum(real, λ) < 5e-8
156
end
157
end
158
159
@testset "Test Jacobian of DG (2D)" begin
160
@timed_testset "TreeMesh: Linear advection" begin
161
trixi_include(@__MODULE__,
162
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
163
"elixir_advection_extended.jl"),
164
tspan = (0.0, 0.0), initial_refinement_level = 2)
165
A, _ = linear_structure(semi)
166
167
J = jacobian_ad_forward(semi)
168
@test Matrix(A) J
169
λ = eigvals(J)
170
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
171
172
J = jacobian_fd(semi)
173
@test Matrix(A) J
174
λ = eigvals(J)
175
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
176
end
177
178
@timed_testset "TreeMesh: Linear advection-diffusion" begin
179
trixi_include(@__MODULE__,
180
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
181
"elixir_advection_diffusion.jl"),
182
tspan = (0.0, 0.0), initial_refinement_level = 2)
183
184
J = jacobian_ad_forward(semi)
185
λ = eigvals(J)
186
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
187
188
J_parabolic = jacobian_ad_forward_parabolic(semi)
189
λ_parabolic = eigvals(J_parabolic)
190
# Parabolic spectrum is real and negative
191
@test maximum(real, λ_parabolic) < 10^(-14)
192
@test maximum(imag, λ_parabolic) < 10^(-14)
193
end
194
195
@timed_testset "TreeMesh: Compressible Euler equations" begin
196
trixi_include(@__MODULE__,
197
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
198
"elixir_euler_density_wave.jl"),
199
tspan = (0.0, 0.0), initial_refinement_level = 1)
200
201
J = jacobian_ad_forward(semi)
202
λ = eigvals(J)
203
@test maximum(real, λ) < 7.0e-7
204
205
J = jacobian_fd(semi)
206
λ = eigvals(J)
207
@test maximum(real, λ) < 7.0e-3
208
209
# This does not work yet because of the indicators...
210
@test_skip begin
211
trixi_include(@__MODULE__,
212
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
213
"elixir_euler_shockcapturing.jl"),
214
tspan = (0.0, 0.0), initial_refinement_level = 1)
215
jacobian_ad_forward(semi)
216
end
217
218
@timed_testset "DGMulti: Euler, weak form" begin
219
equations = CompressibleEulerEquations2D(1.4)
220
initial_condition = initial_condition_density_wave
221
222
solver = DGMulti(polydeg = 5, element_type = Quad(),
223
approximation_type = SBP(),
224
surface_integral = SurfaceIntegralWeakForm(flux_central),
225
volume_integral = VolumeIntegralWeakForm())
226
227
# DGMultiMesh is on [-1, 1]^ndims by default
228
cells_per_dimension = (2, 2)
229
mesh = DGMultiMesh(solver, cells_per_dimension,
230
periodicity = (true, true))
231
232
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition,
233
solver)
234
235
J = jacobian_ad_forward(semi)
236
λ = eigvals(J)
237
@test maximum(real, λ) < 7.0e-7
238
end
239
240
@timed_testset "DGMulti: Euler, SBP & flux differencing" begin
241
equations = CompressibleEulerEquations2D(1.4)
242
initial_condition = initial_condition_density_wave
243
244
solver = DGMulti(polydeg = 5, element_type = Quad(),
245
approximation_type = SBP(),
246
surface_integral = SurfaceIntegralWeakForm(flux_central),
247
volume_integral = VolumeIntegralFluxDifferencing(flux_central))
248
249
# DGMultiMesh is on [-1, 1]^ndims by default
250
cells_per_dimension = (2, 2)
251
mesh = DGMultiMesh(solver, cells_per_dimension,
252
periodicity = (true, true))
253
254
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition,
255
solver)
256
257
J = jacobian_ad_forward(semi)
258
λ = eigvals(J)
259
@test maximum(real, λ) < 7.0e-7
260
end
261
end
262
263
@timed_testset "TreeMesh: Navier-Stokes" begin
264
trixi_include(@__MODULE__,
265
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
266
"elixir_navierstokes_taylor_green_vortex.jl"),
267
tspan = (0.0, 0.0), initial_refinement_level = 2)
268
269
J = jacobian_ad_forward(semi)
270
λ = eigvals(J)
271
@test maximum(real, λ) < 0.2
272
273
J_parabolic = jacobian_ad_forward_parabolic(semi)
274
λ_parabolic = eigvals(J_parabolic)
275
# Parabolic spectrum is real and negative
276
@test maximum(real, λ_parabolic) < 10^(-16)
277
@test maximum(imag, λ_parabolic) < 10^(-15)
278
end
279
280
@timed_testset "TreeMesh: MHD" begin
281
trixi_include(@__MODULE__,
282
joinpath(EXAMPLES_DIR, "tree_2d_dgsem",
283
"elixir_mhd_alfven_wave.jl"),
284
tspan = (0.0, 0.0), initial_refinement_level = 0)
285
@test_nowarn jacobian_ad_forward(semi)
286
end
287
288
@timed_testset "UnstructuredMesh2D: Advection" begin
289
trixi_include(@__MODULE__,
290
joinpath(EXAMPLES_DIR, "unstructured_2d_dgsem",
291
"elixir_advection_basic.jl"),
292
tspan = (0.0, 0.0))
293
294
@test_nowarn jacobian_ad_forward(semi)
295
end
296
297
@timed_testset "TreeMesh: EulerGravity" begin
298
trixi_include(@__MODULE__,
299
joinpath(EXAMPLES_DIR,
300
"paper_self_gravitating_gas_dynamics",
301
"elixir_eulergravity_convergence.jl"),
302
tspan = (0.0, 0.0), initial_refinement_level = 1)
303
J = jacobian_ad_forward(semi)
304
λ = eigvals(J)
305
@test maximum(real, λ) < 1.5
306
end
307
308
@timed_testset "StructuredMesh: Polytropic Euler equations" begin
309
trixi_include(@__MODULE__,
310
joinpath(EXAMPLES_DIR, "structured_2d_dgsem",
311
"elixir_eulerpolytropic_wave.jl"),
312
cells_per_dimension = (6, 6),
313
tspan = (0.0, 0.0))
314
315
J = jacobian_ad_forward(semi)
316
λ = eigvals(J)
317
@test maximum(real, λ) < 0.05
318
end
319
320
@timed_testset "P4estMesh: Navier-Stokes" begin
321
trixi_include(@__MODULE__,
322
joinpath(EXAMPLES_DIR, "p4est_2d_dgsem",
323
"elixir_navierstokes_viscous_shock.jl"),
324
tspan = (0.0, 0.0))
325
326
J = jacobian_ad_forward(semi)
327
λ = eigvals(J)
328
@test maximum(real, λ) < 0.05
329
330
J_parabolic = jacobian_ad_forward_parabolic(semi)
331
λ_parabolic = eigvals(J_parabolic)
332
# Parabolic spectrum is real and negative
333
@test maximum(real, λ_parabolic) < 8e-14
334
@test maximum(imag, λ_parabolic) < 8e-14
335
end
336
337
@timed_testset "T8codeMesh: Advection" begin
338
trixi_include(@__MODULE__,
339
joinpath(EXAMPLES_DIR, "t8code_2d_dgsem",
340
"elixir_advection_unstructured_flag.jl"),
341
tspan = (0.0, 0.0), initial_refinement_level = 0,
342
polydeg = 2)
343
@test_nowarn jacobian_ad_forward(semi)
344
end
345
end
346
347
@timed_testset "Test linear structure (3D)" begin
348
trixi_include(@__MODULE__,
349
joinpath(EXAMPLES_DIR, "tree_3d_dgsem",
350
"elixir_advection_extended.jl"),
351
tspan = (0.0, 0.0), initial_refinement_level = 1)
352
A, b = linear_structure(semi)
353
λ = eigvals(Matrix(A))
354
@test maximum(real, λ) < 10 * sqrt(eps(real(semi)))
355
end
356
357
@timed_testset "Test Jacobian of DG (3D)" begin
358
@timed_testset "TreeMesh: Advection" begin
359
trixi_include(@__MODULE__,
360
joinpath(EXAMPLES_DIR, "tree_3d_dgsem",
361
"elixir_advection_extended.jl"),
362
tspan = (0.0, 0.0), initial_refinement_level = 1)
363
A, _ = linear_structure(semi)
364
365
J = jacobian_ad_forward(semi)
366
@test Matrix(A) J
367
368
J = jacobian_fd(semi)
369
@test Matrix(A) J
370
end
371
372
@timed_testset "StructuredMesh: MHD" begin
373
trixi_include(@__MODULE__,
374
joinpath(EXAMPLES_DIR, "structured_3d_dgsem",
375
"elixir_mhd_alfven_wave.jl"),
376
cells_per_dimension = (2, 2, 2),
377
polydeg = 2,
378
tspan = (0.0, 0.0))
379
380
@test_nowarn jacobian_ad_forward(semi)
381
end
382
383
@timed_testset "P4estMesh: Navier-Stokes" begin
384
trixi_include(@__MODULE__,
385
joinpath(EXAMPLES_DIR, "p4est_3d_dgsem",
386
"elixir_navierstokes_convergence.jl"),
387
initial_refinement_level = 0,
388
tspan = (0.0, 0.0))
389
390
@test_nowarn jacobian_ad_forward(semi)
391
@test_nowarn jacobian_ad_forward_parabolic(semi)
392
end
393
394
@timed_testset "T8CodeMesh: Advection" begin
395
trixi_include(@__MODULE__,
396
joinpath(EXAMPLES_DIR, "t8code_3d_dgsem",
397
"elixir_advection_cubed_sphere.jl"),
398
polydeg = 2,
399
trees_per_face_dimension = 3, layers = 2,
400
tspan = (0.0, 0.0))
401
402
@test_nowarn jacobian_ad_forward(semi)
403
end
404
end
405
406
@testset "AD using ForwardDiff" begin
407
@timed_testset "Euler equations 1D" begin
408
function entropy_at_final_time(k) # k is the wave number of the initial condition
409
equations = CompressibleEulerEquations1D(1.4)
410
mesh = TreeMesh((-1.0,), (1.0,), initial_refinement_level = 3,
411
n_cells_max = 10^4)
412
solver = DGSEM(3, FluxHLL(min_max_speed_naive),
413
VolumeIntegralFluxDifferencing(flux_ranocha))
414
initial_condition = (x, t, equations) -> begin
415
rho = 2 + sinpi(k * sum(x))
416
v1 = 0.1
417
p = 10.0
418
return prim2cons(SVector(rho, v1, p), equations)
419
end
420
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition,
421
solver,
422
uEltype = typeof(k))
423
ode = semidiscretize(semi, (0.0, 1.0))
424
summary_callback = SummaryCallback()
425
analysis_interval = 100
426
analysis_callback = AnalysisCallback(semi, interval = analysis_interval)
427
alive_callback = AliveCallback(analysis_interval = analysis_interval)
428
callbacks = CallbackSet(summary_callback,
429
analysis_callback,
430
alive_callback)
431
sol = solve(ode, SSPRK43(), callback = callbacks)
432
Trixi.integrate(entropy, sol.u[end], semi)
433
end
434
ForwardDiff.derivative(entropy_at_final_time, 1.0) -0.4524664696235628
435
end
436
437
@timed_testset "Linear advection 2D" begin
438
function energy_at_final_time(k) # k is the wave number of the initial condition
439
equations = LinearScalarAdvectionEquation2D(0.2, -0.7)
440
mesh = TreeMesh((-1.0, -1.0), (1.0, 1.0), initial_refinement_level = 3,
441
n_cells_max = 10^4)
442
solver = DGSEM(3, flux_lax_friedrichs)
443
initial_condition = (x, t, equation) -> begin
444
x_trans = Trixi.x_trans_periodic_2d(x - equation.advection_velocity * t)
445
return SVector(sinpi(k * sum(x_trans)))
446
end
447
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition,
448
solver,
449
uEltype = typeof(k))
450
ode = semidiscretize(semi, (0.0, 1.0))
451
summary_callback = SummaryCallback()
452
analysis_interval = 100
453
analysis_callback = AnalysisCallback(semi, interval = analysis_interval)
454
alive_callback = AliveCallback(analysis_interval = analysis_interval)
455
stepsize_callback = StepsizeCallback(cfl = 1.6)
456
callbacks = CallbackSet(summary_callback,
457
analysis_callback,
458
alive_callback,
459
stepsize_callback)
460
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);
461
ode_default_options()..., adaptive = false, dt = 1.0,
462
callback = callbacks)
463
Trixi.integrate(energy_total, sol.u[end], semi)
464
end
465
ForwardDiff.derivative(energy_at_final_time, 1.0) 1.4388628342896945e-5
466
end
467
468
@timed_testset "elixir_euler_ad.jl" begin
469
@test_trixi_include(joinpath(examples_dir(), "special_elixirs",
470
"elixir_euler_ad.jl"))
471
end
472
end
473
end
474
475
# Clean up afterwards: delete Trixi.jl output directory
476
@test_nowarn rm(outdir, recursive = true)
477
478
end #module
479
480