Path: blob/main/examples/structured_1d_dgsem/elixir_burgers_perk3.jl
5586 views
# Convex and ECOS are imported because they are used for finding the optimal time step and optimal1# monomial coefficients in the stability polynomial of PERK time integrators.2using Convex, ECOS34# NLsolve is imported to solve the system of nonlinear equations to find the coefficients5# in the Butcher tableau in the third order PERK time integrator.6using NLsolve78using Trixi910###############################################################################11# semidiscretization of the (inviscid) Burgers equation1213equations = InviscidBurgersEquation1D()1415initial_condition = initial_condition_convergence_test1617# Create DG solver with polynomial degree = 4 and (local) Lax-Friedrichs/Rusanov flux as surface flux18solver = DGSEM(polydeg = 4, surface_flux = flux_lax_friedrichs)1920coordinates_min = (0.0,) # minimum coordinate21coordinates_max = (1.0,) # maximum coordinate22cells_per_dimension = (64,)2324mesh = StructuredMesh(cells_per_dimension, coordinates_min, coordinates_max,25periodicity = true)2627semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;28source_terms = source_terms_convergence_test,29boundary_conditions = boundary_condition_periodic)3031###############################################################################32# ODE solvers, callbacks etc.3334tspan = (0.0, 2.0)35ode = semidiscretize(semi, tspan)3637summary_callback = SummaryCallback()3839analysis_interval = 20040analysis_callback = AnalysisCallback(semi, interval = analysis_interval)4142alive_callback = AliveCallback(analysis_interval = analysis_interval)4344save_solution = SaveSolutionCallback(dt = 0.1,45save_initial_solution = true,46save_final_solution = true,47solution_variables = cons2prim)4849# Construct third order paired explicit Runge-Kutta method with 8 stages for given simulation setup.50# Pass `tspan` to calculate maximum time step allowed for the bisection algorithm used51# in calculating the polynomial coefficients in the ODE algorithm.52ode_algorithm = Trixi.PairedExplicitRK3(8, tspan, semi)5354cfl_number = Trixi.calculate_cfl(ode_algorithm, ode)55# For non-linear problems, the CFL number should be reduced by a safety factor56# as the spectrum changes (in general) over the course of a simulation57stepsize_callback = StepsizeCallback(cfl = 0.85 * cfl_number)5859callbacks = CallbackSet(summary_callback,60analysis_callback, alive_callback, save_solution,61stepsize_callback)6263###############################################################################64# run the simulation65sol = Trixi.solve(ode, ode_algorithm;66dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback67ode_default_options()..., callback = callbacks);686970