Path: blob/main/examples/tree_2d_dgsem/elixir_acoustics_gauss.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the acoustic perturbation equations56v_mean_global = (0.25, 0.25)7c_mean_global = 1.08rho_mean_global = 1.09equations = AcousticPerturbationEquations2D(v_mean_global, c_mean_global, rho_mean_global)1011# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux1213# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of14# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.15# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.16# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.17# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.18# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the19# `StepsizeCallback` (CFL-Condition) and less diffusion.20solver = DGSEM(polydeg = 3, surface_flux = FluxLaxFriedrichs(max_abs_speed_naive))2122coordinates_min = (-1.0, -1.0) # minimum coordinates (min(x), min(y))23coordinates_max = (1.0, 1.0) # maximum coordinates (max(x), max(y))2425# 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 structure2930# A semidiscretization collects data structures and functions for the spatial discretization31semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition_gauss, solver)3233###############################################################################34# ODE solvers, callbacks etc.3536# Create ODE problem with time span from 0.0 to 1.037tspan = (0.0, 1.0)38ode = semidiscretize(semi, tspan)3940# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup41# and resets the timers42summary_callback = SummaryCallback()4344# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results45analysis_callback = AnalysisCallback(semi, interval = 100)4647# The SaveSolutionCallback allows to save the solution to a file in regular intervals48save_solution = SaveSolutionCallback(interval = 100,49solution_variables = cons2prim)5051# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step52stepsize_callback = StepsizeCallback(cfl = 0.5)5354# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver55callbacks = CallbackSet(summary_callback, analysis_callback, save_solution,56stepsize_callback)5758###############################################################################59# run the simulation6061# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks62sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);63dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback64ode_default_options()..., callback = callbacks);656667