Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Added Shu-Osher initialization for 1D compressible Euler with Gauss nodes #1943

Merged
merged 8 commits into from
May 16, 2024
Original file line number Diff line number Diff line change
@@ -0,0 +1,93 @@
using Trixi
using OrdinaryDiffEq

gamma_gas = 1.4
equations = CompressibleEulerEquations1D(gamma_gas)

###############################################################################
# setup the GSBP DG discretization that uses the Gauss operators from
# Chan, Del Rey Fernandez, Carpenter (2019).
# [https://doi.org/10.1137/18M1209234](https://doi.org/10.1137/18M1209234)

# Shu-Osher initial condition for 1D compressible Euler equations
# Example 8 from Shu, Osher (1989).
# [https://doi.org/10.1016/0021-9991(89)90222-2](https://doi.org/10.1016/0021-9991(89)90222-2)
function initial_condition_shu_osher(x, t, equations::CompressibleEulerEquations1D)
x0 = -4

rho_left = 27 / 7
v_left = 4 * sqrt(35) / 9
p_left = 31 / 3

# Replaced v_right = 0 to v_right = 0.1 to avoid positivity issues.
v_right = 0.1
p_right = 1.0

rho = ifelse(x[1] > x0, 1 + 1 / 5 * sin(5 * x[1]), rho_left)
v = ifelse(x[1] > x0, v_right, v_left)
p = ifelse(x[1] > x0, p_right, p_left)

return prim2cons(SVector(rho + 0.1 * rand(), v + 0.1 * rand(), p + 0.1 * rand()),
equations)
end

initial_condition = initial_condition_shu_osher

surface_flux = flux_lax_friedrichs
volume_flux = flux_ranocha

polydeg = 3
basis = DGMultiBasis(Line(), polydeg, approximation_type = GaussSBP())

indicator_sc = IndicatorHennemannGassner(equations, basis,
alpha_max = 0.5,
alpha_min = 0.001,
alpha_smooth = true,
variable = density_pressure)
volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;
volume_flux_dg = volume_flux,
volume_flux_fv = surface_flux)

dg = DGMulti(basis,
surface_integral = SurfaceIntegralWeakForm(surface_flux),
volume_integral = volume_integral)

boundary_condition = BoundaryConditionDirichlet(initial_condition)
boundary_conditions = (; :entire_boundary => boundary_condition)

###############################################################################
# setup the 1D mesh

cells_per_dimension = (64,)
mesh = DGMultiMesh(dg, cells_per_dimension,
coordinates_min = (-5.0,), coordinates_max = (5.0,),
periodicity = false)

###############################################################################
# setup the semidiscretization and ODE problem

semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition,
dg, boundary_conditions = boundary_conditions)

tspan = (0.0, 1.0)
ode = semidiscretize(semi, tspan)

###############################################################################
# setup the callbacks

# prints a summary of the simulation setup and resets the timers
summary_callback = SummaryCallback()

# analyse the solution in regular intervals and prints the results
analysis_callback = AnalysisCallback(semi, interval = 100, uEltype = real(dg))

# handles the re-calculation of the maximum Δt after each time step
stepsize_callback = StepsizeCallback(cfl = 0.1)

# collect all callbacks such that they can be passed to the ODE solver
callbacks = CallbackSet(summary_callback, analysis_callback, stepsize_callback)

# ###############################################################################
# # run the simulation

sol = solve(ode, SSPRK43(), adaptive = true, callback = callbacks, save_everystep = false)
49 changes: 37 additions & 12 deletions test/test_dgmulti_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -52,12 +52,37 @@ end
l2=[
7.853842541289665e-7,
9.609905503440606e-7,
2.832322219966481e-6,
2.832322219966481e-6
] ./ sqrt(2.0),
linf=[
1.5003758788711963e-6,
1.802998748523521e-6,
4.83599270806323e-6,
4.83599270806323e-6
])
# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
let
t = sol.t[end]
u_ode = sol.u[end]
du_ode = similar(u_ode)
@test (@allocated Trixi.rhs!(du_ode, u_ode, semi, t)) < 1000
end
end

@trixi_testset "elixir_euler_shu_osher_gauss_shock_capturing.jl " begin
@test_trixi_include(joinpath(EXAMPLES_DIR,
"elixir_euler_shu_osher_gauss_shock_capturing.jl"),
cells_per_dimension=(64,), tspan=(0.0, 1.0),
# division by sqrt(2.0) corresponds to normalization by the square root of the size of the domain
ranocha marked this conversation as resolved.
Show resolved Hide resolved
JoshuaLampert marked this conversation as resolved.
Show resolved Hide resolved
l2=[
1.7177729727131328,
6.191308529732632,
22.281999765458117
],
linf=[
3.228727516356452,
10.989654948503484,
38.65258452259893
])
# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
Expand Down Expand Up @@ -93,12 +118,12 @@ end
l2=[
6.437827414849647e-6,
2.1840558851820947e-6,
1.3245669629438228e-5,
1.3245669629438228e-5
],
linf=[
2.0715843751295537e-5,
8.519520630301258e-6,
4.2642194098885255e-5,
4.2642194098885255e-5
])
# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
Expand All @@ -122,12 +147,12 @@ end
l2=[
1.8684509287853788e-5,
1.0641411823379635e-5,
5.178010291876143e-5,
5.178010291876143e-5
],
linf=[
6.933493585936645e-5,
3.0277366229292113e-5,
0.0002220020568932668,
0.0002220020568932668
])
show(stdout, semi.solver.basis)
show(stdout, MIME"text/plain"(), semi.solver.basis)
Expand All @@ -145,11 +170,11 @@ end
@test_trixi_include(joinpath(EXAMPLES_DIR, "elixir_euler_fdsbp_periodic.jl"),
l2=[
9.146929178341782e-7, 1.8997616876521201e-6,
3.991417701005622e-6,
3.991417701005622e-6
],
linf=[
1.7321089882393892e-6, 3.3252888869128583e-6,
6.525278767988141e-6,
6.525278767988141e-6
])
show(stdout, semi.solver.basis)
show(stdout, MIME"text/plain"(), semi.solver.basis)
Expand Down Expand Up @@ -186,13 +211,13 @@ end
3.03001101100507e-6,
1.692177335948727e-5,
3.002634351734614e-16,
1.1636653574178203e-15,
1.1636653574178203e-15
],
linf=[
1.2043401988570679e-5,
5.346847010329059e-5,
9.43689570931383e-16,
2.220446049250313e-15,
2.220446049250313e-15
])
# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
Expand All @@ -212,13 +237,13 @@ end
1.633271343738687e-5,
9.575385661756332e-6,
1.2700331443128421e-5,
0.0,
0.0
],
linf=[
7.304984704381567e-5,
5.2365944135601694e-5,
6.469559594934893e-5,
0.0,
0.0
])
# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
Expand Down