Path: blob/main/examples/p4est_2d_dgsem/elixir_eulergravity_convergence.jl
5586 views
using OrdinaryDiffEqLowStorageRK1using Trixi23initial_condition = initial_condition_eoc_test_coupled_euler_gravity45###############################################################################6# semidiscretization of the compressible Euler equations7gamma = 2.08equations_euler = CompressibleEulerEquations2D(gamma)910polydeg = 311solver_euler = DGSEM(polydeg, FluxHLL(min_max_speed_naive))1213coordinates_min = (0.0, 0.0)14coordinates_max = (2.0, 2.0)1516trees_per_dimension = (1, 1)17mesh = P4estMesh(trees_per_dimension, polydeg = 1,18coordinates_min = coordinates_min, coordinates_max = coordinates_max,19initial_refinement_level = 2,20periodicity = true)2122semi_euler = SemidiscretizationHyperbolic(mesh, equations_euler, initial_condition,23solver_euler;24source_terms = source_terms_eoc_test_coupled_euler_gravity,25boundary_conditions = boundary_condition_periodic)2627###############################################################################28# semidiscretization of the hyperbolic diffusion equations29equations_gravity = HyperbolicDiffusionEquations2D()3031# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of32# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.33# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.34# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.35# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.36# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the37# `StepsizeCallback` (CFL-Condition) and less diffusion.38solver_gravity = DGSEM(polydeg, FluxLaxFriedrichs(max_abs_speed_naive))3940semi_gravity = SemidiscretizationHyperbolic(mesh, equations_gravity, initial_condition,41solver_gravity;42source_terms = source_terms_harmonic,43boundary_conditions = boundary_condition_periodic)4445###############################################################################46# combining both semidiscretizations for Euler + self-gravity47parameters = ParametersEulerGravity(background_density = 2.0, # aka rho048# rho0 is (ab)used to add a "+8π" term to the source terms49# for the manufactured solution50gravitational_constant = 1.0, # aka G51cfl = 1.1,52resid_tol = 1.0e-10,53n_iterations_max = 1000,54timestep_gravity = timestep_gravity_erk52_3Sstar!)5556semi = SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)5758###############################################################################59# ODE solvers, callbacks etc.60tspan = (0.0, 0.5)61ode = semidiscretize(semi, tspan)6263summary_callback = SummaryCallback()6465stepsize_callback = StepsizeCallback(cfl = 0.8)6667save_solution = SaveSolutionCallback(interval = 10,68save_initial_solution = true,69save_final_solution = true,70solution_variables = cons2prim)7172save_restart = SaveRestartCallback(interval = 100,73save_final_restart = true)7475analysis_interval = 10076alive_callback = AliveCallback(analysis_interval = analysis_interval)7778analysis_callback = AnalysisCallback(semi_euler, interval = analysis_interval,79save_analysis = true)8081callbacks = CallbackSet(summary_callback, stepsize_callback,82save_restart, save_solution,83analysis_callback, alive_callback)8485###############################################################################86# run the simulation87sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);88dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback89ode_default_options()..., callback = callbacks);9091println("Number of gravity subcycles: ", semi.gravity_counter.ncalls_since_readout)929394