Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

allow periodic FDSBP operators #1570

Merged
merged 4 commits into from
Jul 15, 2023
Merged
Show file tree
Hide file tree
Changes from 3 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -79,7 +79,7 @@ StaticArrayInterface = "1.4"
StaticArrays = "1"
StrideArrays = "0.1.18"
StructArrays = "0.6"
SummationByPartsOperators = "0.5.25"
SummationByPartsOperators = "0.5.41"
TimerOutputs = "0.5"
Triangulate = "2.0"
TriplotBase = "0.1"
Expand Down
3 changes: 2 additions & 1 deletion examples/tree_1d_fdsbp/elixir_advection_upwind.jl
Original file line number Diff line number Diff line change
Expand Up @@ -27,7 +27,8 @@ coordinates_min = -1.0
coordinates_max = 1.0
mesh = TreeMesh(coordinates_min, coordinates_max,
initial_refinement_level = 4,
n_cells_max = 10_000)
n_cells_max = 10_000,
periodicity = true)

semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition_sin, solver)

Expand Down
57 changes: 57 additions & 0 deletions examples/tree_1d_fdsbp/elixir_advection_upwind_periodic.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,57 @@
# !!! warning "Experimental implementation (upwind SBP)"
# This is an experimental feature and may change in future releases.

using OrdinaryDiffEq
using Trixi

###############################################################################
# semidiscretization of the linear scalar advection equation equation

equations = LinearScalarAdvectionEquation1D(1.0)

function initial_condition_sin(x, t, equation::LinearScalarAdvectionEquation1D)
return SVector(sinpi(x[1] - equations.advection_velocity[1] * t))
end

D_upw = upwind_operators(SummationByPartsOperators.periodic_derivative_operator,
accuracy_order = 4,
xmin = -1.0, xmax = 1.0,
N = 64)
flux_splitting = splitting_lax_friedrichs
solver = FDSBP(D_upw,
surface_integral = SurfaceIntegralUpwind(flux_splitting),
volume_integral = VolumeIntegralUpwind(flux_splitting))
sloede marked this conversation as resolved.
Show resolved Hide resolved

coordinates_min = -1.0
coordinates_max = 1.0
mesh = TreeMesh(coordinates_min, coordinates_max,
initial_refinement_level = 0,
n_cells_max = 10_000,
periodicity = true)

semi = SemidiscretizationHyperbolic(mesh, equations, initial_condition_sin, solver)


###############################################################################
# ODE solvers, callbacks etc.

tspan = (0.0, 2.0)
ode = semidiscretize(semi, tspan)

summary_callback = SummaryCallback()

analysis_interval = 1000
analysis_callback = AnalysisCallback(semi, interval=analysis_interval)

alive_callback = AliveCallback(analysis_interval=analysis_interval)

callbacks = CallbackSet(summary_callback,
analysis_callback, alive_callback)


###############################################################################
# run the simulation

