Path: blob/main/examples/p4est_2d_dgsem/elixir_navierstokes_NACA0012airfoil_mach08.jl
2835 views
using OrdinaryDiffEqLowStorageRK1using Trixi23###############################################################################4# semidiscretization of the compressible Euler equations56# Laminar transonic flow around a NACA0012 airfoil.78# This test is taken from the paper of Swanson and Langer. The values for the drag and lift coefficients9# from Case 5 in Table 3 are used to validate the scheme and computation of surface forces.1011# References:12# - Roy Charles Swanson, Stefan Langer (2016)13# Structured and Unstructured Grid Methods (2016)14# [https://ntrs.nasa.gov/citations/20160003623] (https://ntrs.nasa.gov/citations/20160003623)15# - Deep Ray, Praveen Chandrashekar (2017)16# An entropy stable finite volume scheme for the17# two dimensional Navier–Stokes equations on triangular grids18# [DOI:10.1016/j.amc.2017.07.020](https://doi.org/10.1016/j.amc.2017.07.020)1920equations = CompressibleEulerEquations2D(1.4)2122prandtl_number() = 0.7223mu() = 0.003195997496870108824equations_parabolic = CompressibleNavierStokesDiffusion2D(equations, mu = mu(),25Prandtl = prandtl_number())2627rho_inf() = 1.028p_inf() = 2.8529aoa() = 10.0 * pi / 180.0 # 10 degree angle of attack30l_inf() = 1.031mach_inf() = 0.832u_inf(equations) = mach_inf() * sqrt(equations.gamma * p_inf() / rho_inf())33@inline function initial_condition_mach08_flow(x, t,34equations::CompressibleEulerEquations2D)35# set the freestream flow parameters36gamma = equations.gamma37u_inf = mach_inf() * sqrt(gamma * p_inf() / rho_inf())3839v1 = u_inf * cos(aoa())40v2 = u_inf * sin(aoa())4142prim = SVector(rho_inf(), v1, v2, p_inf())43return prim2cons(prim, equations)44end4546initial_condition = initial_condition_mach08_flow4748# Up to version 0.13.0, `max_abs_speed_naive` was used as the default wave speed estimate of49# `const flux_lax_friedrichs = FluxLaxFriedrichs(), i.e., `FluxLaxFriedrichs(max_abs_speed = max_abs_speed_naive)`.50# In the `StepsizeCallback`, though, the less diffusive `max_abs_speeds` is employed which is consistent with `max_abs_speed`.51# Thus, we exchanged in PR#2458 the default wave speed used in the LLF flux to `max_abs_speed`.52# To ensure that every example still runs we specify explicitly `FluxLaxFriedrichs(max_abs_speed_naive)`.53# We remark, however, that the now default `max_abs_speed` is in general recommended due to compliance with the54# `StepsizeCallback` (CFL-Condition) and less diffusion.55surface_flux = FluxLaxFriedrichs(max_abs_speed_naive)5657polydeg = 358solver = DGSEM(polydeg = polydeg, surface_flux = surface_flux)5960mesh_file = Trixi.download("https://gist.githubusercontent.com/Arpit-Babbar/339662b4b46164a016e35c81c66383bb/raw/8bf94f5b426ba907ace87405cfcc1dcc2ef7cbda/NACA0012.inp",61joinpath(@__DIR__, "NACA0012.inp"))6263mesh = P4estMesh{2}(mesh_file, initial_refinement_level = 1)6465# The boundary values across outer boundary are constant but subsonic, so we cannot compute the66# boundary flux from the external information alone. Thus, we use the numerical flux to distinguish67# between inflow and outflow characteristics68@inline function boundary_condition_subsonic_constant(u_inner,69normal_direction::AbstractVector, x,70t,71surface_flux_function,72equations::CompressibleEulerEquations2D)73u_boundary = initial_condition_mach08_flow(x, t, equations)7475return flux_hll(u_inner, u_boundary, normal_direction, equations)76end7778boundary_conditions = (; Left = boundary_condition_subsonic_constant,79Right = boundary_condition_subsonic_constant,80Top = boundary_condition_subsonic_constant,81Bottom = boundary_condition_subsonic_constant,82AirfoilBottom = boundary_condition_slip_wall,83AirfoilTop = boundary_condition_slip_wall)8485velocity_airfoil = NoSlip((x, t, equations_parabolic) -> SVector(0.0, 0.0))8687heat_airfoil = Adiabatic((x, t, equations_parabolic) -> 0.0)8889boundary_conditions_airfoil = BoundaryConditionNavierStokesWall(velocity_airfoil,90heat_airfoil)9192function velocities_initial_condition_mach08_flow(x, t,93equations::CompressibleEulerEquations2D)94u_cons = initial_condition_mach08_flow(x, t, equations)95return SVector(u_cons[2] / u_cons[1], u_cons[3] / u_cons[1])96end9798velocity_bc_square = NoSlip((x, t, equations_parabolic) -> velocities_initial_condition_mach08_flow(x,99t,100equations_parabolic.equations_hyperbolic))101102heat_bc_square = Adiabatic((x, t, equations_parabolic) -> 0.0)103boundary_condition_square = BoundaryConditionNavierStokesWall(velocity_bc_square,104heat_bc_square)105106boundary_conditions_parabolic = (; Left = boundary_condition_square,107Right = boundary_condition_square,108Top = boundary_condition_square,109Bottom = boundary_condition_square,110AirfoilBottom = boundary_conditions_airfoil,111AirfoilTop = boundary_conditions_airfoil)112113semi = SemidiscretizationHyperbolicParabolic(mesh, (equations, equations_parabolic),114initial_condition, solver;115boundary_conditions = (boundary_conditions,116boundary_conditions_parabolic))117118###############################################################################119# ODE solvers120121# Run for a long time to reach a state where forces stabilize up to 3 digits122tspan = (0.0, 10.0)123ode = semidiscretize(semi, tspan)124125# Callbacks126127summary_callback = SummaryCallback()128129analysis_interval = 2000130131force_boundary_names = (:AirfoilBottom, :AirfoilTop)132drag_coefficient = AnalysisSurfaceIntegral(force_boundary_names,133DragCoefficientPressure2D(aoa(), rho_inf(),134u_inf(equations),135l_inf()))136137lift_coefficient = AnalysisSurfaceIntegral(force_boundary_names,138LiftCoefficientPressure2D(aoa(), rho_inf(),139u_inf(equations),140l_inf()))141142drag_coefficient_shear_force = AnalysisSurfaceIntegral(force_boundary_names,143DragCoefficientShearStress2D(aoa(),144rho_inf(),145u_inf(equations),146l_inf()))147148lift_coefficient_shear_force = AnalysisSurfaceIntegral(force_boundary_names,149LiftCoefficientShearStress2D(aoa(),150rho_inf(),151u_inf(equations),152l_inf()))153154analysis_callback = AnalysisCallback(semi, interval = analysis_interval,155output_directory = "out",156save_analysis = true,157analysis_errors = Symbol[],158analysis_integrals = (drag_coefficient,159lift_coefficient,160drag_coefficient_shear_force,161lift_coefficient_shear_force))162163alive_callback = AliveCallback(analysis_interval = analysis_interval)164165save_solution = SaveSolutionCallback(interval = 500,166save_initial_solution = true,167save_final_solution = true,168solution_variables = cons2prim)169170callbacks = CallbackSet(summary_callback, analysis_callback, alive_callback, save_solution)171172###############################################################################173# run the simulation174175sol = solve(ode, RDPK3SpFSAL49(thread = Trixi.True());176abstol = 1e-8, reltol = 1e-8,177ode_default_options()..., callback = callbacks)178179180