Path: blob/main/examples/tree_1d_fdsbp/elixir_euler_convergence.jl
5586 views
# !!! warning "Experimental implementation (upwind SBP)"1# This is an experimental feature and may change in future releases.23using OrdinaryDiffEqSSPRK4using Trixi56###############################################################################7# semidiscretization of the compressible Euler equations89equations = CompressibleEulerEquations1D(1.4)1011initial_condition = initial_condition_convergence_test1213D_upw = upwind_operators(SummationByPartsOperators.Mattsson2017,14derivative_order = 1,15accuracy_order = 4,16xmin = -1.0, xmax = 1.0,17N = 32)18flux_splitting = splitting_steger_warming19solver = FDSBP(D_upw,20surface_integral = SurfaceIntegralUpwind(flux_splitting),21volume_integral = VolumeIntegralUpwind(flux_splitting))2223coordinates_min = 0.024coordinates_max = 2.025mesh = TreeMesh(coordinates_min, coordinates_max,26initial_refinement_level = 1,27n_cells_max = 10_000, periodicity = true)2829semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;30source_terms = source_terms_convergence_test,31boundary_conditions = boundary_condition_periodic)3233###############################################################################34# ODE solvers, callbacks etc.3536tspan = (0.0, 2.0)37ode = semidiscretize(semi, tspan)3839summary_callback = SummaryCallback()4041analysis_interval = 10042analysis_callback = AnalysisCallback(semi, interval = analysis_interval,43extra_analysis_errors = (:l2_error_primitive,44:linf_error_primitive))4546alive_callback = AliveCallback(analysis_interval = analysis_interval)4748save_solution = SaveSolutionCallback(interval = 100,49save_initial_solution = true,50save_final_solution = true,51solution_variables = cons2prim)5253callbacks = CallbackSet(summary_callback,54analysis_callback, alive_callback,55save_solution)5657###############################################################################58# run the simulation5960sol = solve(ode, SSPRK43(); abstol = 1.0e-6, reltol = 1.0e-6,61ode_default_options()..., callback = callbacks)626364