Path: blob/main/examples/tree_3d_dgsem/elixir_eulergravity_convergence.jl
5586 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations5gamma = 2.06equations_euler = CompressibleEulerEquations3D(gamma)78initial_condition = initial_condition_eoc_test_coupled_euler_gravity910polydeg = 311solver_euler = DGSEM(polydeg, FluxHLL(min_max_speed_naive))1213coordinates_min = (0.0, 0.0, 0.0)14coordinates_max = (2.0, 2.0, 2.0)15mesh = TreeMesh(coordinates_min, coordinates_max,16initial_refinement_level = 2,17n_cells_max = 10_000, periodicity = true)1819semi_euler = SemidiscretizationHyperbolic(mesh, equations_euler, initial_condition,20solver_euler;21boundary_conditions = boundary_condition_periodic,22source_terms = source_terms_eoc_test_coupled_euler_gravity)2324###############################################################################25# semidiscretization of the hyperbolic diffusion equations26equations_gravity = HyperbolicDiffusionEquations3D()2728# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of29# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.30# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.31# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.32# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.33# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the34# `StepsizeCallback` (CFL-Condition) and less diffusion.35solver_gravity = DGSEM(polydeg, FluxLaxFriedrichs(max_abs_speed_naive))3637semi_gravity = SemidiscretizationHyperbolic(mesh, equations_gravity, initial_condition,38solver_gravity;39source_terms = source_terms_harmonic,40boundary_conditions = boundary_condition_periodic)4142###############################################################################43# combining both semidiscretizations for Euler + self-gravity44parameters = ParametersEulerGravity(background_density = 2.0, # aka rho045gravitational_constant = 1.0, # aka G46cfl = 1.5,47resid_tol = 1.0e-10,48n_iterations_max = 1000,49timestep_gravity = timestep_gravity_erk52_3Sstar!)5051semi = SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)5253###############################################################################54# ODE solvers, callbacks etc.55tspan = (0.0, 0.5)56ode = semidiscretize(semi, tspan)5758summary_callback = SummaryCallback()5960analysis_interval = 10061analysis_callback = AnalysisCallback(semi_euler, interval = analysis_interval,62save_analysis = true)6364alive_callback = AliveCallback(analysis_interval = analysis_interval)6566save_restart = SaveRestartCallback(interval = 100,67save_final_restart = true)6869save_solution = SaveSolutionCallback(interval = 10,70save_initial_solution = true,71save_final_solution = true,72solution_variables = cons2prim)7374stepsize_callback = StepsizeCallback(cfl = 1.1)7576callbacks = CallbackSet(summary_callback,77analysis_callback,78alive_callback,79save_restart,80save_solution,81stepsize_callback)8283###############################################################################84# run the simulation85sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);86dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback87ode_default_options()..., callback = callbacks);8889println("Number of gravity subcycles: ", semi.gravity_counter.ncalls_since_readout)909192