Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/docs/literate/src/files/non_periodic_boundaries.jl
5591 views
1
#src # Non-periodic boundary conditions
2
3
# # Dirichlet boundary condition
4
# First, let's look at the Dirichlet boundary condition [`BoundaryConditionDirichlet`](@ref).
5
# ````julia
6
# BoundaryConditionDirichlet(boundary_value_function)
7
# ````
8
# In Trixi.jl, this creates a Dirichlet boundary condition where the function `boundary_value_function`
9
# is used to set the values at the boundary. It can be used to create a boundary condition that sets
10
# exact boundary values by passing the exact solution of the equation.
11
12
# It is important to note that standard Dirichlet boundary conditions for hyperbolic PDEs do not
13
# make sense in most cases. However, we are using a special weak form of the Dirichlet boundary
14
# condition, based on the application of the numerical surface flux. The numerical surface flux
15
# takes the solution value from inside the domain and the prescribed value of the outer boundary
16
# state as arguments, and solves an approximate Riemann problem to introduce dissipation (and
17
# hence stabilization) at the boundary. Hence, the performance of the Dirichlet BC depends on the
18
# fidelity of the numerical surface flux.
19
# An easy-to read introductory reference on this topic is the paper by
20
# [Mengaldo et al.](https://doi.org/10.2514/6.2014-2923).
21
22
# The passed boundary value function is called with the same arguments as an initial condition
23
# function, i.e.
24
# ````julia
25
# boundary_value_function(x, t, equations)
26
# ````
27
# where `x` specifies the spatial coordinates, `t` is the current time, and `equations` is the
28
# corresponding system of equations.
29
30
# We want to give a short example for a simulation with such a Dirichlet BC.
31
32
# Consider the one-dimensional linear advection equation with domain $\Omega=[0, 2]$ and a constant
33
# zero initial condition.
34
using OrdinaryDiffEqLowStorageRK
35
using Trixi
36
37
advection_velocity = 1.0
38
equations = LinearScalarAdvectionEquation1D(advection_velocity)
39
40
initial_condition_zero(x, t, equation::LinearScalarAdvectionEquation1D) = SVector(0.0)
41
initial_condition = initial_condition_zero
42
43
using Plots
44
plot(x -> sum(initial_condition(x, 0.0, equations)), label = "initial condition",
45
ylim = (-1.5, 1.5))
46
47
# Using an advection velocity of `1.0` and the (local) Lax-Friedrichs/Rusanov flux
48
# [`FluxLaxFriedrichs`](@ref) as a numerical surface flux, we are able to create an inflow boundary
49
# on the left and an outflow boundary on the right, as the Lax-Friedrichs flux is in this case an
50
# exact characteristics Riemann solver. We note that for more complex PDEs different strategies for
51
# inflow/outflow boundaries are necessary. To define the inflow values, we initialize a `boundary_value_function`.
52
function boundary_condition_sine_sector(x, t, equation::LinearScalarAdvectionEquation1D)
53
if 1 <= t <= 3
54
scalar = sinpi(2 * sum(t - 1))
55
else
56
scalar = zero(t)
57
end
58
return SVector(scalar)
59
end
60
boundary_condition = boundary_condition_sine_sector
61
62
# We set the BC in negative and positive x-direction.
63
boundary_conditions = (x_neg = BoundaryConditionDirichlet(boundary_condition),
64
x_pos = BoundaryConditionDirichlet(boundary_condition))
65
#-
66
solver = DGSEM(polydeg = 3, surface_flux = flux_lax_friedrichs)
67
68
coordinates_min = (0.0,)
69
coordinates_max = (2.0,)
70
71
# For the mesh type `TreeMesh` the parameter `periodicity` must be set to `false` in the
72
# corresponding direction.
73
mesh = TreeMesh(coordinates_min, coordinates_max,
74
initial_refinement_level = 4,
75
n_cells_max = 10_000,
76
periodicity = false)
77
78
semi = SemidiscretizationHyperbolic(mesh, equations,
79
initial_condition,
80
solver,
81
boundary_conditions = boundary_conditions)
82
83
tspan = (0.0, 6.0)
84
ode = semidiscretize(semi, tspan)
85
86
analysis_callback = AnalysisCallback(semi, interval = 100)
87
88
stepsize_callback = StepsizeCallback(cfl = 0.9)
89
90
callbacks = CallbackSet(analysis_callback,
91
stepsize_callback);
92
93
# We define some equidistant nodes for the visualization
94
visnodes = range(tspan[1], tspan[2], length = 300)
95
96
# and run the simulation.
97
sol = solve(ode, CarpenterKennedy2N54(williamson_condition = false);
98
dt = 1, # solve needs some value here but it will be overwritten by the stepsize_callback
99
ode_default_options()..., saveat = visnodes, callback = callbacks);
100
101
using Plots
102
@gif for step in eachindex(sol.u)
103
plot(sol.u[step], semi, ylim = (-1.5, 1.5), legend = true, label = "approximation",
104
title = "time t=$(round(sol.t[step], digits=5))")
105
scatter!([0.0], [sum(boundary_condition(SVector(0.0), sol.t[step], equations))],
106
label = "boundary condition")
107
end
108
109
# # Other available example elixirs with non-trivial BC
110
# Moreover, there are other boundary conditions in Trixi.jl. For instance, you can use the slip wall
111
# boundary condition [`boundary_condition_slip_wall`](@ref).
112
113
# Trixi.jl provides some interesting examples with different combinations of boundary conditions, e.g.
114
# using [`boundary_condition_slip_wall`](@ref) and other self-defined boundary conditions using
115
# [`BoundaryConditionDirichlet`](@ref).
116
117
# For instance, there is a 2D compressible Euler setup for a Mach 3 wind tunnel flow with a forward
118
# 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)
119
# discretized with a [`P4estMesh`](@ref) using adaptive mesh refinement (AMR).
120
# ```@raw html
121
# <!--
122
# Video details
123
# * Source: https://www.youtube.com/watch?v=glAug1aIxio
124
# * Author: Andrew R. Winters (https://liu.se/en/employee/andwi94)
125
# * Obtain responsive code by inserting link on https://embedresponsively.com
126
# -->
127
# <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>
128
# ```
129
# Source: [`Video`](https://www.youtube.com/watch?v=glAug1aIxio) on Trixi.jl's YouTube channel [`Trixi Framework`](https://www.youtube.com/watch?v=WElqqdMhY4A)
130
131
# A double Mach reflection problem for the 2D compressible Euler equations
132
# [`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)
133
# exercises a special boundary conditions along the bottom of the domain that is a mixture of
134
# Dirichlet and slip wall.
135
# ```@raw html
136
# <!--
137
# Video details
138
# * Source: https://www.youtube.com/watch?v=WElqqdMhY4A
139
# * Author: Andrew R. Winters (https://liu.se/en/employee/andwi94)
140
# * Obtain responsive code by inserting link on https://embedresponsively.com
141
# -->
142
# <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>
143
# ```
144
# Source: [`Video`](https://www.youtube.com/watch?v=WElqqdMhY4A) on Trixi.jl's YouTube channel [`Trixi Framework`](https://www.youtube.com/watch?v=WElqqdMhY4A)
145
146
# A channel flow around a cylinder at Mach 3
147
# [`elixir_euler_supersonic_cylinder.jl`](https://github.com/trixi-framework/Trixi.jl/blob/main/examples/p4est_2d_dgsem/elixir_euler_supersonic_cylinder.jl)
148
# contains supersonic Mach 3 inflow at the left portion of the domain and supersonic outflow at the
149
# right portion of the domain. The top and bottom of the channel as well as the cylinder are treated
150
# as Euler slip wall boundaries.
151
# ```@raw html
152
# <!--
153
# Video details
154
# * Source: https://www.youtube.com/watch?v=w0A9X38cSe4
155
# * Author: Andrew R. Winters (https://liu.se/en/employee/andwi94)
156
# * Obtain responsive code by inserting link on https://embedresponsively.com
157
# -->
158
# <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>
159
# ```
160
# Source: [`Video`](https://www.youtube.com/watch?v=w0A9X38cSe4) on Trixi.jl's YouTube channel [`Trixi Framework`](https://www.youtube.com/watch?v=WElqqdMhY4A)
161
162
# Furthermore, Trixi.jl also handles equations that include non-conservative terms.
163
# For such equations, the tuple of conservative and non-conservative surfaces fluxes is passed to the boundary condition,
164
# which then returns a tuple containing the boundary condition values for both the conservative and non-conservative terms.
165
# 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).
166
167
# ## Package versions
168
169
# These results were obtained using the following versions.
170
171
using InteractiveUtils
172
versioninfo()
173
174
using Pkg
175
Pkg.status(["Trixi", "OrdinaryDiffEqLowStorageRK", "Plots"],
176
mode = PKGMODE_MANIFEST)
177
178