Path: blob/main/examples/structured_1d_dgsem/elixir_linearizedeuler_characteristic_system.jl
2802 views
using OrdinaryDiffEqLowStorageRK1using LinearAlgebra: dot2using Trixi34###############################################################################5# semidiscretization of the linearized Euler equations67rho_0 = 1.08v_0 = 1.09c_0 = 1.010equations = LinearizedEulerEquations1D(rho_0, v_0, c_0)1112solver = DGSEM(polydeg = 3, surface_flux = flux_hll)1314coordinates_min = (0.0,) # minimum coordinate15coordinates_max = (1.0,) # maximum coordinate16cells_per_dimension = (64,)1718mesh = StructuredMesh(cells_per_dimension, coordinates_min, coordinates_max,19periodicity = true)2021# For eigensystem of the linearized Euler equations see e.g.22# https://www.nas.nasa.gov/assets/nas/pdf/ams/2018/introtocfd/Intro2CFD_Lecture1_Pulliam_Euler_WaveEQ.pdf23# Linearized Euler: Eigensystem24lin_euler_eigvals = [v_0 - c_0; v_0; v_0 + c_0]25lin_euler_eigvecs = [-rho_0/c_0 1 rho_0/c_0;261 0 1;27-rho_0*c_0 0 rho_0*c_0]28lin_euler_eigvecs_inv = inv(lin_euler_eigvecs)2930# Trace back characteristics.31# See https://metaphor.ethz.ch/x/2019/hs/401-4671-00L/literature/mishra_hyperbolic_pdes.pdf, p.9532function compute_char_initial_pos(x, t)33return SVector(x[1], x[1], x[1]) .- t * lin_euler_eigvals34end3536function compute_primal_sol(char_vars)37return lin_euler_eigvecs * char_vars38end3940# Initial condition is in principle arbitrary, only periodicity is required41function initial_condition_entropy_wave(x, t, equations::LinearizedEulerEquations1D)42# Parameters43alpha = 1.044beta = 150.045center = 0.54647rho_prime = alpha * exp(-beta * (x[1] - center)^2)48v_prime = 0.049p_prime = 0.05051return SVector(rho_prime, v_prime, p_prime)52end5354function initial_condition_char_vars(x, t, equations::LinearizedEulerEquations1D)55# Trace back characteristics56x_char = compute_char_initial_pos(x, t)5758# Employ periodicity59for p in 1:360while x_char[p] < coordinates_min[1]61x_char[p] += coordinates_max[1] - coordinates_min[1]62end63while x_char[p] > coordinates_max[1]64x_char[p] -= coordinates_max[1] - coordinates_min[1]65end66end6768# Set up characteristic variables69w = zeros(3)70t_0 = 0 # Assumes t_0 = 071for p in 1:372u_char = initial_condition_entropy_wave(x_char[p], t_0, equations)73w[p] = dot(lin_euler_eigvecs_inv[p, :], u_char)74end7576return compute_primal_sol(w)77end7879initial_condition = initial_condition_char_vars8081semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver;82boundary_conditions = boundary_condition_periodic)8384###############################################################################85# ODE solvers, callbacks etc.8687tspan = (0.0, 0.3)88ode = semidiscretize(semi, tspan)8990summary_callback = SummaryCallback()9192analysis_interval = 10093analysis_callback = AnalysisCallback(semi, interval = analysis_interval)9495alive_callback = AliveCallback(analysis_interval = analysis_interval)9697stepsize_callback = StepsizeCallback(cfl = 1.0)9899# Create a CallbackSet to collect all callbacks such that they can be passed to the ODE solver100callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback,101stepsize_callback)102103###############################################################################104# run the simulation105106sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);107dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback108ode_default_options()..., callback = callbacks);109110111