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

Add reflecting boundary conditions for Svärd-Kalisch equations #166

Merged
merged 29 commits into from
Dec 11, 2024
Merged
Show file tree
Hide file tree
Changes from 2 commits
Commits
Show all changes
29 commits
Select commit Hold shift + click to select a range
3db9929
add reflecting boundary conditions for Svärd-Kalisch equations
JoshuaLampert Dec 6, 2024
94dad02
Apply suggestions from code review
JoshuaLampert Dec 6, 2024
0d5ca48
add back accidently removed FourierDerivativeOperator
JoshuaLampert Dec 6, 2024
c7f3d4e
also implement upwind operators
JoshuaLampert Dec 8, 2024
9e9d472
reformat
JoshuaLampert Dec 8, 2024
874028b
Apply suggestions from code review
JoshuaLampert Dec 10, 2024
25ad669
Update src/equations/svaerd_kalisch_1d.jl
JoshuaLampert Dec 10, 2024
d8149c6
Merge branch 'main' into reflecting-bc-SK
JoshuaLampert Dec 10, 2024
8fe7401
set stricter time integration tolerances
JoshuaLampert Dec 10, 2024
75ec37f
lower tolerance
JoshuaLampert Dec 10, 2024
8859d06
use Cholesky decomposition
JoshuaLampert Dec 11, 2024
40468f3
Update src/equations/svaerd_kalisch_1d.jl
JoshuaLampert Dec 11, 2024
40f52cb
temporarily add lu decomposition again
JoshuaLampert Dec 11, 2024
f017619
fix scale_by_mass_matrix!
JoshuaLampert Dec 11, 2024
ca345a5
remove LU decomposition again
JoshuaLampert Dec 11, 2024
3fa969f
scale_by_mass_matrix! outside of solve_system_matrix!
JoshuaLampert Dec 11, 2024
82becee
try CI without tolerances
JoshuaLampert Dec 11, 2024
60e800e
reorder imports
JoshuaLampert Dec 11, 2024
08c04e6
fix variable name
JoshuaLampert Dec 11, 2024
94991f3
refactor to reduce code redundancy
JoshuaLampert Dec 11, 2024
320d883
fix
JoshuaLampert Dec 11, 2024
a2d8497
pass boundary_conditions
JoshuaLampert Dec 11, 2024
9bfebe3
format
JoshuaLampert Dec 11, 2024
2d50833
adjust tolerances
JoshuaLampert Dec 11, 2024
209cc76
Apply suggestions from code review
JoshuaLampert Dec 11, 2024
2308121
Update src/equations/equations.jl
JoshuaLampert Dec 11, 2024
3302cd0
begin + 1instead of 2
JoshuaLampert Dec 11, 2024
7358093
set tolerances in tests
JoshuaLampert Dec 11, 2024
7092265
add test for ArgumentError
JoshuaLampert Dec 11, 2024
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
49 changes: 49 additions & 0 deletions examples/svaerd_kalisch_1d/svaerd_kalisch_1d_basic_reflecting.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,49 @@
using OrdinaryDiffEqTsit5
using DispersiveShallowWater
using SummationByPartsOperators: MattssonNordström2004, derivative_operator

###############################################################################
# Semidiscretization of the Svärd-Kalisch equations

equations = SvaerdKalischEquations1D(gravity_constant = 1.0, eta0 = 0.0,
alpha = 0.0, beta = 1 / 3, gamma = 0.0)

initial_condition = initial_condition_manufactured_reflecting
source_terms = source_terms_manufactured_reflecting
boundary_conditions = boundary_condition_reflecting

# create homogeneous mesh
coordinates_min = 0.0
coordinates_max = 1.0
N = 512
mesh = Mesh1D(coordinates_min, coordinates_max, N)

# create solver with SBP operators of accuracy order 4
accuracy_order = 4
D1 = derivative_operator(MattssonNordström2004(),
derivative_order = 1, accuracy_order = accuracy_order,
xmin = mesh.xmin, xmax = mesh.xmax, N = mesh.N)
D2 = derivative_operator(MattssonNordström2004(),
derivative_order = 2, accuracy_order = accuracy_order,
xmin = mesh.xmin, xmax = mesh.xmax, N = mesh.N)
solver = Solver(D1, D2)

