-
Notifications
You must be signed in to change notification settings - Fork 56
/
Copy pathusymqr.jl
323 lines (266 loc) · 12.6 KB
/
usymqr.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
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
# An implementation of USYMQR for the solution of linear system Ax = b.
#
# This method is described in
#
# M. A. Saunders, H. D. Simon, and E. L. Yip
# Two Conjugate-Gradient-Type Methods for Unsymmetric Linear Equations.
# SIAM Journal on Numerical Analysis, 25(4), pp. 927--940, 1988.
#
# A. Buttari, D. Orban, D. Ruiz and D. Titley-Peloquin
# A tridiagonalization method for symmetric saddle-point and quasi-definite systems.
# SIAM Journal on Scientific Computing, 41(5), pp. 409--432, 2019.
#
# A. Montoison and D. Orban
# BiLQ: An Iterative Method for Nonsymmetric Linear Systems with a Quasi-Minimum Error Property.
# SIAM Journal on Matrix Analysis and Applications, 41(3), pp. 1145--1166, 2020.
#
# Alexis Montoison, <[email protected]>
# Montreal, November 2018.
export usymqr, usymqr!
"""
(x, stats) = usymqr(A, b::AbstractVector{FC}, c::AbstractVector{FC};
atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
verbose::Int=0, history::Bool=false, callback=solver->false)
`T` is an `AbstractFloat` such as `Float32`, `Float64` or `BigFloat`.
`FC` is `T` or `Complex{T}`.
(x, stats) = usymqr(A, b, c, x0::AbstractVector; kwargs...)
USYMQR can be warm-started from an initial guess `x0` where `kwargs` are the same keyword arguments as above.
Solve the linear system Ax = b of size m × n using USYMQR.
USYMQR is based on the orthogonal tridiagonalization process and requires two initial nonzero vectors `b` and `c`.
The vector `c` is only used to initialize the process and a default value can be `b` or `Aᴴb` depending on the shape of `A`.
The residual norm ‖b - Ax‖ monotonously decreases in USYMQR.
It's considered as a generalization of MINRES.
It can also be applied to under-determined and over-determined problems.
USYMQR finds the minimum-norm solution if problems are inconsistent.
The callback is called as `callback(solver)` and should return `true` if the main loop should terminate,
and `false` otherwise.
#### Input arguments
* `A`: a linear operator that models a matrix of dimension m × n;
* `b`: a vector of length m;
* `c`: a vector of length n.
#### Optional argument
* `x0`: a vector of length n that represents an initial guess of the solution x.
#### Output arguments
* `x`: a dense vector of length n;
* `stats`: statistics collected on the run in a [`SimpleStats`](@ref) structure.
#### References
* M. A. Saunders, H. D. Simon, and E. L. Yip, [*Two Conjugate-Gradient-Type Methods for Unsymmetric Linear Equations*](https://doi.org/10.1137/0725052), SIAM Journal on Numerical Analysis, 25(4), pp. 927--940, 1988.
* A. Buttari, D. Orban, D. Ruiz and D. Titley-Peloquin, [*A tridiagonalization method for symmetric saddle-point and quasi-definite systems*](https://doi.org/10.1137/18M1194900), SIAM Journal on Scientific Computing, 41(5), pp. 409--432, 2019.
* A. Montoison and D. Orban, [*BiLQ: An Iterative Method for Nonsymmetric Linear Systems with a Quasi-Minimum Error Property*](https://doi.org/10.1137/19M1290991), SIAM Journal on Matrix Analysis and Applications, 41(3), pp. 1145--1166, 2020.
"""
function usymqr end
function usymqr(A, b :: AbstractVector{FC}, c :: AbstractVector{FC}, x0 :: AbstractVector; kwargs...) where FC <: FloatOrComplex
solver = UsymqrSolver(A, b)
usymqr!(solver, A, b, c, x0; kwargs...)
return (solver.x, solver.stats)
end
function usymqr(A, b :: AbstractVector{FC}, c :: AbstractVector{FC}; kwargs...) where FC <: FloatOrComplex
solver = UsymqrSolver(A, b)
usymqr!(solver, A, b, c; kwargs...)
return (solver.x, solver.stats)
end
"""
solver = usymqr!(solver::UsymqrSolver, A, b, c; kwargs...)
solver = usymqr!(solver::UsymqrSolver, A, b, c, x0; kwargs...)
where `kwargs` are keyword arguments of [`usymqr`](@ref).
See [`UsymqrSolver`](@ref) for more details about the `solver`.
"""
function usymqr! end
function usymqr!(solver :: UsymqrSolver{T,FC,S}, A, b :: AbstractVector{FC}, c :: AbstractVector{FC},
x0 :: AbstractVector; kwargs...) where {T <: AbstractFloat, FC <: FloatOrComplex{T}, S <: DenseVector{FC}}
warm_start!(solver, x0)
usymqr!(solver, A, b, c; kwargs...)
return solver
end
function usymqr!(solver :: UsymqrSolver{T,FC,S}, A, b :: AbstractVector{FC}, c :: AbstractVector{FC};
atol :: T=√eps(T), rtol :: T=√eps(T),
itmax :: Int=0, verbose :: Int=0, history :: Bool=false,
callback = solver -> false) where {T <: AbstractFloat, FC <: FloatOrComplex{T}, S <: DenseVector{FC}}
m, n = size(A)
length(b) == m || error("Inconsistent problem size")
length(c) == n || error("Inconsistent problem size")
(verbose > 0) && @printf("USYMQR: system of %d equations in %d variables\n", m, n)
# Check type consistency
eltype(A) == FC || error("eltype(A) ≠ $FC")
ktypeof(b) <: S || error("ktypeof(b) is not a subtype of $S")
ktypeof(c) <: S || error("ktypeof(c) is not a subtype of $S")
# Compute the adjoint of A
Aᴴ = A'
# Set up workspace.
vₖ₋₁, vₖ, q, Δx, x, p = solver.vₖ₋₁, solver.vₖ, solver.q, solver.Δx, solver.x, solver.p
wₖ₋₂, wₖ₋₁, uₖ₋₁, uₖ, stats = solver.wₖ₋₂, solver.wₖ₋₁, solver.uₖ₋₁, solver.uₖ, solver.stats
warm_start = solver.warm_start
rNorms, AᴴrNorms = stats.residuals, stats.Aresiduals
reset!(stats)
r₀ = warm_start ? q : b
if warm_start
mul!(r₀, A, Δx)
@kaxpby!(n, one(FC), b, -one(FC), r₀)
end
# Initial solution x₀ and residual norm ‖r₀‖.
x .= zero(FC)
rNorm = @knrm2(m, r₀)
history && push!(rNorms, rNorm)
if rNorm == 0
stats.niter = 0
stats.solved = true
stats.inconsistent = false
stats.status = "x = 0 is a zero-residual solution"
solver.warm_start = false
return solver
end
iter = 0
itmax == 0 && (itmax = m+n)
ε = atol + rtol * rNorm
κ = zero(T)
(verbose > 0) && @printf("%5s %7s %7s\n", "k", "‖rₖ‖", "‖Aᴴrₖ₋₁‖")
kdisplay(iter, verbose) && @printf("%5d %7.1e %7s\n", iter, rNorm, "✗ ✗ ✗ ✗")
βₖ = @knrm2(m, r₀) # β₁ = ‖v₁‖ = ‖r₀‖
γₖ = @knrm2(n, c) # γ₁ = ‖u₁‖ = ‖c‖
vₖ₋₁ .= zero(FC) # v₀ = 0
uₖ₋₁ .= zero(FC) # u₀ = 0
vₖ .= r₀ ./ βₖ # v₁ = (b - Ax₀) / β₁
uₖ .= c ./ γₖ # u₁ = c / γ₁
cₖ₋₂ = cₖ₋₁ = cₖ = one(T) # Givens cosines used for the QR factorization of Tₖ₊₁.ₖ
sₖ₋₂ = sₖ₋₁ = sₖ = zero(FC) # Givens sines used for the QR factorization of Tₖ₊₁.ₖ
wₖ₋₂ .= zero(FC) # Column k-2 of Wₖ = Uₖ(Rₖ)⁻¹
wₖ₋₁ .= zero(FC) # Column k-1 of Wₖ = Uₖ(Rₖ)⁻¹
ζbarₖ = βₖ # ζbarₖ is the last component of z̅ₖ = (Qₖ)ᴴβ₁e₁
# Stopping criterion.
solved = rNorm ≤ ε
inconsistent = false
tired = iter ≥ itmax
status = "unknown"
user_requested_exit = false
while !(solved || tired || inconsistent || user_requested_exit)
# Update iteration index.
iter = iter + 1
# Continue the SSY tridiagonalization process.
# AUₖ = VₖTₖ + βₖ₊₁vₖ₊₁(eₖ)ᵀ = Vₖ₊₁Tₖ₊₁.ₖ
# AᴴVₖ = Uₖ(Tₖ)ᴴ + γₖ₊₁uₖ₊₁(eₖ)ᵀ = Uₖ₊₁(Tₖ.ₖ₊₁)ᴴ
mul!(q, A , uₖ) # Forms vₖ₊₁ : q ← Auₖ
mul!(p, Aᴴ, vₖ) # Forms uₖ₊₁ : p ← Aᴴvₖ
@kaxpy!(m, -γₖ, vₖ₋₁, q) # q ← q - γₖ * vₖ₋₁
@kaxpy!(n, -βₖ, uₖ₋₁, p) # p ← p - βₖ * uₖ₋₁
αₖ = @kdot(m, vₖ, q) # αₖ = ⟨vₖ,q⟩
@kaxpy!(m, - αₖ , vₖ, q) # q ← q - αₖ * vₖ
@kaxpy!(n, -conj(αₖ), uₖ, p) # p ← p - ᾱₖ * uₖ
βₖ₊₁ = @knrm2(m, q) # βₖ₊₁ = ‖q‖
γₖ₊₁ = @knrm2(n, p) # γₖ₊₁ = ‖p‖
# Update the QR factorization of Tₖ₊₁.ₖ = Qₖ [ Rₖ ].
# [ Oᵀ ]
# [ α₁ γ₂ 0 • • • 0 ] [ δ₁ λ₁ ϵ₁ 0 • • 0 ]
# [ β₂ α₂ γ₃ • • ] [ 0 δ₂ λ₂ • • • ]
# [ 0 • • • • • ] [ • • δ₃ • • • • ]
# [ • • • • • • • ] = Qₖ [ • • • • • 0 ]
# [ • • • • • 0 ] [ • • • • ϵₖ₋₂]
# [ • • • • γₖ ] [ • • • λₖ₋₁]
# [ • • βₖ αₖ ] [ 0 • • • • 0 δₖ ]
# [ 0 • • • • 0 βₖ₊₁] [ 0 • • • • • 0 ]
#
# If k = 1, we don't have any previous reflexion.
# If k = 2, we apply the last reflexion.
# If k ≥ 3, we only apply the two previous reflexions.
# Apply previous Givens reflections Qₖ₋₂.ₖ₋₁
if iter ≥ 3
# [cₖ₋₂ sₖ₋₂] [0 ] = [ ϵₖ₋₂ ]
# [s̄ₖ₋₂ -cₖ₋₂] [γₖ] [λbarₖ₋₁]
ϵₖ₋₂ = sₖ₋₂ * γₖ
λbarₖ₋₁ = -cₖ₋₂ * γₖ
end
# Apply previous Givens reflections Qₖ₋₁.ₖ
if iter ≥ 2
iter == 2 && (λbarₖ₋₁ = γₖ)
# [cₖ₋₁ sₖ₋₁] [λbarₖ₋₁] = [λₖ₋₁ ]
# [s̄ₖ₋₁ -cₖ₋₁] [ αₖ ] [δbarₖ]
λₖ₋₁ = cₖ₋₁ * λbarₖ₋₁ + sₖ₋₁ * αₖ
δbarₖ = conj(sₖ₋₁) * λbarₖ₋₁ - cₖ₋₁ * αₖ
end
# Compute and apply current Givens reflection Qₖ.ₖ₊₁
iter == 1 && (δbarₖ = αₖ)
# [cₖ sₖ] [δbarₖ] = [δₖ]
# [s̄ₖ -cₖ] [βₖ₊₁ ] [0 ]
(cₖ, sₖ, δₖ) = sym_givens(δbarₖ, βₖ₊₁)
# Update z̅ₖ₊₁ = Qₖ.ₖ₊₁ [ z̄ₖ ]
# [ 0 ]
#
# [cₖ sₖ] [ζbarₖ] = [ ζₖ ]
# [s̄ₖ -cₖ] [ 0 ] [ζbarₖ₊₁]
ζₖ = cₖ * ζbarₖ
ζbarₖ₊₁ = conj(sₖ) * ζbarₖ
# Compute the direction wₖ, the last column of Wₖ = Uₖ(Rₖ)⁻¹ ⟷ (Rₖ)ᵀ(Wₖ)ᵀ = (Uₖ)ᵀ.
# w₁ = u₁ / δ₁
if iter == 1
wₖ = wₖ₋₁
@kaxpy!(n, one(FC), uₖ, wₖ)
@. wₖ = wₖ / δₖ
end
# w₂ = (u₂ - λ₁w₁) / δ₂
if iter == 2
wₖ = wₖ₋₂
@kaxpy!(n, -λₖ₋₁, wₖ₋₁, wₖ)
@kaxpy!(n, one(FC), uₖ, wₖ)
@. wₖ = wₖ / δₖ
end
# wₖ = (uₖ - λₖ₋₁wₖ₋₁ - ϵₖ₋₂wₖ₋₂) / δₖ
if iter ≥ 3
@kscal!(n, -ϵₖ₋₂, wₖ₋₂)
wₖ = wₖ₋₂
@kaxpy!(n, -λₖ₋₁, wₖ₋₁, wₖ)
@kaxpy!(n, one(FC), uₖ, wₖ)
@. wₖ = wₖ / δₖ
end
# Compute solution xₖ.
# xₖ ← xₖ₋₁ + ζₖ * wₖ
@kaxpy!(n, ζₖ, wₖ, x)
# Compute ‖rₖ‖ = |ζbarₖ₊₁|.
rNorm = abs(ζbarₖ₊₁)
history && push!(rNorms, rNorm)
# Compute ‖Aᴴrₖ₋₁‖ = |ζbarₖ| * √(|δbarₖ|² + |λbarₖ|²).
AᴴrNorm = abs(ζbarₖ) * √(abs2(δbarₖ) + abs2(cₖ₋₁ * γₖ₊₁))
history && push!(AᴴrNorms, AᴴrNorm)
# Compute uₖ₊₁ and uₖ₊₁.
@. vₖ₋₁ = vₖ # vₖ₋₁ ← vₖ
@. uₖ₋₁ = uₖ # uₖ₋₁ ← uₖ
if βₖ₊₁ ≠ zero(T)
@. vₖ = q / βₖ₊₁ # βₖ₊₁vₖ₊₁ = q
end
if γₖ₊₁ ≠ zero(T)
@. uₖ = p / γₖ₊₁ # γₖ₊₁uₖ₊₁ = p
end
# Update directions for x.
if iter ≥ 2
@kswap(wₖ₋₂, wₖ₋₁)
end
# Update sₖ₋₂, cₖ₋₂, sₖ₋₁, cₖ₋₁, ζbarₖ, γₖ, βₖ.
if iter ≥ 2
sₖ₋₂ = sₖ₋₁
cₖ₋₂ = cₖ₋₁
end
sₖ₋₁ = sₖ
cₖ₋₁ = cₖ
ζbarₖ = ζbarₖ₊₁
γₖ = γₖ₊₁
βₖ = βₖ₊₁
# Update stopping criterion.
iter == 1 && (κ = atol + rtol * AᴴrNorm)
user_requested_exit = callback(solver) :: Bool
solved = rNorm ≤ ε
inconsistent = !solved && AᴴrNorm ≤ κ
tired = iter ≥ itmax
kdisplay(iter, verbose) && @printf("%5d %7.1e %7.1e\n", iter, rNorm, AᴴrNorm)
end
(verbose > 0) && @printf("\n")
tired && (status = "maximum number of iterations exceeded")
solved && (status = "solution good enough given atol and rtol")
user_requested_exit && (status = "user-requested exit")
# Update x
warm_start && @kaxpy!(n, one(FC), Δx, x)
solver.warm_start = false
# Update stats
stats.niter = iter
stats.solved = solved
stats.inconsistent = inconsistent
stats.status = status
return solver
end