Path: blob/main/examples/unstructured_2d_dgsem/elixir_euler_basic.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations56equations = CompressibleEulerEquations2D(1.4)78initial_condition = initial_condition_convergence_test9source_terms = source_terms_convergence_test1011boundary_condition_convergence_test = BoundaryConditionDirichlet(initial_condition)12boundary_conditions = Dict(:Slant => boundary_condition_convergence_test,13:Bezier => boundary_condition_convergence_test,14:Right => boundary_condition_convergence_test,15:Bottom => boundary_condition_convergence_test,16:Top => boundary_condition_convergence_test)1718###############################################################################19# Get the DG approximation space2021# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of22# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.23# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.24# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.25# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.26# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the27# `StepsizeCallback` (CFL-Condition) and less diffusion.28solver = DGSEM(polydeg = 8, surface_flux = FluxLaxFriedrichs(max_abs_speed_naive))2930###############################################################################31# Get the curved quad mesh from a file (downloads the file if not available locally)32mesh_file = Trixi.download("https://gist.githubusercontent.com/andrewwinters5000/52056f1487853fab63b7f4ed7f171c80/raw/9d573387dfdbb8bce2a55db7246f4207663ac07f/mesh_trixi_unstructured_mesh_docs.mesh",33joinpath(@__DIR__, "mesh_trixi_unstructured_mesh_docs.mesh"))3435mesh = UnstructuredMesh2D(mesh_file)3637###############################################################################38# create the semi discretization object3940semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,41source_terms = source_terms,42boundary_conditions = boundary_conditions)4344###############################################################################45# ODE solvers, callbacks etc.4647tspan = (0.0, 3.0)48ode = semidiscretize(semi, tspan)4950summary_callback = SummaryCallback()5152analysis_interval = 10053analysis_callback = AnalysisCallback(semi, interval = analysis_interval)5455alive_callback = AliveCallback(analysis_interval = analysis_interval)5657save_restart = SaveRestartCallback(interval = 50,58save_final_restart = true)5960save_solution = SaveSolutionCallback(interval = 10,61save_initial_solution = true,62save_final_solution = true)6364stepsize_callback = StepsizeCallback(cfl = 0.9)6566callbacks = CallbackSet(summary_callback,67analysis_callback,68alive_callback,69save_restart,70save_solution,71stepsize_callback)7273###############################################################################74# run the simulation7576sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);77dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback78ode_default_options()..., callback = callbacks);798081