Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/examples/p4est_3d_dgsem/elixir_mhd_shockcapturing_amr.jl
5586 views
1
using OrdinaryDiffEqLowStorageRK
2
using Trixi
3
4
###############################################################################
5
# semidiscretization of the compressible ideal GLM-MHD equations
6
7
equations = IdealGlmMhdEquations3D(1.4)
8
9
"""
10
initial_condition_blast_wave(x, t, equations::IdealGlmMhdEquations3D)
11
12
Weak magnetic blast wave setup taken from Section 6.1 of the paper:
13
- A. M. Rueda-Ramírez, S. Hennemann, F. J. Hindenlang, A. R. Winters, G. J. Gassner (2021)
14
An entropy stable nodal discontinuous Galerkin method for the resistive MHD
15
equations. Part II: Subcell finite volume shock capturing
16
[doi: 10.1016/j.jcp.2021.110580](https://doi.org/10.1016/j.jcp.2021.110580)
17
"""
18
function initial_condition_blast_wave(x, t, equations::IdealGlmMhdEquations3D)
19
# Center of the blast wave is selected for the domain [0, 3]^3
20
inicenter = (1.5, 1.5, 1.5)
21
x_norm = x[1] - inicenter[1]
22
y_norm = x[2] - inicenter[2]
23
z_norm = x[3] - inicenter[3]
24
r = sqrt(x_norm^2 + y_norm^2 + z_norm^2)
25
26
delta_0 = 0.1
27
r_0 = 0.3
28
lambda = exp(5.0 / delta_0 * (r - r_0))
29
30
prim_inner = SVector(1.2, 0.1, 0.0, 0.1, 0.9, 1.0, 1.0, 1.0, 0.0)
31
prim_outer = SVector(1.2, 0.2, -0.4, 0.2, 0.3, 1.0, 1.0, 1.0, 0.0)
32
prim_vars = (prim_inner + lambda * prim_outer) / (1.0 + lambda)
33
34
return prim2cons(prim_vars, equations)
35
end
36
initial_condition = initial_condition_blast_wave
37
38
# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of
39
# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.
40
# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.
41
# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.
42
# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.
43
# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the
44
# `StepsizeCallback` (CFL-Condition) and less diffusion.
45
surface_flux = (FluxLaxFriedrichs(max_abs_speed_naive), flux_nonconservative_powell)
46
volume_flux = (flux_hindenlang_gassner, flux_nonconservative_powell)
47
polydeg = 3
48
basis = LobattoLegendreBasis(polydeg)
49
indicator_sc = IndicatorHennemannGassner(equations, basis,
50
alpha_max = 0.5,
51
alpha_min = 0.001,
52
alpha_smooth = true,
53
variable = density_pressure)
54
volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;
55
volume_flux_dg = volume_flux,
56
volume_flux_fv = surface_flux)
57
58
solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux,
59
volume_integral = volume_integral)
60
61
# Mapping as described in https://arxiv.org/abs/2012.12040 but with slightly less warping.
62
# The mapping will be interpolated at tree level, and then refined without changing
63
# the geometry interpolant.
64
function mapping(xi_, eta_, zeta_)
65
# Transform input variables between -1 and 1 onto [0,3]
66
xi = 1.5 * xi_ + 1.5
67
eta = 1.5 * eta_ + 1.5
68
zeta = 1.5 * zeta_ + 1.5
69
70
y = eta +
71
3 / 11 * (cos(1.5 * pi * (2 * xi - 3) / 3) *
72
cos(0.5 * pi * (2 * eta - 3) / 3) *
73
cos(0.5 * pi * (2 * zeta - 3) / 3))
74
75
x = xi +
76
3 / 11 * (cos(0.5 * pi * (2 * xi - 3) / 3) *
77
cos(2 * pi * (2 * y - 3) / 3) *
78
cos(0.5 * pi * (2 * zeta - 3) / 3))
79
80
z = zeta +
81
3 / 11 * (cos(0.5 * pi * (2 * x - 3) / 3) *
82
cos(pi * (2 * y - 3) / 3) *
83
cos(0.5 * pi * (2 * zeta - 3) / 3))
84
85
return SVector(x, y, z)
86
end
87
88
trees_per_dimension = (2, 2, 2)
89
mesh = P4estMesh(trees_per_dimension,
90
polydeg = 3,
91
mapping = mapping,
92
initial_refinement_level = 2,
93
periodicity = true)
94
95
# create the semi discretization object
96
semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;
97
boundary_conditions = boundary_condition_periodic)
98
99
###############################################################################
100
# ODE solvers, callbacks etc.
101
102
tspan = (0.0, 0.5)
103
ode = semidiscretize(semi, tspan)
104
105
summary_callback = SummaryCallback()
106
107
analysis_interval = 100
108
analysis_callback = AnalysisCallback(semi, interval = analysis_interval)
109
110
alive_callback = AliveCallback(analysis_interval = analysis_interval)
111
112
amr_indicator = IndicatorLöhner(semi,
113
variable = density_pressure)
114
amr_controller = ControllerThreeLevel(semi, amr_indicator,
115
base_level = 2,
116
max_level = 4, max_threshold = 0.15)
117
amr_callback = AMRCallback(semi, amr_controller,
118
interval = 5,
119
adapt_initial_condition = true,
120
adapt_initial_condition_only_refine = true)
121
122
cfl = 1.4
123
stepsize_callback = StepsizeCallback(cfl = cfl)
124
125
glm_speed_callback = GlmSpeedCallback(glm_scale = 0.5, cfl = cfl)
126
127
callbacks = CallbackSet(summary_callback,
128
analysis_callback,
129
alive_callback,
130
amr_callback,
131
stepsize_callback,
132
glm_speed_callback)
133
134
###############################################################################
135
# run the simulation
136
137
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);
138
dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback
139
ode_default_options()..., callback = callbacks);
140
141