forked from Jutho/KrylovKit.jl
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathlinsolve.jl
139 lines (130 loc) · 6.66 KB
/
linsolve.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
# Test CG complete
@testset "CG small problem ($mode)" for mode in (:vector, :inplace, :outplace)
scalartypes = mode === :vector ? (Float32, Float64, ComplexF32, ComplexF64) :
(ComplexF64,)
@testset for T in scalartypes
A = rand(T, (n, n))
A = sqrt(A * A')
b = rand(T, n)
alg = CG(; maxiter=2n, tol=tolerance(T) * norm(b), verbosity=2) # because of loss of orthogonality, we choose maxiter = 2n
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode));
ishermitian=true, isposdef=true, maxiter=2n,
krylovdim=1, rtol=tolerance(T),
verbosity=1)
@test info.converged > 0
@test unwrapvec(b) ≈ A * unwrapvec(x)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)), x;
ishermitian=true, isposdef=true, maxiter=2n,
krylovdim=1, rtol=tolerance(T))
@test info.numops == 1
A = rand(T, (n, n))
A = sqrt(A * A')
α₀ = rand(real(T)) + 1
α₁ = rand(real(T))
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)),
wrapvec(zerovector(b), Val(mode)), alg, α₀, α₁)
@test unwrapvec(b) ≈ (α₀ * I + α₁ * A) * unwrapvec(x)
@test info.converged > 0
end
end
# Test CG complete
@testset "CG large problem ($mode)" for mode in (:vector, :inplace, :outplace)
scalartypes = mode === :vector ? (Float32, Float64, ComplexF32, ComplexF64) :
(ComplexF64,)
@testset for T in scalartypes
A = rand(T, (N, N))
A = sqrt(sqrt(A * A')) / N
b = rand(T, N)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode));
isposdef=true, maxiter=1, krylovdim=N,
rtol=tolerance(T))
@test unwrapvec(b) ≈ A * unwrapvec(x) + unwrapvec(info.residual)
α₀ = rand(real(T)) + 1
α₁ = rand(real(T))
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)),
α₀, α₁;
isposdef=true, maxiter=1, krylovdim=N,
rtol=tolerance(T))
@test unwrapvec(b) ≈ (α₀ * I + α₁ * A) * unwrapvec(x) + unwrapvec(info.residual)
end
end
# Test GMRES complete
@testset "GMRES full factorization ($mode)" for mode in (:vector, :inplace, :outplace)
scalartypes = mode === :vector ? (Float32, Float64, ComplexF32, ComplexF64) :
(ComplexF64,)
@testset for T in scalartypes
A = rand(T, (n, n)) .- one(T) / 2
b = rand(T, n)
alg = GMRES(; krylovdim=n, maxiter=2, tol=tolerance(T) * norm(b), verbosity=2)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode));
krylovdim=n, maxiter=2,
rtol=tolerance(T), verbosity=1)
@test info.converged > 0
@test unwrapvec(b) ≈ A * unwrapvec(x)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)), x;
krylovdim=n, maxiter=2,
rtol=tolerance(T))
@test info.numops == 1
A = rand(T, (n, n))
α₀ = rand(T)
α₁ = -rand(T)
x, info = @constinferred(linsolve(A, b, zerovector(b), alg, α₀, α₁))
@test unwrapvec(b) ≈ (α₀ * I + α₁ * A) * unwrapvec(x)
@test info.converged > 0
end
end
# Test GMRES with restart
@testset "GMRES with restarts ($mode)" for mode in (:vector, :inplace, :outplace)
scalartypes = mode === :vector ? (Float32, Float64, ComplexF32, ComplexF64) :
(ComplexF64,)
@testset for T in scalartypes
A = rand(T, (N, N)) .- one(T) / 2
A = I - T(9 / 10) * A / maximum(abs, eigvals(A))
b = rand(T, N)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode));
krylovdim=3 * n,
maxiter=50, rtol=tolerance(T))
@test unwrapvec(b) ≈ A * unwrapvec(x) + unwrapvec(info.residual)
A = rand(T, (N, N)) .- one(T) / 2
α₀ = maximum(abs, eigvals(A))
α₁ = -rand(T)
α₁ *= T(9) / T(10) / abs(α₁)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)), α₀,
α₁; krylovdim=3 * n,
maxiter=50, rtol=tolerance(T))
@test unwrapvec(b) ≈ (α₀ * I + α₁ * A) * unwrapvec(x) + unwrapvec(info.residual)
end
end
# Test BICGStab
@testset "BiCGStab ($mode)" for mode in (:vector, :inplace, :outplace)
scalartypes = mode === :vector ? (Float32, Float64, ComplexF32, ComplexF64) :
(ComplexF64,)
@testset for T in scalartypes
A = rand(T, (n, n)) .- one(T) / 2
A = I - T(9 / 10) * A / maximum(abs, eigvals(A))
b = rand(T, n)
alg = BiCGStab(; maxiter=4n, tol=tolerance(T) * norm(b), verbosity=2)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)),
wrapvec(zerovector(b), Val(mode)), alg)
@test info.converged > 0
@test unwrapvec(b) ≈ A * unwrapvec(x)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)), x,
alg)
@test info.numops == 1
A = rand(T, (N, N)) .- one(T) / 2
b = rand(T, N)
α₀ = maximum(abs, eigvals(A))
α₁ = -rand(T)
α₁ *= T(9) / T(10) / abs(α₁)
alg = BiCGStab(; maxiter=2, tol=tolerance(T) * norm(b), verbosity=1)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)),
wrapvec(zerovector(b), Val(mode)), alg, α₀,
α₁)
@test unwrapvec(b) ≈ (α₀ * I + α₁ * A) * unwrapvec(x) + unwrapvec(info.residual)
alg = BiCGStab(; maxiter=10 * N, tol=tolerance(T) * norm(b), verbosity=0)
x, info = @constinferred linsolve(wrapop(A, Val(mode)), wrapvec(b, Val(mode)), x,
alg, α₀, α₁)
@test info.converged > 0
@test unwrapvec(b) ≈ (α₀ * I + α₁ * A) * unwrapvec(x)
end
end