Path: blob/main/examples/tree_2d_dgsem/elixir_euler_blob_amr.jl
5585 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations5gamma = 5 / 36equations = CompressibleEulerEquations2D(gamma)78"""9initial_condition_blob(x, t, equations::CompressibleEulerEquations2D)1011The blob test case taken from12- Agertz et al. (2006)13Fundamental differences between SPH and grid methods14[arXiv: astro-ph/0610051](https://arxiv.org/abs/astro-ph/0610051)15"""16function initial_condition_blob(x, t, equations::CompressibleEulerEquations2D)17# blob test case, see Agertz et al. https://arxiv.org/pdf/astro-ph/0610051.pdf18# other reference: https://doi.org/10.1111/j.1365-2966.2007.12183.x19# change discontinuity to tanh20# typical domain is rectangular, we change it to a square21# resolution 128^2, 256^222# domain size is [-20.0,20.0]^223# gamma = 5/3 for this test case24RealT = eltype(x)25R = 1 # radius of the blob26# background density27dens0 = 128Chi = convert(RealT, 10) # density contrast29# reference time of characteristic growth of KH instability equal to 1.030tau_kh = 131tau_cr = tau_kh / convert(RealT, 1.6) # crushing time32# determine background velocity33velx0 = 2 * R * sqrt(Chi) / tau_cr34vely0 = 035Ma0 = convert(RealT, 2.7) # background flow Mach number Ma=v/c36c = velx0 / Ma0 # sound speed37# use perfect gas assumption to compute background pressure via the sound speed c^2 = gamma * pressure/density38p0 = c * c * dens0 / equations.gamma39# initial center of the blob40inicenter = [-15, 0]41x_rel = x - inicenter42r = sqrt(x_rel[1]^2 + x_rel[2]^2)43# steepness of the tanh transition zone44slope = 245# density blob46dens = dens0 +47(Chi - 1) * 0.5f0 * (1 + (tanh(slope * (r + R)) - (tanh(slope * (r - R)) + 1)))48# velocity blob is zero49velx = velx0 -50velx0 * 0.5f0 * (1 + (tanh(slope * (r + R)) - (tanh(slope * (r - R)) + 1)))51return prim2cons(SVector(dens, velx, vely0, p0), equations)52end53initial_condition = initial_condition_blob5455# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of56# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.57# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.58# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.59# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.60# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the61# `StepsizeCallback` (CFL-Condition) and less diffusion.62surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)63volume_flux = flux_ranocha64basis = LobattoLegendreBasis(4)6566indicator_sc = IndicatorHennemannGassner(equations, basis,67alpha_max = 0.4,68alpha_min = 0.0001,69alpha_smooth = true,70variable = pressure)7172volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;73volume_flux_dg = volume_flux,74volume_flux_fv = surface_flux)7576solver = DGSEM(basis, surface_flux, volume_integral)7778coordinates_min = (-20.0, -20.0)79coordinates_max = (20.0, 20.0)8081mesh = TreeMesh(coordinates_min, coordinates_max,82initial_refinement_level = 6,83n_cells_max = 100_000, periodicity = true)8485semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;86boundary_conditions = boundary_condition_periodic)8788###############################################################################89# ODE solvers, callbacks etc.9091tspan = (0.0, 8.0)92ode = semidiscretize(semi, tspan)9394summary_callback = SummaryCallback()9596analysis_interval = 10097analysis_callback = AnalysisCallback(semi, interval = analysis_interval)9899alive_callback = AliveCallback(analysis_interval = analysis_interval)100101save_solution = SaveSolutionCallback(interval = 100,102save_initial_solution = true,103save_final_solution = true,104solution_variables = cons2prim)105106amr_indicator = IndicatorHennemannGassner(semi,107alpha_max = 1.0,108alpha_min = 0.0001,109alpha_smooth = false,110variable = Trixi.density)111amr_controller = ControllerThreeLevelCombined(semi, amr_indicator, indicator_sc,112base_level = 4,113med_level = 0, med_threshold = 0.0003, # med_level = current level114max_level = 7, max_threshold = 0.003,115max_threshold_secondary = indicator_sc.alpha_max)116amr_callback = AMRCallback(semi, amr_controller,117interval = 1,118adapt_initial_condition = true,119adapt_initial_condition_only_refine = true)120121stepsize_callback = StepsizeCallback(cfl = 0.25)122123callbacks = CallbackSet(summary_callback,124analysis_callback, alive_callback,125save_solution,126amr_callback, stepsize_callback)127128###############################################################################129# run the simulation130131sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);132dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback133ode_default_options()..., callback = callbacks);134135136