Path: blob/main/examples/tree_2d_dgsem/elixir_euler_colliding_flow.jl
5586 views
using OrdinaryDiffEqSSPRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations5gamma = 1.001 # almost isothermal when gamma reaches 16equations = CompressibleEulerEquations2D(gamma)78# This is a hand made colliding flow setup without reference. Features Mach=70 inflow from both9# sides, with relative low temperature, such that pressure keeps relatively small10# Computed with gamma close to 1, to simulate isothermal gas11function initial_condition_colliding_flow_astro(x, t,12equations::CompressibleEulerEquations2D)13# change discontinuity to tanh14# resolution 128^2 elements (refined close to the interface) and polydeg=3 (total of 512^2 DOF)15# domain size is [-64,+64]^216RealT = eltype(x)17@unpack gamma = equations18# the quantities are chosen such, that they are as close as possible to the astro examples19# keep in mind, that in the astro example, the physical units are weird (parsec, mega years, ...)20rho = convert(RealT, 0.0247)21c = convert(RealT, 0.2)22p = c^2 / gamma * rho23vel = convert(RealT, 13.907432274789372)24slope = 125v1 = -vel * tanh(slope * x[1])26# add small initial disturbance to the field, but only close to the interface27if abs(x[1]) < 1028v1 = v1 * (1 + RealT(0.01) * sinpi(x[2]))29end30v2 = 031return prim2cons(SVector(rho, v1, v2, p), equations)32end33initial_condition = initial_condition_colliding_flow_astro3435boundary_conditions = (;36x_neg = BoundaryConditionDirichlet(initial_condition_colliding_flow_astro),37x_pos = BoundaryConditionDirichlet(initial_condition_colliding_flow_astro),38y_neg = boundary_condition_periodic,39y_pos = boundary_condition_periodic)4041# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of42# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.43# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.44# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.45# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.46# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the47# `StepsizeCallback` (CFL-Condition) and less diffusion.48surface_flux = FluxLaxFriedrichs(max_abs_speed_naive) # HLLC needs more shock capturing (alpha_max)49volume_flux = flux_ranocha # works with Chandrashekar flux as well50polydeg = 351basis = LobattoLegendreBasis(polydeg)5253# shock capturing necessary for this tough example, however alpha_max = 0.5 is fine54indicator_sc = IndicatorHennemannGassner(equations, basis,55alpha_max = 0.5,56alpha_min = 0.0001,57alpha_smooth = true,58variable = density_pressure)59volume_integral = VolumeIntegralShockCapturingHG(indicator_sc;60volume_flux_dg = volume_flux,61volume_flux_fv = surface_flux)62solver = DGSEM(basis, surface_flux, volume_integral)6364coordinates_min = (-64.0, -64.0)65coordinates_max = (64.0, 64.0)6667# only refinement in a patch. Needs x=-17/+17 to trigger refinement due to coarse base mesh68refinement_patches = ((type = "box", coordinates_min = (-17, -64),69coordinates_max = (17, 64)),70(type = "box", coordinates_min = (-17, -64),71coordinates_max = (17, 64)),72(type = "box", coordinates_min = (-17, -64),73coordinates_max = (17, 64)),74(type = "box", coordinates_min = (-17, -64),75coordinates_max = (17, 64))76#(type="box", coordinates_min=(-17, -64), coordinates_max=(17, 64)), # very high resolution, takes about 1000s on 2 cores77)78mesh = TreeMesh(coordinates_min, coordinates_max,79initial_refinement_level = 3,80refinement_patches = refinement_patches,81periodicity = (false, true),82n_cells_max = 100_000)83semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;84boundary_conditions = boundary_conditions)8586###############################################################################87# ODE solvers, callbacks etc.8889tspan = (0.0, 25.0)90ode = semidiscretize(semi, tspan)9192summary_callback = SummaryCallback()9394analysis_interval = 100095analysis_callback = AnalysisCallback(semi, interval = analysis_interval)9697alive_callback = AliveCallback(analysis_interval = analysis_interval)9899save_solution = SaveSolutionCallback(interval = 1000,100save_initial_solution = true,101save_final_solution = true,102solution_variables = cons2prim)103104callbacks = CallbackSet(summary_callback,105analysis_callback, alive_callback,106save_solution)107108# positivity limiter necessary for this tough example109stage_limiter! = PositivityPreservingLimiterZhangShu(thresholds = (5.0e-6, 5.0e-6),110variables = (Trixi.density, pressure))111112###############################################################################113# run the simulation114# use adaptive time stepping based on error estimates, time step roughly dt = 5e-3115sol = solve(ode, SSPRK43(stage_limiter!);116ode_default_options()..., callback = callbacks);117118119