Skip to content

Commit ae01dfe

Browse files
authored
Make orthogonalization method chooseable for gmres (#293)
1 parent e8b795c commit ae01dfe

File tree

4 files changed

+24
-15
lines changed

4 files changed

+24
-15
lines changed

docs/src/linear_systems/gmres.md

+2-2
Original file line numberDiff line numberDiff line change
@@ -13,9 +13,9 @@ gmres!
1313

1414
The implementation pre-allocates a matrix $V$ of size `n` by `restart` whose columns form an orthonormal basis for the Krylov subspace. This allows BLAS2 operations when updating the solution vector $x$. The Hessenberg matrix is also pre-allocated.
1515

16-
Modified Gram-Schmidt is used to orthogonalize the columns of $V$.
16+
By default, modified Gram-Schmidt is used to orthogonalize the columns of $V$, since it is numerically more stable than classical Gram-Schmidt. Modified Gram-Schmidt is however inherently sequential, and if stability is not a concern, classical Gram-Schmidt can be used, which is implemented using BLAS2 operations. As a compromise the "DGKS criterion" can be used, which conditionally applies classical Gram-Schmidt repeatedly to stabilize it, and is typically one to two times slower than classical Gram-Schmidt.
1717

1818
The computation of the residual norm is implemented in a non-standard way, namely keeping track of a vector $\gamma$ in the null-space of $H_k^*$, which is the adjoint of the $(k + 1) \times k$ Hessenberg matrix $H_k$ at the $k$th iteration. Only when $x$ needs to be updated is the Hessenberg matrix mutated with Givens rotations.
1919

2020
!!! tip
21-
GMRES can be used as an [iterator](@ref Iterators). This makes it possible to access the Hessenberg matrix and Krylov basis vectors during the iterations.
21+
GMRES can be used as an [iterator](@ref Iterators). This makes it possible to access the Hessenberg matrix and Krylov basis vectors during the iterations.

src/gmres.jl

+14-6
Original file line numberDiff line numberDiff line change
@@ -28,7 +28,7 @@ Residual(order, T::Type) = Residual{T, real(T)}(
2828
one(real(T))
2929
)
3030

31-
mutable struct GMRESIterable{preclT, precrT, solT, rhsT, vecT, arnoldiT <: ArnoldiDecomp, residualT <: Residual, resT <: Real}
31+
mutable struct GMRESIterable{preclT, precrT, solT, rhsT, vecT, arnoldiT <: ArnoldiDecomp, residualT <: Residual, resT <: Real, orthmethT}
3232
Pl::preclT
3333
Pr::precrT
3434
x::solT
@@ -44,6 +44,8 @@ mutable struct GMRESIterable{preclT, precrT, solT, rhsT, vecT, arnoldiT <: Arnol
4444
maxiter::Int
4545
tol::resT
4646
β::resT
47+
48+
orth_meth::orthmethT
4749
end
4850

4951
converged(g::GMRESIterable) = g.residual.current g.tol
@@ -66,7 +68,8 @@ function iterate(g::GMRESIterable, iteration::Int=start(g))
6668
g.arnoldi.H[g.k + 1, g.k] = orthogonalize_and_normalize!(
6769
view(g.arnoldi.V, :, 1 : g.k),
6870
view(g.arnoldi.V, :, g.k + 1),
69-
view(g.arnoldi.H, 1 : g.k, g.k)
71+
view(g.arnoldi.H, 1 : g.k, g.k),
72+
g.orth_meth
7073
)
7174

7275
# Implicitly computes the residual
@@ -109,7 +112,8 @@ function gmres_iterable!(x, A, b;
109112
reltol::Real = sqrt(eps(real(eltype(b)))),
110113
restart::Int = min(20, size(A, 2)),
111114
maxiter::Int = size(A, 2),
112-
initially_zero::Bool = false)
115+
initially_zero::Bool = false,
116+
orth_meth::OrthogonalizationMethod = ModifiedGramSchmidt())
113117
T = eltype(x)
114118

115119
# Approximate solution
@@ -126,7 +130,8 @@ function gmres_iterable!(x, A, b;
126130

127131
GMRESIterable(Pl, Pr, x, b, Ax,
128132
arnoldi, residual,
129-
mv_products, restart, 1, maxiter, tolerance, residual.current
133+
mv_products, restart, 1, maxiter, tolerance, residual.current,
134+
orth_meth
130135
)
131136
end
132137

@@ -163,6 +168,7 @@ Solves the problem ``Ax = b`` with restarted GMRES.
163168
- `Pr`: right preconditioner;
164169
- `log::Bool`: keep track of the residual norm in each iteration;
165170
- `verbose::Bool`: print convergence information during the iterations.
171+
- `orth_meth::OrthogonalizationMethod = ModifiedGramSchmidt()`: orthogonalization method (ModifiedGramSchmidt(), ClassicalGramSchmidt(), DGKS())
166172
167173
# Return values
168174
@@ -184,15 +190,17 @@ function gmres!(x, A, b;
184190
maxiter::Int = size(A, 2),
185191
log::Bool = false,
186192
initially_zero::Bool = false,
187-
verbose::Bool = false)
193+
verbose::Bool = false,
194+
orth_meth::OrthogonalizationMethod = ModifiedGramSchmidt())
188195
history = ConvergenceHistory(partial = !log, restart = restart)
189196
history[:abstol] = abstol
190197
history[:reltol] = reltol
191198
log && reserve!(history, :resnorm, maxiter)
192199

193200
iterable = gmres_iterable!(x, A, b; Pl = Pl, Pr = Pr,
194201
abstol = abstol, reltol = reltol, maxiter = maxiter,
195-
restart = restart, initially_zero = initially_zero)
202+
restart = restart, initially_zero = initially_zero,
203+
orth_meth = orth_meth)
196204

