Path: blob/main/docs/literate/src/files/non_periodic_boundaries.jl
5591 views
#src # Non-periodic boundary conditions12# # Dirichlet boundary condition3# First, let's look at the Dirichlet boundary condition [`BoundaryConditionDirichlet`](@ref).4# ````julia5# BoundaryConditionDirichlet(boundary_value_function)6# ````7# In Trixi.jl, this creates a Dirichlet boundary condition where the function `boundary_value_function`8# is used to set the values at the boundary. It can be used to create a boundary condition that sets9# exact boundary values by passing the exact solution of the equation.1011# It is important to note that standard Dirichlet boundary conditions for hyperbolic PDEs do not12# make sense in most cases. However, we are using a special weak form of the Dirichlet boundary13# condition, based on the application of the numerical surface flux. The numerical surface flux14# takes the solution value from inside the domain and the prescribed value of the outer boundary15# state as arguments, and solves an approximate Riemann problem to introduce dissipation (and16# hence stabilization) at the boundary. Hence, the performance of the Dirichlet BC depends on the17# fidelity of the numerical surface flux.18# An easy-to read introductory reference on this topic is the paper by19# [Mengaldo et al.](https://doi.org/10.2514/6.2014-2923).2021# The passed boundary value function is called with the same arguments as an initial condition22# function, i.e.23# ````julia24# boundary_value_function(x, t, equations)25# ````26# where `x` specifies the spatial coordinates, `t` is the current time, and `equations` is the27# corresponding system of equations.2829# We want to give a short example for a simulation with such a Dirichlet BC.3031# Consider the one-dimensional linear advection equation with domain $\Omega=[0, 2]$ and a constant32# zero initial condition.33using OrdinaryDiffEqLowStorageRK34using Trixi3536advection_velocity = 1.037equations = LinearScalarAdvectionEquation1D(advection_velocity)3839initial_condition_zero(x, t, equation::LinearScalarAdvectionEquation1D) = SVector(0.0)40initial_condition = initial_condition_zero4142using Plots43plot(x -> sum(initial_condition(x, 0.0, equations)), label = "initial condition",44ylim = (-1.5, 1.5))4546# Using an advection velocity of `1.0` and the (local) Lax-Friedrichs/Rusanov flux47# [`FluxLaxFriedrichs`](@ref) as a numerical surface flux, we are able to create an inflow boundary48# on the left and an outflow boundary on the right, as the Lax-Friedrichs flux is in this case an49# exact characteristics Riemann solver. We note that for more complex PDEs different strategies for50# inflow/outflow boundaries are necessary. To define the inflow values, we initialize a `boundary_value_function`.51function boundary_condition_sine_sector(x, t, equation::LinearScalarAdvectionEquation1D)52if 1 <= t <= 353scalar = sinpi(2 * sum(t - 1))54else55scalar = zero(t)56end57return SVector(scalar)58end59boundary_condition = boundary_condition_sine_sector6061# We set the BC in negative and positive x-direction.62boundary_conditions = (x_neg = BoundaryConditionDirichlet(boundary_condition),63x_pos = BoundaryConditionDirichlet(boundary_condition))64#-65solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)6667coordinates_min = (0.0,)68coordinates_max = (2.0,)6970# For the mesh type `TreeMesh` the parameter `periodicity` must be set to `false` in the71# corresponding direction.72mesh = TreeMesh(coordinates_min, coordinates_max,73initial_refinement_level = 4,74n_cells_max = 10_000,75periodicity = false)7677semi = SemidiscretizationHyperbolic(mesh, equations,78initial_condition,79solver,80boundary_conditions = boundary_conditions)8182tspan = (0.0, 6.0)83ode = semidiscretize(semi, tspan)8485analysis_callback = AnalysisCallback(semi, interval = 100)8687stepsize_callback = StepsizeCallback(cfl = 0.9)8889callbacks = CallbackSet(analysis_callback,90stepsize_callback);9192# We define some equidistant nodes for the visualization93visnodes = range(tspan[1], tspan[2], length = 300)9495# and run the simulation.96sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);97dt = 1, # solve needs some value here but it will be overwritten by the stepsize_callback98ode_default_options()..., saveat = visnodes, callback = callbacks);99100using Plots101@gif for step in eachindex(sol.u)102plot(sol.u[step], semi, ylim = (-1.5, 1.5), legend = true, label = "approximation",103title = "time t=$(round(sol.t[step], digits=5))")104scatter!([0.0], [sum(boundary_condition(SVector(0.0), sol.t[step], equations))],105label = "boundary condition")106end107108# # Other available example elixirs with non-trivial BC109# Moreover, there are other boundary conditions in Trixi.jl. For instance, you can use the slip wall110# boundary condition [`boundary_condition_slip_wall`](@ref).111112# Trixi.jl provides some interesting examples with different combinations of boundary conditions, e.g.113# using [`boundary_condition_slip_wall`](@ref) and other self-defined boundary conditions using114# [`BoundaryConditionDirichlet`](@ref).115116# For instance, there is a 2D compressible Euler setup for a Mach 3 wind tunnel flow with a forward117# facing step in the elixir [`elixir_euler_forward_step_amr.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/p4est_2d_dgsem/elixir_euler_forward_step_amr.jl)118# discretized with a [`P4estMesh`](@ref) using adaptive mesh refinement (AMR).119# ```@raw html120# <!--121# Video details122# * Source: https://www.youtube.com/watch?v=glAug1aIxio123# * Author: Andrew R. Winters (https://liu.se/en/employee/andwi94)124# * Obtain responsive code by inserting link on https://embedresponsively.com125# -->126# <style>.embed-container { position: relative; padding-bottom: 56.25%; height: 0; overflow: hidden; max-width: 100%; } .embed-container iframe, .embed-container object, .embed-container embed { position: absolute; top: 0; left: 0; width: 100%; height: 100%; }</style><div class='embed-container'><iframe src='https://www.youtube-nocookie.com/embed/glAug1aIxio' frameborder='0' allowfullscreen></iframe></div>127# ```128# Source: [`Video`](https://www.youtube.com/watch?v=glAug1aIxio) on Trixi.jl's YouTube channel [`Trixi Framework`](https://www.youtube.com/watch?v=WElqqdMhY4A)129130# A double Mach reflection problem for the 2D compressible Euler equations131# [`elixir_euler_double_mach_amr.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/p4est_2d_dgsem/elixir_euler_double_mach_amr.jl)132# exercises a special boundary conditions along the bottom of the domain that is a mixture of133# Dirichlet and slip wall.134# ```@raw html135# <!--136# Video details137# * Source: https://www.youtube.com/watch?v=WElqqdMhY4A138# * Author: Andrew R. Winters (https://liu.se/en/employee/andwi94)139# * Obtain responsive code by inserting link on https://embedresponsively.com140# -->141# <style>.embed-container { position: relative; padding-bottom: 56.25%; height: 0; overflow: hidden; max-width: 100%; } .embed-container iframe, .embed-container object, .embed-container embed { position: absolute; top: 0; left: 0; width: 100%; height: 100%; }</style><div class='embed-container'><iframe src='https://www.youtube-nocookie.com/embed/WElqqdMhY4A' frameborder='0' allowfullscreen></iframe></div>142# ```143# Source: [`Video`](https://www.youtube.com/watch?v=WElqqdMhY4A) on Trixi.jl's YouTube channel [`Trixi Framework`](https://www.youtube.com/watch?v=WElqqdMhY4A)144145# A channel flow around a cylinder at Mach 3146# [`elixir_euler_supersonic_cylinder.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/p4est_2d_dgsem/elixir_euler_supersonic_cylinder.jl)147# contains supersonic Mach 3 inflow at the left portion of the domain and supersonic outflow at the148# right portion of the domain. The top and bottom of the channel as well as the cylinder are treated149# as Euler slip wall boundaries.150# ```@raw html151# <!--152# Video details153# * Source: https://www.youtube.com/watch?v=w0A9X38cSe4154# * Author: Andrew R. Winters (https://liu.se/en/employee/andwi94)155# * Obtain responsive code by inserting link on https://embedresponsively.com156# -->157# <style>.embed-container { position: relative; padding-bottom: 56.25%; height: 0; overflow: hidden; max-width: 100%; } .embed-container iframe, .embed-container object, .embed-container embed { position: absolute; top: 0; left: 0; width: 100%; height: 100%; }</style><div class='embed-container'><iframe src='https://www.youtube-nocookie.com/embed/w0A9X38cSe4' frameborder='0' allowfullscreen></iframe></div>158# ```159# Source: [`Video`](https://www.youtube.com/watch?v=w0A9X38cSe4) on Trixi.jl's YouTube channel [`Trixi Framework`](https://www.youtube.com/watch?v=WElqqdMhY4A)160161# Furthermore, Trixi.jl also handles equations that include non-conservative terms.162# For such equations, the tuple of conservative and non-conservative surfaces fluxes is passed to the boundary condition,163# which then returns a tuple containing the boundary condition values for both the conservative and non-conservative terms.164# For instance, a 2D ideal compressible GLM-MHD setup with reflective walls can be found in the elixir ['elixir_mhd_reflective_wall.jl](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/dgmulti_2d/elixir_mhd_reflective_wall.jl).165166# ## Package versions167168# These results were obtained using the following versions.169170using InteractiveUtils171versioninfo()172173using Pkg174Pkg.status(["Trixi", "OrdinaryDiffEqLowStorageRK", "Plots"],175mode = PKGMODE_MANIFEST)176177178