# semidiscretization holds all the necessary data structures for the spatial discretization
semi = Semidiscretization(mesh, equations, initial_condition, solver,
boundary_conditions = boundary_conditions,
source_terms = source_terms)

###############################################################################
# Create `ODEProblem` and run the simulation
tspan = (0.0, 1.0)
ode = semidiscretize(semi, tspan)
summary_callback = SummaryCallback()
analysis_callback = AnalysisCallback(semi; interval = 100,
extra_analysis_errors = (:conservation_error,),
extra_analysis_integrals = (waterheight_total,
velocity, entropy))
callbacks = CallbackSet(analysis_callback, summary_callback)

saveat = range(tspan..., length = 100)
sol = solve(ode, Tsit5(), abstol = 1e-7, reltol = 1e-7,
save_everystep = false, callback = callbacks, saveat = saveat)
3 changes: 2 additions & 1 deletion src/DispersiveShallowWater.jl
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,8 @@ using BandedMatrices: BandedMatrix
using DiffEqBase: DiffEqBase, SciMLBase, terminate!
using FastBroadcast: @..
using Interpolations: Interpolations, linear_interpolation
using LinearAlgebra: mul!, ldiv!, I, Diagonal, Symmetric, diag, lu, cholesky, cholesky!,
using LinearAlgebra: mul!, ldiv!, I, Diagonal, Symmetric, diag, lu, lu!, cholesky,
cholesky!,
issuccess
using PolynomialBases: PolynomialBases
using Printf: @printf, @sprintf
Expand Down
12 changes: 4 additions & 8 deletions src/equations/bbm_bbm_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -338,8 +338,7 @@ function create_cache(mesh, equations::BBMBBMEquations1D{BathymetryFlat},

# homogeneous Neumann boundary conditions
if solver.D1 isa DerivativeOperator ||
solver.D1 isa UniformCoupledOperator ||
solver.D1 isa FourierDerivativeOperator
solver.D1 isa UniformCoupledOperator
JoshuaLampert marked this conversation as resolved.
Show resolved Hide resolved
D1_b = BandedMatrix(solver.D1)
invImD2n = lu(I + 1 / 6 * D^2 * inv(M) * D1_b' * PdM * D1_b)
elseif solver.D1 isa UpwindOperators
Expand Down Expand Up @@ -382,8 +381,7 @@ function create_cache(mesh, equations::BBMBBMEquations1D,

# homogeneous Neumann boundary conditions
if solver.D1 isa DerivativeOperator ||
solver.D1 isa UniformCoupledOperator ||
solver.D1 isa FourierDerivativeOperator
solver.D1 isa UniformCoupledOperator
D1_b = BandedMatrix(solver.D1)
invImD2n = lu(I + 1 / 6 * inv(M) * D1_b' * PdM * K * D1_b)
elseif solver.D1 isa UpwindOperators
Expand Down Expand Up @@ -434,8 +432,7 @@ function rhs!(dq, q, t, mesh, equations::BBMBBMEquations1D, initial_condition,
end
# energy and mass conservative semidiscretization
if solver.D1 isa PeriodicDerivativeOperator ||
solver.D1 isa UniformPeriodicCoupledOperator ||
solver.D1 isa FourierDerivativeOperator
solver.D1 isa UniformPeriodicCoupledOperator
@trixi_timeit timer() "deta hyperbolic" begin
mul!(deta, solver.D1, tmp1)
end
Expand Down Expand Up @@ -501,8 +498,7 @@ function rhs!(dq, q, t, mesh, equations::BBMBBMEquations1D, initial_condition,
end
# energy and mass conservative semidiscretization
if solver.D1 isa DerivativeOperator ||
solver.D1 isa UniformCoupledOperator ||
solver.D1 isa FourierDerivativeOperator
solver.D1 isa UniformCoupledOperator
@trixi_timeit timer() "deta hyperbolic" begin
mul!(deta, solver.D1, tmp1)
end
Expand Down
7 changes: 7 additions & 0 deletions src/equations/equations.jl
Original file line number Diff line number Diff line change
Expand Up @@ -581,6 +581,13 @@ function solve_system_matrix!(dv, system_matrix, rhs,
return nothing
end

# Svärd-Kalisch equations for reflecting boundary conditions uses LU factorization
function solve_system_matrix!(dv, system_matrix, ::SvaerdKalischEquations1D, cache)
(; factorization) = cache
lu!(factorization, system_matrix)
ldiv!(factorization, dv)
end

function solve_system_matrix!(dv, system_matrix, ::Union{BBMEquation1D, BBMBBMEquations1D})
ldiv!(system_matrix, dv)
end
Expand Down
179 changes: 164 additions & 15 deletions src/equations/svaerd_kalisch_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -140,6 +140,79 @@ function source_terms_manufactured(q, x, t, equations::SvaerdKalischEquations1D)
return SVector(dq1, dq2, zero(dq1))
end

"""
initial_condition_manufactured_reflecting(x, t, equations::SvaerdKalischEquations1D, mesh)

A smooth manufactured solution for reflecting boundary conditions in combination
with [`source_terms_manufactured_reflecting`](@ref).
"""
function initial_condition_manufactured_reflecting(x, t,
equations::SvaerdKalischEquations1D,
mesh)
eta = exp(2 * t) * cospi(x)
v = exp(t) * x * sinpi(x)
b = -5 - 2 * cospi(2 * x)
D = equations.eta0 - b
return SVector(eta, v, D)
end

"""
source_terms_manufactured_reflecting(q, x, t, equations::SvaerdKalischEquations1D, mesh)

A smooth manufactured solution for reflecting boundary conditions in combination
with [`initial_condition_manufactured_reflecting`](@ref).
"""
function source_terms_manufactured_reflecting(q, x, t, equations::SvaerdKalischEquations1D)
g = gravity_constant(equations)
eta0 = still_water_surface(q, equations)
beta = equations.beta
a1 = sinpi(2 * x)
a2 = cospi(2 * x)
a9 = exp(t)
a11 = eta0 + 2.0 * a2 + 5.0
a13 = 0.2 * eta0 + 0.4 * a2 + 1
a22 = sinpi(x)
a23 = cospi(x)
a24 = exp(2 * t)
a25 = a24 * a23
a26 = a25 + 2.0 * a2 + 5.0
a27 = a9 * a22
a28 = pi * x * a9 * a23 + a27
a29 = -pi * a24 * a22 - 4.0 * pi * a1
a30 = a26 * a24
a31 = a26 * a27
dq1 = x * a29 * a27 + pi * x * a26 * a9 * a23 + a31 + 2 * a25
dq2 = 100.0 * pi * beta * a28 * a13^2 * a1 + 40.0 * pi * beta * a28 * a13 * a11 * a1 -
25.0 * beta * (-pi^2 * x * a27 + 2 * pi * a9 * a23) * a13^2 * a11 -
pi * g * a30 * a22 + x^2 * a29 * a24 * a22^2 / 2 + pi * x^2 * a30 * a22 * a23 +
x * a28 * a31 / 2 + x * a30 * a22^2 + x * a31 -
x * (x * a29 * a27 + pi * x * a26 * a9 * a23 + a31) * a27 / 2

return SVector(dq1, dq2, zero(dq1))
end

# For periodic boundary conditions
function assemble_system_matrix!(cache, h, D1,
::SvaerdKalischEquations1D)
(; M_h, minus_MD1betaD1) = cache

@.. M_h = h
scale_by_mass_matrix!(M_h, D1)

# Floating point errors accumulate a bit and the system matrix
# is not necessarily perfectly symmetric but only up to
# round-off errors. We wrap it here to avoid issues with the
# factorization.
return Symmetric(Diagonal(M_h) + minus_MD1betaD1)
end

# For reflecting boundary conditions
function assemble_system_matrix!(cache, h,
JoshuaLampert marked this conversation as resolved.
Show resolved Hide resolved
::SvaerdKalischEquations1D)
(; D1betaD1d) = cache
return Diagonal((@view h[2:(end - 1)])) - D1betaD1d
end

function create_cache(mesh, equations::SvaerdKalischEquations1D,
solver, initial_condition,
::BoundaryConditionPeriodic,
Expand Down Expand Up @@ -185,15 +258,15 @@ function create_cache(mesh, equations::SvaerdKalischEquations1D,
minus_MD1betaD1 = D1mat' * Diagonal(M_beta) * D1mat
system_matrix = let cache = (; M_h, minus_MD1betaD1)
assemble_system_matrix!(cache, h,
D1, D1mat, equations)
JoshuaLampert marked this conversation as resolved.
Show resolved Hide resolved
D1, equations)
end
elseif D1 isa PeriodicUpwindOperators
D1_central = D1.central
D1mat_minus = sparse(D1.minus)
minus_MD1betaD1 = D1mat_minus' * Diagonal(M_beta) * D1mat_minus
system_matrix = let cache = (; M_h, minus_MD1betaD1)
assemble_system_matrix!(cache, h,
D1, D1mat_minus, equations)
D1, equations)
end
else
throw(ArgumentError("unknown type of first-derivative operator: $(typeof(D1))"))
Expand All @@ -211,18 +284,48 @@ function create_cache(mesh, equations::SvaerdKalischEquations1D,
return cache
end

function assemble_system_matrix!(cache, h, D1, D1mat,
::SvaerdKalischEquations1D)
(; M_h, minus_MD1betaD1) = cache

@.. M_h = h
scale_by_mass_matrix!(M_h, D1)

# Floating point errors accumulate a bit and the system matrix
# is not necessarily perfectly symmetric but only up to
# round-off errors. We wrap it here to avoid issues with the
# factorization.
return Symmetric(Diagonal(M_h) + minus_MD1betaD1)
# Reflecting boundary conditions assume alpha = gamma = 0
function create_cache(mesh, equations::SvaerdKalischEquations1D,
solver, initial_condition,
::BoundaryConditionReflecting,
RealT, uEltype)
if !isapprox(equations.alpha, 0, atol = 1e-14) ||
!isapprox(equations.gamma, 0, atol = 1e-14)
throw(ArgumentError("Reflecting boundary conditions for Svärd-Kalisch equations only implemented for alpha = gamma = 0"))
end
JoshuaLampert marked this conversation as resolved.
Show resolved Hide resolved
D1 = solver.D1
N = nnodes(mesh)
# Assume D is independent of time and compute D evaluated at mesh points once.
D = Array{RealT}(undef, N)
x = grid(solver)
for i in eachnode(solver)
D[i] = still_waterdepth(initial_condition(x[i], 0.0, equations, mesh), equations)
end
h = ones(RealT, N)
hv = zero(h)
b = zero(h)
eta_x = zero(h)
v_x = zero(h)
h_v_x = zero(h)
hv2_x = zero(h)
beta_hat = equations.beta * D .^ 3
JoshuaLampert marked this conversation as resolved.
Show resolved Hide resolved
if D1 isa DerivativeOperator ||
D1 isa UniformCoupledOperator
D1_central = D1
D1mat = sparse(D1_central)
Pd = BandedMatrix((-1 => fill(one(real(mesh)), N - 2),), (N, N - 2))
D1betaD1d = sparse(D1mat * Diagonal(beta_hat) * D1mat * Pd)[2:(end - 1), :]
system_matrix = let cache = (; D1betaD1d)
assemble_system_matrix!(cache, h, equations)
end
factorization = lu(system_matrix)
else
throw(ArgumentError("unknown type of first-derivative operator: $(typeof(D1))"))
end
factorization = lu(system_matrix)
cache = (; factorization, D1betaD1d, D, h, hv, b, eta_x, v_x,
h_v_x, hv2_x, beta_hat, D1_central, D1)
return cache
end

# Discretization that conserves
Expand Down Expand Up @@ -315,7 +418,7 @@ function rhs!(dq, q, t, mesh, equations::SvaerdKalischEquations1D,
@trixi_timeit timer() "source terms" calc_sources!(dq, q, t, source_terms, equations,
solver)
@trixi_timeit timer() "assemble system matrix" begin
system_matrix = assemble_system_matrix!(cache, h, D1, D1_central, equations)
system_matrix = assemble_system_matrix!(cache, h, D1, equations)
end
@trixi_timeit timer() "solving elliptic system" begin
tmp1 .= dv
Expand All @@ -326,6 +429,52 @@ function rhs!(dq, q, t, mesh, equations::SvaerdKalischEquations1D,
return nothing
end

function rhs!(dq, q, t, mesh, equations::SvaerdKalischEquations1D,
initial_condition, ::BoundaryConditionReflecting, source_terms,
solver, cache)
(; D, h, hv, b, eta_x, v_x, h_v_x, hv2_x, tmp1, D1_central, D1) = cache

g = gravity_constant(equations)
eta, v = q.x
deta, dv, dD = dq.x
fill!(dD, zero(eltype(dD)))

@trixi_timeit timer() "deta hyperbolic" begin
@.. b = equations.eta0 - D
@.. h = eta - b
@.. hv = h * v

mul!(deta, D1_central, -hv)
end

# split form
@trixi_timeit timer() "dv hyperbolic" begin
mul!(eta_x, D1_central, eta)
mul!(v_x, D1_central, v)
mul!(h_v_x, D1_central, hv)
@.. tmp1 = hv * v
mul!(hv2_x, D1_central, tmp1)
@.. dv = -(0.5 * (hv2_x + hv * v_x - v * h_v_x) +
g * h * eta_x)
end

# no split form
# dv[:] = -(D1_central * (hv .* v) - v .* (D1_central * hv)+
# g * h .* eta_x)

JoshuaLampert marked this conversation as resolved.
Show resolved Hide resolved
@trixi_timeit timer() "source terms" calc_sources!(dq, q, t, source_terms, equations,
solver)
@trixi_timeit timer() "assemble system matrix" begin
system_matrix = assemble_system_matrix!(cache, h, equations)
end
@trixi_timeit timer() "solving elliptic system" begin
solve_system_matrix!((@view dv[2:(end - 1)]), system_matrix, equations, cache)
dv[1] = dv[end] = zero(eltype(dv))
end

return nothing
end

# The modified entropy/energy takes the whole `q` for every point in space
"""
DispersiveShallowWater.energy_total_modified!(e, q_global, equations::SvaerdKalischEquations1D, cache)
Expand Down
36 changes: 36 additions & 0 deletions test/test_svaerd_kalisch_1d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -17,6 +17,42 @@ end
@test_allocations(semi, sol, allocs=90_000)
end

@testitem "svaerd_kalisch_1d_basic_reflecting" setup=[
Setup,
SvaerdKalischEquations1D,
AdditionalImports
] begin
@test_trixi_include(joinpath(EXAMPLES_DIR, "svaerd_kalisch_1d_basic_reflecting.jl"),
tspan=(0.0, 1.0),
l2=[4.067039325687744e-5 6.705583873243185e-8 0.0],
linf=[0.00028373482419530305 1.4602668701318988e-7 0.0],
cons_error=[1.0484869654379814e-9 0.5469460930247622 0.0],
change_waterheight=1.0484869654379814e-9,
change_velocity=0.5469460930247622,
change_entropy=14.059432069002527)

@test_allocations(semi, sol, allocs=650_000)

# # test upwind operators
# D1 = upwind_operators(Mattsson2017; derivative_order = 1,
# accuracy_order = accuracy_order, xmin = mesh.xmin,
# xmax = mesh.xmax,
# N = mesh.N)
# D2 = sparse(D1.plus) * sparse(D1.minus)
# solver = Solver(D1, D2)
# @test_trixi_include(joinpath(EXAMPLES_DIR, "svaerd_kalisch_1d_basic_reflecting.jl"),
# tspan=(0.0, 1.0),
# solver=solver,
# l2=[0.002345799818043513 3.254313503127441e-8 0.0],
# linf=[0.05625950533062252 6.815531732318192e-7 0.0],
# cons_error=[1.6607871307809518e-9 0.5469460993745239 0.0],
# change_waterheight=-1.6607871307809518e-9,
# change_velocity=0.5469460993745239,
# change_entropy=132.10938489083918)

# @test_allocations(semi, sol, allocs=2_000)
end

@testitem "svaerd_kalisch_1d_dingemans" setup=[
Setup,
SvaerdKalischEquations1D,
Expand Down
Loading