Path: blob/main/examples/p4est_2d_dgsem/elixir_linearizedeuler_gaussian_source.jl
2055 views
using OrdinaryDiffEqLowStorageRK1using Trixi23# Based on the TreeMesh example `elixir_acoustics_gaussian_source.jl`.4# The acoustic perturbation equations have been replaced with the linearized Euler5# equations and instead of the Cartesian `TreeMesh` a rotated `P4estMesh` is used67# Oscillating Gaussian-shaped source terms8function source_terms_gauss(u, x, t, equations::LinearizedEulerEquations2D)9r = 0.110A = 1.011f = 2.01213# Velocity sources14s2 = 0.015s3 = 0.016# Density and pressure source17s1 = s4 = exp(-(x[1]^2 + x[2]^2) / (2 * r^2)) * A * sin(2 * pi * f * t)1819return SVector(s1, s2, s3, s4)20end2122function initial_condition_zero(x, t, equations::LinearizedEulerEquations2D)23SVector(0.0, 0.0, 0.0, 0.0)24end2526###############################################################################27# semidiscretization of the linearized Euler equations2829# Create a domain that is a 30° rotated version of [-3, 3]^230c = cospi(2 * 30.0 / 360.0)31s = sinpi(2 * 30.0 / 360.0)32rot_mat = Trixi.SMatrix{2, 2}([c -s; s c])33mapping(xi, eta) = rot_mat * SVector(3.0 * xi, 3.0 * eta)3435# Mean density and speed of sound are slightly off from 1.0 to allow proper verification of36# curved LEE implementation using this elixir (some things in the LEE cancel if both are 1.0)37equations = LinearizedEulerEquations2D(v_mean_global = Tuple(rot_mat * SVector(-0.5, 0.25)),38c_mean_global = 1.02, rho_mean_global = 1.01)3940initial_condition = initial_condition_zero4142# Create DG solver with polynomial degree = 3 and upwind flux as surface flux43solver = DGSEM(polydeg = 3, surface_flux = flux_godunov)4445# Create a uniformly refined mesh with periodic boundaries46trees_per_dimension = (4, 4)47mesh = P4estMesh(trees_per_dimension, polydeg = 1,48mapping = mapping,49periodicity = true, initial_refinement_level = 2)5051# A semidiscretization collects data structures and functions for the spatial discretization52semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,53source_terms = source_terms_gauss)5455###############################################################################56# ODE solvers, callbacks etc.5758# Create ODE problem with time span from 0.0 to 2.059tspan = (0.0, 2.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 = 100)7172# The StepsizeCallback handles the re-calculation of the maximum Δt after each time step73stepsize_callback = StepsizeCallback(cfl = 0.5)7475# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver76callbacks = CallbackSet(summary_callback, analysis_callback, save_solution,77stepsize_callback)7879###############################################################################80# run the simulation8182# OrdinaryDiffEq's `solve` method evolves the solution in time and executes the passed callbacks83sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);84dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback85ode_default_options()..., callback = callbacks);868788