Path: blob/main/examples/structured_1d_dgsem/elixir_burgers_perk3.jl
2055 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)2526semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,27source_terms = source_terms_convergence_test)2829###############################################################################30# ODE solvers, callbacks etc.3132tspan = (0.0, 2.0)33ode = semidiscretize(semi, tspan)3435summary_callback = SummaryCallback()3637analysis_interval = 20038analysis_callback = AnalysisCallback(semi, interval = analysis_interval)3940alive_callback = AliveCallback(analysis_interval = analysis_interval)4142save_solution = SaveSolutionCallback(dt = 0.1,43save_initial_solution = true,44save_final_solution = true,45solution_variables = cons2prim)4647# Construct third order paired explicit Runge-Kutta method with 8 stages for given simulation setup.48# Pass `tspan` to calculate maximum time step allowed for the bisection algorithm used49# in calculating the polynomial coefficients in the ODE algorithm.50ode_algorithm = Trixi.PairedExplicitRK3(8, tspan, semi)5152cfl_number = Trixi.calculate_cfl(ode_algorithm, ode)53# For non-linear problems, the CFL number should be reduced by a safety factor54# as the spectrum changes (in general) over the course of a simulation55stepsize_callback = StepsizeCallback(cfl = 0.85 * cfl_number)5657callbacks = CallbackSet(summary_callback,58analysis_callback, alive_callback, save_solution,59stepsize_callback)6061###############################################################################62# run the simulation63sol = Trixi.solve(ode, ode_algorithm;64dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback65ode_default_options()..., callback = callbacks);666768