sol = solve(ode, RDPK3SpFSAL49(); abstol=1.0e-6, reltol=1.0e-6,
ode_default_options()..., callback=callbacks);
summary_callback() # print the timer summary
2 changes: 2 additions & 0 deletions src/solvers/fdsbp_tree/fdsbp.jl
Original file line number Diff line number Diff line change
Expand Up @@ -27,6 +27,8 @@ The other arguments have the same meaning as in [`DG`](@ref) or [`DGSEM`](@ref).
"""
const FDSBP = DG{Basis} where {Basis <: AbstractDerivativeOperator}

const PeriodicFDSBP = FDSBP{Basis} where {Basis <: AbstractPeriodicDerivativeOperator}

function FDSBP(D_SBP::AbstractDerivativeOperator; surface_integral, volume_integral)
# `nothing` is passed as `mortar`
return DG(D_SBP, nothing, surface_integral, volume_integral)
Expand Down
28 changes: 26 additions & 2 deletions src/solvers/fdsbp_tree/fdsbp_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -165,6 +165,14 @@ function calc_surface_integral!(du, u, mesh::TreeMesh{1},
return nothing
end

# Periodic FDSBP operators need to use a single element without boundaries
function calc_surface_integral!(du, u, mesh::TreeMesh1D,
equations, surface_integral::SurfaceIntegralStrongForm,
sloede marked this conversation as resolved.
Show resolved Hide resolved
dg::PeriodicFDSBP, cache)
@assert nelements(dg, cache) == 1
return nothing
end

# Specialized interface flux computation because the upwind solver does
# not require a standard numerical flux (Riemann solver). The flux splitting
# already separates the solution information into right-traveling and
Expand Down Expand Up @@ -239,13 +247,25 @@ function calc_surface_integral!(du, u, mesh::TreeMesh{1},
return nothing
end

# Periodic FDSBP operators need to use a single element without boundaries
function calc_surface_integral!(du, u, mesh::TreeMesh1D,
equations, surface_integral::SurfaceIntegralUpwind,
dg::PeriodicFDSBP, cache)
@assert nelements(dg, cache) == 1
return nothing
end

# AnalysisCallback

function integrate_via_indices(func::Func, u,
mesh::TreeMesh{1}, equations,
dg::FDSBP, cache, args...; normalize = true) where {Func}
# TODO: FD. This is rather inefficient right now and allocates...
weights = diag(SummationByPartsOperators.mass_matrix(dg.basis))
M = SummationByPartsOperators.mass_matrix(dg.basis)
if M isa UniformScaling
M = M(nnodes(dg))
end
weights = diag(M)

# Initialize integral with zeros of the right shape
integral = zero(func(u, 1, 1, equations, dg, args...))
Expand All @@ -271,7 +291,11 @@ function calc_error_norms(func, u, t, analyzer,
mesh::TreeMesh{1}, equations, initial_condition,
dg::FDSBP, cache, cache_analysis)
# TODO: FD. This is rather inefficient right now and allocates...
weights = diag(SummationByPartsOperators.mass_matrix(dg.basis))
M = SummationByPartsOperators.mass_matrix(dg.basis)
if M isa UniformScaling
M = M(nnodes(dg))
end
weights = diag(M)
@unpack node_coordinates = cache.elements

# Set up data structures
Expand Down
28 changes: 26 additions & 2 deletions src/solvers/fdsbp_tree/fdsbp_2d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -201,6 +201,14 @@ function calc_surface_integral!(du, u, mesh::TreeMesh{2},
return nothing
end

# Periodic FDSBP operators need to use a single element without boundaries
function calc_surface_integral!(du, u, mesh::TreeMesh2D,
equations, surface_integral::SurfaceIntegralStrongForm,
dg::PeriodicFDSBP, cache)
@assert nelements(dg, cache) == 1
return nothing
end

# Specialized interface flux computation because the upwind solver does
# not require a standard numerical flux (Riemann solver). The flux splitting
# already separates the solution information into right-traveling and
Expand Down Expand Up @@ -295,12 +303,24 @@ function calc_surface_integral!(du, u, mesh::TreeMesh{2},
return nothing
end

# Periodic FDSBP operators need to use a single element without boundaries
function calc_surface_integral!(du, u, mesh::TreeMesh2D,
equations, surface_integral::SurfaceIntegralUpwind,
dg::PeriodicFDSBP, cache)
@assert nelements(dg, cache) == 1
return nothing
end

# AnalysisCallback
function integrate_via_indices(func::Func, u,
mesh::TreeMesh{2}, equations,
dg::FDSBP, cache, args...; normalize = true) where {Func}
# TODO: FD. This is rather inefficient right now and allocates...
weights = diag(SummationByPartsOperators.mass_matrix(dg.basis))
M = SummationByPartsOperators.mass_matrix(dg.basis)
if M isa UniformScaling
M = M(nnodes(dg))
end
weights = diag(M)

# Initialize integral with zeros of the right shape
integral = zero(func(u, 1, 1, 1, equations, dg, args...))
Expand All @@ -326,7 +346,11 @@ function calc_error_norms(func, u, t, analyzer,
mesh::TreeMesh{2}, equations, initial_condition,
dg::FDSBP, cache, cache_analysis)
# TODO: FD. This is rather inefficient right now and allocates...
weights = diag(SummationByPartsOperators.mass_matrix(dg.basis))
M = SummationByPartsOperators.mass_matrix(dg.basis)
if M isa UniformScaling
M = M(nnodes(dg))
end
weights = diag(M)
@unpack node_coordinates = cache.elements

# Set up data structures
Expand Down
28 changes: 26 additions & 2 deletions src/solvers/fdsbp_tree/fdsbp_3d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -237,6 +237,14 @@ function calc_surface_integral!(du, u, mesh::TreeMesh{3},
return nothing
end

# Periodic FDSBP operators need to use a single element without boundaries
function calc_surface_integral!(du, u, mesh::TreeMesh3D,
equations, surface_integral::SurfaceIntegralStrongForm,
dg::PeriodicFDSBP, cache)
@assert nelements(dg, cache) == 1
return nothing
end

# Specialized interface flux computation because the upwind solver does
# not require a standard numerical flux (Riemann solver). The flux splitting
# already separates the solution information into right-traveling and
Expand Down Expand Up @@ -346,13 +354,25 @@ function calc_surface_integral!(du, u, mesh::TreeMesh{3},
return nothing
end

# Periodic FDSBP operators need to use a single element without boundaries
function calc_surface_integral!(du, u, mesh::TreeMesh3D,
equations, surface_integral::SurfaceIntegralUpwind,
dg::PeriodicFDSBP, cache)
@assert nelements(dg, cache) == 1
return nothing
end

# AnalysisCallback

function integrate_via_indices(func::Func, u,
mesh::TreeMesh{3}, equations,
dg::FDSBP, cache, args...; normalize = true) where {Func}
# TODO: FD. This is rather inefficient right now and allocates...
weights = diag(SummationByPartsOperators.mass_matrix(dg.basis))
M = SummationByPartsOperators.mass_matrix(dg.basis)
if M isa UniformScaling
M = M(nnodes(dg))
end
weights = diag(M)

# Initialize integral with zeros of the right shape
integral = zero(func(u, 1, 1, 1, 1, equations, dg, args...))
Expand All @@ -378,7 +398,11 @@ function calc_error_norms(func, u, t, analyzer,
mesh::TreeMesh{3}, equations, initial_condition,
dg::FDSBP, cache, cache_analysis)
# TODO: FD. This is rather inefficient right now and allocates...
weights = diag(SummationByPartsOperators.mass_matrix(dg.basis))
M = SummationByPartsOperators.mass_matrix(dg.basis)
if M isa UniformScaling
M = M(nnodes(dg))
end
weights = diag(M)
@unpack node_coordinates = cache.elements

# Set up data structures
Expand Down
15 changes: 15 additions & 0 deletions test/test_tree_1d_fdsbp.jl
Original file line number Diff line number Diff line change
Expand Up @@ -23,6 +23,21 @@ EXAMPLES_DIR = pkgdir(Trixi, "examples", "tree_1d_fdsbp")
@test (@allocated Trixi.rhs!(du_ode, u_ode, semi, t)) < 1000
end
end

@trixi_testset "elixir_advection_upwind_periodic.jl" begin
@test_trixi_include(joinpath(EXAMPLES_DIR, "elixir_advection_upwind_periodic.jl"),
l2 = [1.1672962783692568e-5],
linf = [1.650514414558435e-5])

# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
let
t = sol.t[end]
u_ode = sol.u[end]
du_ode = similar(u_ode)
@test (@allocated Trixi.rhs!(du_ode, u_ode, semi, t)) < 1000
end
end
end

@testset "Inviscid Burgers" begin
Expand Down
18 changes: 18 additions & 0 deletions test/test_tree_2d_fdsbp.jl
Original file line number Diff line number Diff line change
Expand Up @@ -23,6 +23,24 @@ EXAMPLES_DIR = pkgdir(Trixi, "examples", "tree_2d_fdsbp")
@test (@allocated Trixi.rhs!(du_ode, u_ode, semi, t)) < 1000
end
end

@trixi_testset "elixir_advection_extended.jl with periodic operators" begin
@test_trixi_include(joinpath(EXAMPLES_DIR, "elixir_advection_extended.jl"),
l2 = [1.1239649404463432e-5],
linf = [1.5895264629195438e-5],
D_SBP = SummationByPartsOperators.periodic_derivative_operator(
derivative_order = 1, accuracy_order = 4, xmin = 0.0, xmax = 1.0, N = 40),
initial_refinement_level = 0)

# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
let
t = sol.t[end]
u_ode = sol.u[end]
du_ode = similar(u_ode)
@test (@allocated Trixi.rhs!(du_ode, u_ode, semi, t)) < 1000
end
end
end

@testset "Compressible Euler" begin
Expand Down
23 changes: 22 additions & 1 deletion test/test_tree_3d_fdsbp.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,7 @@ include("test_trixi.jl")

EXAMPLES_DIR = pkgdir(Trixi, "examples", "tree_3d_fdsbp")

@testset "Compressible Euler" begin
@testset "Linear scalar advection" begin
@trixi_testset "elixir_advection_extended.jl" begin
@test_trixi_include(joinpath(EXAMPLES_DIR, "elixir_advection_extended.jl"),
l2 = [0.005355755365412444],
Expand All @@ -23,6 +23,27 @@ EXAMPLES_DIR = pkgdir(Trixi, "examples", "tree_3d_fdsbp")
end
end

@trixi_testset "elixir_advection_extended.jl with periodic operators" begin
@test_trixi_include(joinpath(EXAMPLES_DIR, "elixir_advection_extended.jl"),
l2 = [1.3819894522373702e-8],
linf = [3.381866298113323e-8],
D_SBP = SummationByPartsOperators.periodic_derivative_operator(
derivative_order = 1, accuracy_order = 4, xmin = 0.0, xmax = 1.0, N = 10),
initial_refinement_level = 0,
tspan = (0.0, 5.0))

# Ensure that we do not have excessive memory allocations
# (e.g., from type instabilities)
let
t = sol.t[end]
u_ode = sol.u[end]
du_ode = similar(u_ode)
@test (@allocated Trixi.rhs!(du_ode, u_ode, semi, t)) < 1000
end
end
end

@testset "Compressible Euler" begin
@trixi_testset "elixir_euler_convergence.jl" begin
@test_trixi_include(joinpath(EXAMPLES_DIR, "elixir_euler_convergence.jl"),
l2 = [2.247522803543667e-5, 2.2499169224681058e-5, 2.24991692246826e-5, 2.2499169224684707e-5, 5.814121361417382e-5],
Expand Down