Skip to content

Commit 94aad1d

Browse files
committed
Add tests for TensorAlgebra.MatrixAlgebra
1 parent e68a135 commit 94aad1d

File tree

2 files changed

+149
-1
lines changed

2 files changed

+149
-1
lines changed

src/MatrixAlgebra.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -128,7 +128,7 @@ for (factorize, orth_f) in ((:factorize, :(MatrixAlgebra.orth)), (:factorize!, :
128128
else
129129
throw(ArgumentError("`orth=$orth` not supported."))
130130
end
131-
return f(A; kwargs...)
131+
return f(A; side=orth, kwargs...)
132132
end
133133
end
134134
end

test/test_matrixalgebra.jl

Lines changed: 148 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,148 @@
1+
using LinearAlgebra: I, diag
2+
using TensorAlgebra.MatrixAlgebra: MatrixAlgebra
3+
using Test: @test, @testset
4+
5+
elts = (Float32, Float64, ComplexF32, ComplexF64)
6+
7+
@testset "elt=$elt" for elt in elts
8+
A = randn(elt, 3, 2)
9+
for positive in (false, true)
10+
for (Q, R) in (MatrixAlgebra.qr(A; positive), MatrixAlgebra.qr(A; full=false, positive))
11+
@test A Q * R
12+
@test size(Q) == size(A)
13+
@test size(R) == (size(A, 2), size(A, 2))
14+
@test Q' * Q I
15+
@test Q * Q' I
16+
if positive
17+
@test all((0), real(diag(R)))
18+
@test all((0), imag(diag(R)))
19+
end
20+
end
21+
end
22+
23+
A = randn(elt, 3, 2)
24+
for positive in (false, true)
25+
Q, R = MatrixAlgebra.qr(A; full=true, positive)
26+
@test A Q * R
27+
@test size(Q) == (size(A, 1), size(A, 1))
28+
@test size(R) == size(A)
29+
@test Q' * Q I
30+
@test Q * Q' I
31+
if positive
32+
@test all((0), real(diag(R)))
33+
@test all((0), imag(diag(R)))
34+
end
35+
end
36+
37+
A = randn(elt, 2, 3)
38+
for positive in (false, true)
39+
for (L, Q) in (MatrixAlgebra.lq(A; positive), MatrixAlgebra.lq(A; full=false, positive))
40+
@test A L * Q
41+
@test size(L) == (size(A, 1), size(A, 1))
42+
@test size(Q) == size(A)
43+
@test Q * Q' I
44+
@test Q' * Q I
45+
if positive
46+
@test all((0), real(diag(L)))
47+
@test all((0), imag(diag(L)))
48+
end
49+
end
50+
end
51+
52+
A = randn(elt, 3, 2)
53+
for positive in (false, true)
54+
L, Q = MatrixAlgebra.lq(A; full=true, positive)
55+
@test A L * Q
56+
@test size(L) == size(A)
57+
@test size(Q) == (size(A, 2), size(A, 2))
58+
@test Q * Q' I
59+
@test Q' * Q I
60+
if positive
61+
@test all((0), real(diag(L)))
62+
@test all((0), imag(diag(L)))
63+
end
64+
end
65+
66+
A = randn(elt, 3, 2)
67+
for (W, C) in (MatrixAlgebra.orth(A), MatrixAlgebra.orth(A; side=:left))
68+
@test A W * C
69+
@test size(W) == size(A)
70+
@test size(C) == (size(A, 2), size(A, 2))
71+
@test W' * W I
72+
@test W * W' I
73+
end
74+
75+
A = randn(elt, 2, 3)
76+
C, W = MatrixAlgebra.orth(A; side=:right)
77+
@test A C * W
78+
@test size(C) == (size(A, 1), size(A, 1))
79+
@test size(W) == size(A)
80+
@test W * W' I
81+
@test W' * W I
82+
83+
A = randn(elt, 3, 2)
84+
for (W, P) in (MatrixAlgebra.polar(A), MatrixAlgebra.polar(A; side=:left))
85+
@test A W * P
86+
@test size(W) == size(A)
87+
@test size(P) == (size(A, 2), size(A, 2))
88+
@test W' * W I
89+
@test W * W' I
90+
@test isposdef(P)
91+
end
92+
93+
A = randn(elt, 2, 3)
94+
P, W = MatrixAlgebra.polar(A; side=:right)
95+
@test A P * W
96+
@test size(P) == (size(A, 1), size(A, 1))
97+
@test size(W) == size(A)
98+
@test W * W' I
99+
@test W' * W I
100+
@test isposdef(P)
101+
102+
A = randn(elt, 3, 2)
103+
for (W, C) in (MatrixAlgebra.factorize(A), MatrixAlgebra.factorize(A; orth=:left))
104+
@test A W * C
105+
@test size(W) == size(A)
106+
@test size(C) == (size(A, 2), size(A, 2))
107+
@test W' * W I
108+
@test W * W' I
109+
end
110+
111+
A = randn(elt, 2, 3)
112+
C, W = MatrixAlgebra.factorize(A; orth=:right)
113+
@test A C * W
114+
@test size(C) == (size(A, 1), size(A, 1))
115+
@test size(W) == size(A)
116+
@test W * W' I
117+
@test W' * W I
118+
119+
A = randn(elt, 3, 3)
120+
D, V = MatrixAlgebra.eigen(A)
121+
@test A * V V * D
122+
@test MatrixAlgebra.eigvals(A) diag(D)
123+
124+
A = randn(elt, 3, 2)
125+
for (U, S, V) in (MatrixAlgebra.svd(A), MatrixAlgebra.svd(A; full=false))
126+
@test A U * S * V
127+
@test size(U) == size(A)
128+
@test size(S) == (size(A, 2), size(A, 2))
129+
@test size(V) == (size(A, 2), size(A, 2))
130+
@test U' * U I
131+
@test U * U' I
132+
@test V * V' I
133+
@test V' * V I
134+
@test MatrixAlgebra.svdvals(A) diag(S)
135+
end
136+
137+
A = randn(elt, 3, 2)
138+
U, S, V = MatrixAlgebra.svd(A; full=true)
139+
@test A U * S * V
140+
@test size(U) == (size(A, 1), size(A, 1))
141+
@test size(S) == size(A)
142+
@test size(V) == (size(A, 2), size(A, 2))
143+
@test U' * U I
144+
@test U * U' I
145+
@test V * V' I
146+
@test V' * V I
147+
@test MatrixAlgebra.svdvals(A) diag(S)
148+
end

0 commit comments

Comments
 (0)