block_minres and block_gmres work on GPUs with Julia 1.11. Version 11.2.0 or later of GPUArrays.jl is also required.
If you want to use block_minres and block_gmres on previous Julia versions, you can overload the function Krylov.copy_triangle with the following code:
using KernelAbstractions, Krylov
@kernel function copy_triangle_kernel!(dest, src)
i, j = @index(Global, NTuple)
if j >= i
@inbounds dest[i, j] = src[i, j]
end
end
function Krylov.copy_triangle(Q::AbstractMatrix{FC}, R::AbstractMatrix{FC}, k::Int) where FC <: Krylov.FloatOrComplex
backend = get_backend(Q)
ndrange = (k, k)
copy_triangle_kernel!(backend)(R, Q; ndrange=ndrange)
KernelAbstractions.synchronize(backend)
endBlock-MINRES
Krylov.block_minres — Function
(X, stats) = block_minres(A, b::AbstractMatrix{FC};
M=I, ldiv::Bool=false,
atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(X, stats) = block_minres(A, B, X0::AbstractMatrix; kwargs...)Block-MINRES can be warm-started from an initial guess X0 where kwargs are the same keyword arguments as above.
Solve the Hermitian linear system AX = B of size n with p right-hand sides using block-MINRES.
Interface
To easily switch between block Krylov methods, use the generic interface krylov_solve with method = :block_minres.
For an in-place variant that reuses memory across solves, see block_minres!.
Input arguments
A: a linear operator that models a Hermitian matrix of dimensionn;B: a matrix of sizen × p.
Optional argument
X0: a matrix of sizen × pthat represents an initial guess of the solutionX.
Keyword arguments
M: linear operator that models a Hermitian positive-definite matrix of sizenused for centered preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2 * div(n,p);timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the block-Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
X: a dense matrix of sizen × p;stats: statistics collected on the run in aSimpleStatsstructure.
Krylov.block_minres! — Function
workspace = block_minres!(workspace::BlockMinresWorkspace, B; kwargs...)
workspace = block_minres!(workspace::BlockMinresWorkspace, B, X0; kwargs...)In these calls, kwargs are keyword arguments of block_minres.
See BlockMinresWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :block_minres to allocate the workspace, and krylov_solve! to run the block Krylov method in-place.
Block-GMRES
Krylov.block_gmres — Function
(X, stats) = block_gmres(A, b::AbstractMatrix{FC};
memory::Int=5, M=I, N=I, ldiv::Bool=false,
restart::Bool=false, reorthogonalization::Bool=false,
atol::T=√eps(T), rtol::T=√eps(T), itmax::Int=0,
timemax::Float64=Inf, verbose::Int=0, history::Bool=false,
callback=workspace->false, iostream::IO=kstdout)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
(X, stats) = block_gmres(A, B, X0::AbstractMatrix; kwargs...)Block-GMRES 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 n with p right-hand sides using block-GMRES.
Interface
To easily switch between block Krylov methods, use the generic interface krylov_solve with method = :block_gmres.
For an in-place variant that reuses memory across solves, see block_gmres!.
Input arguments
A: a linear operator that models a matrix of dimensionn;B: a matrix of sizen × p.
Optional argument
X0: a matrix of sizen × pthat represents an initial guess of the solutionX.
Keyword arguments
memory: ifrestart = true, the restarted version block-GMRES(k) is used withk = memory. Ifrestart = false, the parametermemoryshould be used as a hint of the number of iterations to limit dynamic memory allocations. Additional storage will be allocated if the number of iterations exceedsmemory;M: linear operator that models a nonsingular matrix of sizenused for left preconditioning;N: linear operator that models a nonsingular matrix of sizenused for right preconditioning;ldiv: define whether the preconditioners useldiv!ormul!;restart: restart the method aftermemoryiterations;reorthogonalization: reorthogonalize the new matrices of the block-Krylov basis against all previous matrix;atol: absolute stopping tolerance based on the residual norm;rtol: relative stopping tolerance based on the residual norm;itmax: the maximum number of iterations. Ifitmax=0, the default number of iterations is set to2 * div(n,p);timemax: the time limit in seconds;verbose: additional details can be displayed if verbose mode is enabled (verbose > 0). Information will be displayed everyverboseiterations;history: collect additional statistics on the run such as residual norms;callback: function or functor called ascallback(workspace)that returnstrueif the block-Krylov method should terminate, andfalseotherwise;iostream: stream to which output is logged.
Output arguments
X: a dense matrix of sizen × p;stats: statistics collected on the run in aSimpleStatsstructure.
Krylov.block_gmres! — Function
workspace = block_gmres!(workspace::BlockGmresWorkspace, B; kwargs...)
workspace = block_gmres!(workspace::BlockGmresWorkspace, B, X0; kwargs...)In these calls, kwargs are keyword arguments of block_gmres. The keyword argument memory is the only exception. It is only supported by block_gmres and is required to create a BlockGmresWorkspace. It cannot be changed later.
See BlockGmresWorkspace for instructions on how to create the workspace.
For a more generic interface, you can use krylov_workspace with method = :block_gmres to allocate the workspace, and krylov_solve! to run the block Krylov method in-place.