Path: blob/main/examples/tree_2d_dgsem/elixir_acoustics_gaussian_source.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23# Oscillating Gaussian-shaped source terms4function source_terms_gauss(u, x, t, equations::AcousticPerturbationEquations2D)5RealT = eltype(u)6r = convert(RealT, 0.1)7A = 18f = 2910# Velocity sources11s1 = 012s2 = 013# Pressure source14s3 = exp(-(x[1]^2 + x[2]^2) / (2 * r^2)) * A * sinpi(2 * f * t)1516# Mean sources17s4 = s5 = s6 = s7 = 01819return SVector(s1, s2, s3, s4, s5, s6, s7)20end2122###############################################################################23# semidiscretization of the acoustic perturbation equations2425equations = AcousticPerturbationEquations2D(v_mean_global = (-0.5, 0.25),26c_mean_global = 1.0,27rho_mean_global = 1.0)2829initial_condition = initial_condition_constant3031# Create DG solver with polynomial degree = 3 and (local) Lax-Friedrichs/Rusanov flux as surface flux3233# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of34# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.35# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.36# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.37# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.38# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the39# `StepsizeCallback` (CFL-Condition) and less diffusion.40solver = DGSEM(polydeg = 3, surface_flux = FluxLaxFriedrichs(max_abs_speed_naive))4142coordinates_min = (-3.0, -3.0) # minimum coordinates (min(x), min(y))43coordinates_max = (3.0, 3.0) # maximum coordinates (max(x), max(y))4445# Create a uniformly refined mesh with periodic boundaries46mesh = TreeMesh(coordinates_min, coordinates_max,47initial_refinement_level = 4,48n_cells_max = 30_000) # set maximum capacity of tree data structure4950# A semidiscretization collects data structures and functions for the spatial discretization51semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,52source_terms = source_terms_gauss)5354###############################################################################55# ODE solvers, callbacks etc.5657# Create ODE problem with time span from 0.0 to 2.058tspan = (0.0, 2.0)59ode = semidiscretize(semi, tspan)6061# At the beginning of the main loop, the SummaryCallback prints a summary of the simulation setup62# and resets the timers63summary_callback = SummaryCallback()6465# The AnalysisCallback allows to analyse the solution in regular intervals and prints the results66analysis_callback = AnalysisCallback(semi, interval = 100)6768# The SaveSolutionCallback allows to save the solution to a file in regular intervals69save_solution = SaveSolutionCallback(interval = 100,70solution_variables = cons2prim)7172# The TimeSeriesCallback records the solution at the given points over time73time_series = TimeSeriesCallback(semi, [(0.0, 0.0), (-1.0, 0.5)])7475# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step76stepsize_callback = StepsizeCallback(cfl = 0.5)7778# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver79callbacks = CallbackSet(summary_callback, analysis_callback, save_solution, time_series,80stepsize_callback)8182###############################################################################83# run the simulation8485# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks86sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);87dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback88ode_default_options()..., callback = callbacks);899091