Path: blob/main/examples/tree_3d_dgsem/elixir_eulergravity_convergence.jl
2055 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)1819semi_euler = SemidiscretizationHyperbolic(mesh, equations_euler, initial_condition,20solver_euler,21source_terms = source_terms_eoc_test_coupled_euler_gravity)2223###############################################################################24# semidiscretization of the hyperbolic diffusion equations25equations_gravity = HyperbolicDiffusionEquations3D()2627# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of28# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.29# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.30# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.31# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.32# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the33# `StepsizeCallback` (CFL-Condition) and less diffusion.34solver_gravity = DGSEM(polydeg, FluxLaxFriedrichs(max_abs_speed_naive))3536semi_gravity = SemidiscretizationHyperbolic(mesh, equations_gravity, initial_condition,37solver_gravity,38source_terms = source_terms_harmonic)3940###############################################################################41# combining both semidiscretizations for Euler + self-gravity42parameters = ParametersEulerGravity(background_density = 2.0, # aka rho043gravitational_constant = 1.0, # aka G44cfl = 1.5,45resid_tol = 1.0e-10,46n_iterations_max = 1000,47timestep_gravity = timestep_gravity_erk52_3Sstar!)4849semi = SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)5051###############################################################################52# ODE solvers, callbacks etc.53tspan = (0.0, 0.5)54ode = semidiscretize(semi, tspan)5556summary_callback = SummaryCallback()5758analysis_interval = 10059analysis_callback = AnalysisCallback(semi_euler, interval = analysis_interval,60save_analysis = true)6162alive_callback = AliveCallback(analysis_interval = analysis_interval)6364save_restart = SaveRestartCallback(interval = 100,65save_final_restart = true)6667save_solution = SaveSolutionCallback(interval = 10,68save_initial_solution = true,69save_final_solution = true,70solution_variables = cons2prim)7172stepsize_callback = StepsizeCallback(cfl = 1.1)7374callbacks = CallbackSet(summary_callback,75analysis_callback,76alive_callback,77save_restart,78save_solution,79stepsize_callback)8081###############################################################################82# run the simulation83sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);84dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback85ode_default_options()..., callback = callbacks);8687println("Number of gravity subcycles: ", semi.gravity_counter.ncalls_since_readout)888990