Path: blob/main/examples/p4est_2d_dgsem/elixir_euler_subsonic_constant.jl
2055 views
using OrdinaryDiffEqSSPRK1using Trixi2using LinearAlgebra: norm34###############################################################################5## Semidiscretization of the compressible Euler equations67equations = CompressibleEulerEquations2D(1.4)8polydeg = 39solver = DGSEM(polydeg = polydeg, surface_flux = flux_lax_friedrichs)1011@inline function initial_condition_subsonic(x_, t, equations::CompressibleEulerEquations2D)12rho, v1, v2, p = (0.5313, 0.0, 0.0, 0.4)1314prim = SVector(rho, v1, v2, p)15return prim2cons(prim, equations)16end1718initial_condition = initial_condition_subsonic1920# Calculate the boundary flux from the inner state while using the pressure from the outer state21# when the flow is subsonic (which is always the case in this example).2223# If the naive approach of only using the inner state is used, the errors increase with the24# increase of refinement level, see https://github.com/trixi-framework/Trixi.jl/issues/253025# These errors arise from the corner points in this test.2627# See the reference below for a discussion on inflow/outflow boundary conditions. The subsonic28# outflow boundary conditions are discussed in Section 2.3.29#30# - Jan-Reneé Carlson (2011)31# Inflow/Outflow Boundary Conditions with Application to FUN3D.32# [NASA TM 20110022658](https://ntrs.nasa.gov/citations/20110022658)33@inline function boundary_condition_outflow_general(u_inner,34normal_direction::AbstractVector, x, t,35surface_flux_function,36equations::CompressibleEulerEquations2D)3738# This would be for the general case where we need to check the magnitude of the local Mach number39norm_ = norm(normal_direction)40# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later41normal = normal_direction / norm_4243# Rotate the internal solution state44u_local = Trixi.rotate_to_x(u_inner, normal, equations)4546# Compute the primitive variables47rho_local, v_normal, v_tangent, p_local = cons2prim(u_local, equations)4849# Compute local Mach number50a_local = sqrt(equations.gamma * p_local / rho_local)51Mach_local = abs(v_normal / a_local)52if Mach_local <= 1.0 # The `if` is not needed in this elixir but kept for generality53# In general, `p_local` need not be available from the initial condition54p_local = pressure(initial_condition_subsonic(x, t, equations), equations)55end5657# Create the `u_surface` solution state where the local pressure is possibly set from an external value58prim = SVector(rho_local, v_normal, v_tangent, p_local)59u_boundary = prim2cons(prim, equations)60u_surface = Trixi.rotate_from_x(u_boundary, normal, equations)6162# Compute the flux using the appropriate mixture of internal / external solution states63return flux(u_surface, normal_direction, equations)64end6566boundary_conditions = Dict(:x_neg => boundary_condition_outflow_general,67:x_pos => boundary_condition_outflow_general,68:y_neg => boundary_condition_outflow_general,69:y_pos => boundary_condition_outflow_general)7071coordinates_min = (0.0, 0.0)72coordinates_max = (1.0, 1.0)7374trees_per_dimension = (1, 1)7576mesh = P4estMesh(trees_per_dimension, polydeg = polydeg,77coordinates_min = coordinates_min, coordinates_max = coordinates_max,78initial_refinement_level = 3,79periodicity = false)8081semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition, solver,82boundary_conditions = boundary_conditions)8384###############################################################################85## ODE solvers, callbacks etc.8687tspan = (0.0, 0.25)88ode = semidiscretize(semi, tspan)8990summary_callback = SummaryCallback()9192analysis_interval = 100093analysis_callback = AnalysisCallback(semi, interval = analysis_interval)9495alive_callback = AliveCallback(analysis_interval = 100)9697save_solution = SaveSolutionCallback(interval = 100,98save_initial_solution = true,99save_final_solution = true,100solution_variables = cons2prim)101102stepsize_callback = StepsizeCallback(cfl = 0.5)103104callbacks = CallbackSet(summary_callback,105analysis_callback, alive_callback,106save_solution,107stepsize_callback)108109###############################################################################110## Run the simulation111sol = solve(ode, SSPRK54();112dt = 1.0, # solve needs some value here but it will be overwritten by the stepsize_callback113save_everystep = false, callback = callbacks);114115116