Path: blob/main/examples/tree_1d_dgsem/elixir_eulergravity_convergence.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations5gamma = 2.06equations_euler = CompressibleEulerEquations1D(gamma)78initial_condition = initial_condition_eoc_test_coupled_euler_gravity910polydeg = 311solver_euler = DGSEM(polydeg, flux_hll)1213coordinates_min = 0.014coordinates_max = 2.015mesh = TreeMesh(coordinates_min, coordinates_max,16initial_refinement_level = 2,17n_cells_max = 10_000)1819semi_euler = SemidiscretizationHyperbolic(mesh, equations_euler, initial_condition,20solver_euler)2122###############################################################################23# semidiscretization of the hyperbolic diffusion equations24equations_gravity = HyperbolicDiffusionEquations1D()2526# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of27# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.28# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.29# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.30# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.31# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the32# `StepsizeCallback` (CFL-Condition) and less diffusion.33solver_gravity = DGSEM(polydeg, FluxLaxFriedrichs(max_abs_speed_naive))3435semi_gravity = SemidiscretizationHyperbolic(mesh, equations_gravity, initial_condition,36solver_gravity,37source_terms = source_terms_harmonic)3839###############################################################################40# combining both semidiscretizations for Euler + self-gravity41parameters = ParametersEulerGravity(background_density = 2.0, # aka rho042gravitational_constant = 1.0, # aka G43cfl = 1.5,44resid_tol = 1.0e-10,45n_iterations_max = 1000,46timestep_gravity = timestep_gravity_erk52_3Sstar!)4748semi = SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)4950###############################################################################51# ODE solvers, callbacks etc.52tspan = (0.0, 0.5)53ode = semidiscretize(semi, tspan)5455summary_callback = SummaryCallback()5657analysis_interval = 10058analysis_callback = AnalysisCallback(semi_euler, interval = analysis_interval,59save_analysis = true)6061alive_callback = AliveCallback(analysis_interval = analysis_interval)6263save_restart = SaveRestartCallback(interval = 100,64save_final_restart = true)6566save_solution = SaveSolutionCallback(interval = 10,67save_initial_solution = true,68save_final_solution = true,69solution_variables = cons2prim)7071stepsize_callback = StepsizeCallback(cfl = 1.1)7273callbacks = CallbackSet(summary_callback,74analysis_callback,75alive_callback,76save_restart,77save_solution,78stepsize_callback)7980###############################################################################81# run the simulation82sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);83dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback84ode_default_options()..., callback = callbacks);8586println("Number of gravity subcycles: ", semi.gravity_counter.ncalls_since_readout)878889