Path: blob/main/examples/t8code_2d_dgsem/elixir_eulergravity_convergence.jl
5585 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)1718mesh = T8codeMesh(trees_per_dimension, polydeg = 1,19coordinates_min = coordinates_min, coordinates_max = coordinates_max,20initial_refinement_level = 2,21periodicity = true)2223semi_euler = SemidiscretizationHyperbolic(mesh, equations_euler, initial_condition,24solver_euler;25source_terms = source_terms_eoc_test_coupled_euler_gravity,26boundary_conditions = boundary_condition_periodic)2728###############################################################################29# semidiscretization of the hyperbolic diffusion equations30equations_gravity = HyperbolicDiffusionEquations2D()3132# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of33# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.34# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.35# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.36# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.37# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the38# `StepsizeCallback` (CFL-Condition) and less diffusion.39solver_gravity = DGSEM(polydeg, FluxLaxFriedrichs(max_abs_speed_naive))4041semi_gravity = SemidiscretizationHyperbolic(mesh, equations_gravity, initial_condition,42solver_gravity;43source_terms = source_terms_harmonic,44boundary_conditions = boundary_condition_periodic)4546###############################################################################47# combining both semidiscretizations for Euler + self-gravity48parameters = ParametersEulerGravity(background_density = 2.0, # aka rho049# rho0 is (ab)used to add a "+8π" term to the source terms50# for the manufactured solution51gravitational_constant = 1.0, # aka G52cfl = 1.1,53resid_tol = 1.0e-10,54n_iterations_max = 1000,55timestep_gravity = timestep_gravity_erk52_3Sstar!)5657semi = SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)5859###############################################################################60# ODE solvers, callbacks etc.61tspan = (0.0, 0.5)62ode = semidiscretize(semi, tspan)6364summary_callback = SummaryCallback()6566stepsize_callback = StepsizeCallback(cfl = 0.8)6768save_solution = SaveSolutionCallback(interval = 10,69save_initial_solution = true,70save_final_solution = true,71solution_variables = cons2prim)7273save_restart = SaveRestartCallback(interval = 100,74save_final_restart = true)7576analysis_interval = 10077alive_callback = AliveCallback(analysis_interval = analysis_interval)7879analysis_callback = AnalysisCallback(semi_euler, interval = analysis_interval,80save_analysis = true)8182callbacks = CallbackSet(summary_callback, stepsize_callback,83save_restart, save_solution,84analysis_callback, alive_callback)8586###############################################################################87# run the simulation88sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);89dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback90ode_default_options()..., callback = callbacks);9192println("Number of gravity subcycles: ", semi.gravity_counter.ncalls_since_readout)939495