-
Notifications
You must be signed in to change notification settings - Fork 12
/
matrix_operators_test.jl
175 lines (155 loc) · 6.56 KB
/
matrix_operators_test.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
using Test
using LinearAlgebra
using SummationByPartsOperators
using Optim: Optim # to enable loading the function space operator optimization code
# check construction of interior part of upwind operators
@testset "Check against some upwind operators" begin
N = 14
xmin_construction = 0.5
xmax_construction = 1.0
xmin_application = -0.25
xmax_application = 0.75
for acc_order in 2:6
Dm_bounded = derivative_operator(Mattsson2017(:minus ), 1, acc_order,
xmin_construction, xmax_construction, N)
Dp_bounded = derivative_operator(Mattsson2017(:plus ), 1, acc_order,
xmin_construction, xmax_construction, N)
Dc_bounded = derivative_operator(Mattsson2017(:central), 1, acc_order,
xmin_construction, xmax_construction, N)
nodes = collect(grid(Dc_bounded))
weights = diag(mass_matrix(Dc_bounded))
Dm = MatrixDerivativeOperator(xmin_application, xmax_application,
nodes, weights, Matrix(Dm_bounded), acc_order,
source_of_coefficients(Dm_bounded))
Dp = MatrixDerivativeOperator(xmin_application, xmax_application,
nodes, weights, Matrix(Dp_bounded), acc_order,
source_of_coefficients(Dp_bounded))
Dc = MatrixDerivativeOperator(xmin_application, xmax_application,
nodes, weights, Matrix(Dc_bounded), acc_order,
source_of_coefficients(Dc_bounded))
D = UpwindOperators(Dm, Dc, Dp)
for compact in (true, false)
show(IOContext(devnull, :compact=>compact), D)
show(IOContext(devnull, :compact=>compact), Dm)
show(IOContext(devnull, :compact=>compact), Dp)
show(IOContext(devnull, :compact=>compact), Dc)
summary(IOContext(devnull, :compact=>compact), D)
summary(IOContext(devnull, :compact=>compact), Dm)
summary(IOContext(devnull, :compact=>compact), Dp)
summary(IOContext(devnull, :compact=>compact), Dc)
end
M = mass_matrix(D)
@test M == mass_matrix(Dm)
@test M == mass_matrix(Dp)
@test M == mass_matrix(Dc)
x = grid(D)
@test x == grid(Dm)
@test x == grid(Dp)
@test x == grid(Dc)
@test issymmetric(Dm) == false
@test issymmetric(Dp) == false
@test issymmetric(Dc) == false
D_reference = upwind_operators(Mattsson2017;
derivative_order = 1, accuracy_order = acc_order,
xmin = xmin_application, xmax = xmax_application,
N = N)
@test x ≈ grid(D_reference)
@test M ≈ mass_matrix(D_reference)
u = sinpi.(x)
@test D.minus * u ≈ D_reference.minus * u
@test D.plus * u ≈ D_reference.plus * u
@test D.central * u ≈ D_reference.central * u
du = copy(u)
du_reference = copy(u)
α = 1.23
mul!(du, D.minus, u, α)
mul!(du_reference, D_reference.minus, u, α)
@test du ≈ du_reference
mul!(du, D.plus, u, α)
mul!(du_reference, D_reference.plus, u, α)
@test du ≈ du_reference
mul!(du, D.central, u, α)
mul!(du_reference, D_reference.central, u, α)
@test du ≈ du_reference
du = copy(u)
du_reference = copy(u)
β = 4.56
mul!(du, D.minus, u, α, β)
mul!(du_reference, D_reference.minus, u, α, β)
@test du ≈ du_reference
mul!(du, D.plus, u, α, β)
mul!(du_reference, D_reference.plus, u, α, β)
@test du ≈ du_reference
mul!(du, D.central, u, α, β)
mul!(du_reference, D_reference.central, u, α, β)
@test du ≈ du_reference
@test integrate(abs2, u, Dm) ≈ integrate(abs2, u, D_reference.minus)
@test integrate(abs2, u, Dp) ≈ integrate(abs2, u, D_reference.plus)
@test integrate(abs2, u, Dc) ≈ integrate(abs2, u, D_reference.central)
u = sinpi.(x)
u_reference = copy(u)
SummationByPartsOperators.scale_by_mass_matrix!(u, Dm)
SummationByPartsOperators.scale_by_mass_matrix!(u_reference, D_reference.minus)
@test u ≈ u_reference
SummationByPartsOperators.scale_by_inverse_mass_matrix!(u, Dm)
SummationByPartsOperators.scale_by_inverse_mass_matrix!(u_reference, D_reference.minus)
@test u ≈ u_reference
@test SummationByPartsOperators.get_weight(Dm, 1) == left_boundary_weight(D)
@test SummationByPartsOperators.get_weight(Dm, N) == right_boundary_weight(D)
@test SummationByPartsOperators.lower_bandwidth(Dm) == size(Dm, 1) - 1
@test SummationByPartsOperators.upper_bandwidth(Dm) == size(Dm, 1) - 1
end
end
if VERSION >= v"1.9"
@testset "Function space operators" begin
N = 5
x_min = -1.0
x_max = 1.0
nodes = collect(range(x_min, x_max, length=N))
source = GlaubitzNordströmÖffner2023()
for compact in (true, false)
show(IOContext(devnull, :compact=>compact), source)
end
B = zeros(N, N)
B[1, 1] = -1.0
B[N, N] = 1.0
let basis_functions = [x -> x^i for i in 0:3]
D = function_space_operator(basis_functions, nodes, source)
# Only first-derivative operators are implemented yet
@test_throws ArgumentError function_space_operator(basis_functions, nodes, source; derivative_order = 2)
@test grid(D) ≈ nodes
@test all(isapprox.(D * ones(N), zeros(N); atol = 1e-13))
@test D * nodes ≈ ones(N)
@test D * (nodes .^ 2) ≈ 2 * nodes
@test D * (nodes .^ 3) ≈ 3 * (nodes .^ 2)
M = mass_matrix(D)
@test M * D.D + D.D' * M ≈ B
end
let basis_functions = [one, identity, exp]
D = function_space_operator(basis_functions, nodes, source)
@test grid(D) ≈ nodes
@test all(isapprox.(D * ones(N), zeros(N); atol = 1e-13))
@test D * nodes ≈ ones(N)
@test D * exp.(nodes) ≈ exp.(nodes)
M = mass_matrix(D)
@test M * D.D + D.D' * M ≈ B
end
# test non-equidistant nodes generated by `nodes = [0.0, rand(8)..., 1.0]`
nodes = [0.0, 0.01585580467018155, 0.18010381213204507, 0.270467434432868,
0.37699483985320303, 0.5600831197666554, 0.5698824835924449, 0.623949064816263,
0.8574665549914025, 1.0]
N = length(nodes)
B = zeros(N, N)
B[1, 1] = -1.0
B[N, N] = 1.0
let basis_functions = [one, identity, exp]
D = function_space_operator(basis_functions, nodes, source)
@test grid(D) ≈ nodes
@test all(isapprox.(D * ones(N), zeros(N); atol = 1e-11))
@test D * nodes ≈ ones(N)
@test D * exp.(nodes) ≈ exp.(nodes)
M = mass_matrix(D)
@test M * D.D + D.D' * M ≈ B
end
end
end