Path: blob/main/examples/tree_2d_dgsem/elixir_acoustics_gauss_wall.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the acoustic perturbation equations56equations = AcousticPerturbationEquations2D(v_mean_global = (0.5, 0.0), c_mean_global = 1.0,7rho_mean_global = 1.0)89# Create DG solver with polynomial degree = 5 and (local) Lax-Friedrichs/Rusanov flux as surface flux10# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of11# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.12# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.13# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.14# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.15# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the16# `StepsizeCallback` (CFL-Condition) and less diffusion.17solver = DGSEM(polydeg = 5, surface_flux = FluxLaxFriedrichs(max_abs_speed_naive))1819coordinates_min = (-100.0, 0.0) # minimum coordinates (min(x), min(y))20coordinates_max = (100.0, 200.0) # maximum coordinates (max(x), max(y))2122# Create a uniformly refined mesh with periodic boundaries23mesh = TreeMesh(coordinates_min, coordinates_max,24initial_refinement_level = 4,25n_cells_max = 100_000,26periodicity = false)2728"""29initial_condition_gauss_wall(x, t, equations::AcousticPerturbationEquations2D)3031A Gaussian pulse, used in the `gauss_wall` example elixir in combination with32[`boundary_condition_wall`](@ref). Uses the global mean values from `equations`.33"""34function initial_condition_gauss_wall(x, t, equations::AcousticPerturbationEquations2D)35RealT = eltype(x)36v1_prime = 037v2_prime = 038p_prime = exp(-log(convert(RealT, 2)) * (x[1]^2 + (x[2] - 25)^2) / 25)3940prim = SVector(v1_prime, v2_prime, p_prime, global_mean_vars(equations)...)4142return prim2cons(prim, equations)43end44initial_condition = initial_condition_gauss_wall4546# A semidiscretization collects data structures and functions for the spatial discretization47semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,48boundary_conditions = boundary_condition_wall)4950###############################################################################51# ODE solvers, callbacks etc.5253# Create ODE problem with time span from 0.0 to 30.054tspan = (0.0, 30.0)55ode = semidiscretize(semi, tspan)5657# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup58# and resets the timers59summary_callback = SummaryCallback()6061# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results62analysis_callback = AnalysisCallback(semi, interval = 100)6364# The SaveSolutionCallback allows to save the solution to a file in regular intervals65save_solution = SaveSolutionCallback(interval = 100, solution_variables = cons2state)6667# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step68stepsize_callback = StepsizeCallback(cfl = 0.7)6970# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver71callbacks = CallbackSet(summary_callback, analysis_callback, save_solution,72stepsize_callback)7374###############################################################################75# run the simulation7677# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks78sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);79dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback80ode_default_options()..., callback = callbacks)818283