BiLQR

Krylov.bilqrFunction
(x, y, stats) = bilqr(A, b::AbstractVector{FC}, c::AbstractVector{FC};
                      atol::T=√eps(T), rtol::T=√eps(T), transfer_to_bicg::Bool=true,
                      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}.

Combine BiLQ and QMR to solve adjoint systems.

[0  A] [y] = [b]
[Aᴴ 0] [x]   [c]

The relation bᴴc ≠ 0 must be satisfied. BiLQ is used for solving primal system Ax = b. QMR is used for solving dual system Aᴴy = c.

An option gives the possibility of transferring from the BiLQ point to the BiCG point, when it exists. The transfer is based on the residual norm.

BiLQR can be warm-started from initial guesses x0 and y0 with

(x, y, stats) = bilqr(A, b, c, x0, y0; kwargs...)

where kwargs are the same keyword arguments as above.

The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.

Reference

source
Krylov.bilqr!Function
solver = bilqr!(solver::BilqrSolver, A, b, c; kwargs...)
solver = bilqr!(solver::BilqrSolver, A, b, c, x0, y0; kwargs...)

where kwargs are keyword arguments of bilqr.

See BilqrSolver for more details about the solver.

source

TriLQR

Krylov.trilqrFunction
(x, y, stats) = trilqr(A, b::AbstractVector{FC}, c::AbstractVector{FC};
                       atol::T=√eps(T), rtol::T=√eps(T), transfer_to_usymcg::Bool=true,
                       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}.

Combine USYMLQ and USYMQR to solve adjoint systems.

[0  A] [y] = [b]
[Aᴴ 0] [x]   [c]

USYMLQ is used for solving primal system Ax = b. USYMQR is used for solving dual system Aᴴy = c.

An option gives the possibility of transferring from the USYMLQ point to the USYMCG point, when it exists. The transfer is based on the residual norm.

TriLQR can be warm-started from initial guesses x0 and y0 with

(x, y, stats) = trilqr(A, b, c, x0, y0; kwargs...)

where kwargs are the same keyword arguments as above.

The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.

Reference

source
Krylov.trilqr!Function
solver = trilqr!(solver::TrilqrSolver, A, b, c; kwargs...)
solver = trilqr!(solver::TrilqrSolver, A, b, c, x0, y0; kwargs...)

where kwargs are keyword arguments of trilqr.

See TrilqrSolver for more details about the solver.

source