Path: blob/main/examples/tree_1d_dgsem/elixir_advection_perk2.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, ECOS3using Trixi45###############################################################################6# semidiscretization of the linear advection equation78advection_velocity = 1.09equations = LinearScalarAdvectionEquation1D(advection_velocity)1011# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux12solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)1314coordinates_min = -1.0 # minimum coordinate15coordinates_max = 1.0 # maximum coordinate1617# Create a uniformly refined mesh with periodic boundaries18mesh = TreeMesh(coordinates_min, coordinates_max,19initial_refinement_level = 4,20n_cells_max = 30_000, periodicity = true) # set maximum capacity of tree data structure2122# A semidiscretization collects data structures and functions for the spatial discretization23semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition_convergence_test,24solver;25boundary_conditions = boundary_condition_periodic)2627###############################################################################28# ODE solvers, callbacks etc.2930# Create ODE problem with time span from 0.0 to 20.031tspan = (0.0, 20.0)32ode = semidiscretize(semi, tspan)3334# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup35# and resets the timers36summary_callback = SummaryCallback()3738# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results39analysis_interval = 10040analysis_callback = AnalysisCallback(semi, interval = analysis_interval)4142# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step43stepsize_callback = StepsizeCallback(cfl = 2.5)4445alive_callback = AliveCallback(alive_interval = analysis_interval)4647save_solution = SaveSolutionCallback(dt = 0.1,48save_initial_solution = true,49save_final_solution = true,50solution_variables = cons2prim)5152# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver53callbacks = CallbackSet(summary_callback,54alive_callback,55save_solution,56analysis_callback,57stepsize_callback)5859###############################################################################60# run the simulation6162# Construct second order paired explicit Runge-Kutta method with 6 stages for given simulation setup.63# Pass `tspan` to calculate maximum time step allowed for the bisection algorithm used64# in calculating the polynomial coefficients in the ODE algorithm.65ode_algorithm = Trixi.PairedExplicitRK2(6, tspan, semi)6667sol = Trixi.solve(ode, ode_algorithm;68dt = 1.0, # Manual time step value, will be overwritten by the stepsize_callback when it is specified.69ode_default_options()..., callback = callbacks);707172