Path: blob/main/examples/unstructured_2d_dgsem/elixir_acoustics_gauss_wall.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the acoustic perturbation equations56equations = AcousticPerturbationEquations2D(v_mean_global = (0.0, -0.5),7c_mean_global = 1.0,8rho_mean_global = 1.0)910# Create DG solver with polynomial degree = 4 and (local) Lax-Friedrichs/Rusanov flux1112# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of13# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.14# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.15# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.16# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.17# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the18# `StepsizeCallback` (CFL-Condition) and less diffusion.19solver = DGSEM(polydeg = 4, surface_flux = FluxLaxFriedrichs(max_abs_speed_naive))2021# Create unstructured quadrilateral mesh from a file22mesh_file = Trixi.download("https://gist.githubusercontent.com/andrewwinters5000/3c79baad6b4d73bb26ec6420b5d16f45/raw/22aefc4ec2107cf0bffc40e81dfbc52240c625b1/mesh_five_circles_in_circle.mesh",23joinpath(@__DIR__, "mesh_five_circles_in_circle.mesh"))2425mesh = UnstructuredMesh2D(mesh_file)2627"""28initial_condition_gauss_wall(x, t, equations::AcousticPerturbationEquations2D)2930A Gaussian pulse, used in the `gauss_wall` example elixir in combination with31[`boundary_condition_wall`](@ref). Uses the global mean values from `equations`.32"""33function initial_condition_gauss_wall(x, t, equations::AcousticPerturbationEquations2D)34v1_prime = 0.035v2_prime = 0.036p_prime = exp(-log(2) * (x[1]^2 + (x[2] - 25)^2) / 25)3738prim = SVector(v1_prime, v2_prime, p_prime, global_mean_vars(equations)...)3940return prim2cons(prim, equations)41end42initial_condition = initial_condition_gauss_wall4344boundary_conditions = Dict(:OuterCircle => boundary_condition_slip_wall,45:InnerCircle1 => boundary_condition_slip_wall,46:InnerCircle2 => boundary_condition_slip_wall,47:InnerCircle3 => boundary_condition_slip_wall,48:InnerCircle4 => boundary_condition_slip_wall,49:InnerCircle5 => boundary_condition_slip_wall)5051# A semidiscretization collects data structures and functions for the spatial discretization52semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,53boundary_conditions = boundary_conditions)5455###############################################################################56# ODE solvers, callbacks etc.5758# Create ODE problem with time span from 0.0 to 300.059tspan = (0.0, 300.0)60ode = semidiscretize(semi, tspan)6162# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup63# and resets the timers64summary_callback = SummaryCallback()6566# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results67analysis_callback = AnalysisCallback(semi, interval = 100)6869# The SaveSolutionCallback allows to save the solution to a file in regular intervals70save_solution = SaveSolutionCallback(interval = 50, solution_variables = cons2state)7172# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver73callbacks = CallbackSet(summary_callback, analysis_callback, save_solution)7475###############################################################################76# run the simulation7778# use a Runge-Kutta method with automatic (error based) time step size control79sol = solve(ode, RDPK3SpFSAL49(); abstol = 1.0e-6, reltol = 1.0e-6,80ode_default_options()..., callback = callbacks);818283