Path: blob/main/examples/tree_2d_dgsem/elixir_advection_extended.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the linear advection equation56advection_velocity = (0.2, -0.7)7equations = LinearScalarAdvectionEquation2D(advection_velocity)89initial_condition = initial_condition_convergence_test1011# you can either use a single function to impose the BCs weakly in all12# 1*ndims == 2 directions or you can pass a tuple containing BCs for13# each direction14# Note: "boundary_condition_periodic" indicates that it is a periodic boundary and can be omitted on15# fully periodic domains. Here, however, it is included to allow easy override during testing16boundary_conditions = boundary_condition_periodic1718# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux19solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)2021coordinates_min = (-1.0, -1.0) # minimum coordinates (min(x), min(y))22coordinates_max = (1.0, 1.0) # maximum coordinates (max(x), max(y))2324# Create a uniformly refined mesh with periodic boundaries25mesh = TreeMesh(coordinates_min, coordinates_max,26initial_refinement_level = 4,27n_cells_max = 30_000, # set maximum capacity of tree data structure28periodicity = true)2930# A semidiscretization collects data structures and functions for the spatial discretization31semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,32boundary_conditions = boundary_conditions)3334###############################################################################35# ODE solvers, callbacks etc.3637# Create ODE problem with time span from 0.0 to 1.038tspan = (0.0, 1.0)39ode = semidiscretize(semi, tspan)4041# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup42# and resets the timers43summary_callback = SummaryCallback()4445# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results46analysis_interval = 10047analysis_callback = AnalysisCallback(semi, interval = analysis_interval,48extra_analysis_integrals = (entropy, energy_total))4950# The AliveCallback prints short status information in regular intervals51alive_callback = AliveCallback(analysis_interval = analysis_interval)5253# The SaveRestartCallback allows to save a file from which a Trixi.jl simulation can be restarted54save_restart = SaveRestartCallback(interval = 40,55save_final_restart = true)5657# The SaveSolutionCallback allows to save the solution to a file in regular intervals58save_solution = SaveSolutionCallback(interval = 100,59save_initial_solution = true,60save_final_solution = true,61solution_variables = cons2prim)6263# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step64stepsize_callback = StepsizeCallback(cfl = 1.6)6566# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver67callbacks = CallbackSet(summary_callback,68analysis_callback, alive_callback,69save_restart, save_solution,70stepsize_callback)7172###############################################################################73# run the simulation7475# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks76alg = CarpenterKennedy2N54(williamson_condition = false)77sol = solve(ode, alg;78dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback79callback = callbacks,80ode_default_options()...);818283