Path: blob/main/examples/tree_1d_dgsem/elixir_advection_extended.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi2using Plots # For visualization callback34###############################################################################5# semidiscretization of the linear advection equation67advection_velocity = 1.08equations = LinearScalarAdvectionEquation1D(advection_velocity)910initial_condition = initial_condition_convergence_test1112# you can either use a single function to impose the BCs weakly in all13# 1*ndims == 2 directions or you can pass a tuple containing BCs for14# each direction15# Note: "boundary_condition_periodic" indicates that it is a periodic boundary and can be omitted on16# fully periodic domains. Here, however, it is included to allow easy override during testing17boundary_conditions = boundary_condition_periodic1819# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux20solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)2122coordinates_min = -1.0 # minimum coordinate23coordinates_max = 1.0 # maximum coordinate2425# Create a uniformly refined mesh with periodic boundaries26mesh = TreeMesh(coordinates_min, coordinates_max,27initial_refinement_level = 4,28n_cells_max = 30_000, # set maximum capacity of tree data structure29periodicity = true)3031# A semidiscretization collects data structures and functions for the spatial discretization32semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,33boundary_conditions = boundary_conditions)3435###############################################################################36# ODE solvers, callbacks etc.3738# Create ODE problem with time span from 0.0 to 1.039tspan = (0.0, 1.0)40ode = semidiscretize(semi, tspan)4142# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup43# and resets the timers44summary_callback = SummaryCallback()4546# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results47analysis_interval = 10048analysis_callback = AnalysisCallback(semi, interval = analysis_interval,49extra_analysis_integrals = (entropy, energy_total))5051# The AliveCallback prints short status information in regular intervals52alive_callback = AliveCallback(analysis_interval = analysis_interval)5354# The SaveRestartCallback allows to save a file from which a Trixi.jl simulation can be restarted55save_restart = SaveRestartCallback(interval = 100,56save_final_restart = true)5758# The SaveSolutionCallback allows to save the solution to a file in regular intervals59save_solution = SaveSolutionCallback(interval = 100,60save_initial_solution = true,61save_final_solution = true,62solution_variables = cons2prim)6364# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step65stepsize_callback = StepsizeCallback(cfl = 1.6)6667# Enable in-situ visualization with a new plot generated at every time step68visualization = VisualizationCallback(semi; interval = 1)6970# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver71callbacks = CallbackSet(summary_callback,72analysis_callback, alive_callback,73save_restart,74save_solution, visualization,75stepsize_callback)7677###############################################################################78# run the simulation7980# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks81sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);82dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback83ode_default_options()..., callback = callbacks);848586