Path: blob/main/examples/p4est_2d_dgsem/elixir_eulergravity_convergence.jl
2055 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)2021semi_euler = SemidiscretizationHyperbolic(mesh, equations_euler, initial_condition,22solver_euler,23source_terms = source_terms_eoc_test_coupled_euler_gravity)2425###############################################################################26# semidiscretization of the hyperbolic diffusion equations27equations_gravity = HyperbolicDiffusionEquations2D()2829# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of30# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.31# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.32# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.33# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.34# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the35# `StepsizeCallback` (CFL-Condition) and less diffusion.36solver_gravity = DGSEM(polydeg, FluxLaxFriedrichs(max_abs_speed_naive))3738semi_gravity = SemidiscretizationHyperbolic(mesh, equations_gravity, initial_condition,39solver_gravity,40source_terms = source_terms_harmonic)4142###############################################################################43# combining both semidiscretizations for Euler + self-gravity44parameters = ParametersEulerGravity(background_density = 2.0, # aka rho045# rho0 is (ab)used to add a "+8π" term to the source terms46# for the manufactured solution47gravitational_constant = 1.0, # aka G48cfl = 1.1,49resid_tol = 1.0e-10,50n_iterations_max = 1000,51timestep_gravity = timestep_gravity_erk52_3Sstar!)5253semi = SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)5455###############################################################################56# ODE solvers, callbacks etc.57tspan = (0.0, 0.5)58ode = semidiscretize(semi, tspan)5960summary_callback = SummaryCallback()6162stepsize_callback = StepsizeCallback(cfl = 0.8)6364save_solution = SaveSolutionCallback(interval = 10,65save_initial_solution = true,66save_final_solution = true,67solution_variables = cons2prim)6869save_restart = SaveRestartCallback(interval = 100,70save_final_restart = true)7172analysis_interval = 10073alive_callback = AliveCallback(analysis_interval = analysis_interval)7475analysis_callback = AnalysisCallback(semi_euler, interval = analysis_interval,76save_analysis = true)7778callbacks = CallbackSet(summary_callback, stepsize_callback,79save_restart, save_solution,80analysis_callback, alive_callback)8182###############################################################################83# run the simulation84sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);85dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback86ode_default_options()..., callback = callbacks);8788println("Number of gravity subcycles: ", semi.gravity_counter.ncalls_since_readout)899091