197205
verbose && @printf("=== gmres ===\n%4s\t%4s\t%7s\n","rest","iter","resnorm")
198206

src/orthogonalize.jl

+6-5
Original file line numberDiff line numberDiff line change
@@ -7,9 +7,10 @@ struct ClassicalGramSchmidt <: OrthogonalizationMethod end
77
struct ModifiedGramSchmidt <: OrthogonalizationMethod end
88

99
# Default to MGS, good enough for solving linear systems.
10-
@inline orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}) where {T} = orthogonalize_and_normalize!(V, w, h, ModifiedGramSchmidt)
10+
@inline orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}) where {T} =
11+
orthogonalize_and_normalize!(V, w, h, ModifiedGramSchmidt())
1112

12-
function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}, ::Type{DGKS}) where {T}
13+
function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}, ::DGKS) where {T}
1314
# Orthogonalize using BLAS-2 ops
1415
mul!(h, adjoint(V), w)
1516
mul!(w, V, h, -one(T), one(T))
@@ -37,7 +38,7 @@ function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T},
3738
nrm
3839
end
3940

40-
function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}, ::Type{ClassicalGramSchmidt}) where {T}
41+
function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}, ::ClassicalGramSchmidt) where {T}
4142
# Orthogonalize using BLAS-2 ops
4243
mul!(h, adjoint(V), w)
4344
mul!(w, V, h, -one(T), one(T))
@@ -49,7 +50,7 @@ function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T},
4950
nrm
5051
end
5152

52-
function orthogonalize_and_normalize!(V::StridedVector{Vector{T}}, w::StridedVector{T}, h::StridedVector{T}, ::Type{ModifiedGramSchmidt}) where {T}
53+
function orthogonalize_and_normalize!(V::StridedVector{Vector{T}}, w::StridedVector{T}, h::StridedVector{T}, ::ModifiedGramSchmidt) where {T}
5354
# Orthogonalize using BLAS-1 ops
5455
for i = 1 : length(V)
5556
h[i] = dot(V[i], w)
@@ -63,7 +64,7 @@ function orthogonalize_and_normalize!(V::StridedVector{Vector{T}}, w::StridedVec
6364
nrm
6465
end
6566

66-
function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}, ::Type{ModifiedGramSchmidt}) where {T}
67+
function orthogonalize_and_normalize!(V::StridedMatrix{T}, w::StridedVector{T}, h::StridedVector{T}, ::ModifiedGramSchmidt) where {T}
6768
# Orthogonalize using BLAS-1 ops and column views.
6869
for i = 1 : size(V, 2)
6970
column = view(V, :, i)

test/orthogonalize.jl

+2-2
Original file line numberDiff line numberDiff line change
@@ -34,7 +34,7 @@ m = 3
3434
end
3535

3636
# Assuming V is a matrix
37-
@testset "Using $method" for method = (DGKS, ClassicalGramSchmidt, ModifiedGramSchmidt)
37+
@testset "Using $method" for method = (DGKS(), ClassicalGramSchmidt(), ModifiedGramSchmidt())
3838

3939
# Projection size
4040
h = zeros(T, m)
@@ -55,7 +55,7 @@ m = 3
5555

5656
# Orthogonalize w in-place
5757
w = copy(w_original)
58-
nrm = orthogonalize_and_normalize!(V_vec, w, h, ModifiedGramSchmidt)
58+
nrm = orthogonalize_and_normalize!(V_vec, w, h, ModifiedGramSchmidt())
5959

6060
is_orthonormalized(w, h, nrm)
6161
end

0 commit comments

Comments
 (0)