Skip to content

Commit 720304d

Browse files
Use LinearSolve's internal methods for preconditioners (#588)
1 parent 402c27a commit 720304d

File tree

3 files changed

+15
-27
lines changed

3 files changed

+15
-27
lines changed

CHANGELOG.md

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -9,6 +9,7 @@ and this project adheres to [Semantic Versioning](https://semver.org/spec/v2.0.0
99

1010
- Change default solver detection in `eigensolve` when using `sigma` keyword argument (shift-inverse algorithm). If the operator is a `SparseMatrixCSC`, the default solver is `UMFPACKFactorization`, otherwise it is automatically chosen by LinearSolve.jl, depending on the type of the operator. ([#580])
1111
- Add keyword argument `assume_hermitian` to `liouvillian`. This allows users to disable the assumption that the Hamiltonian is Hermitian. ([#581])
12+
- Use LinearSolve's internal methods for preconditioners in `SteadyStateLinearSolver`. ([#588])
1213

1314
## [v0.38.1]
1415
Release date: 2025-10-27
@@ -362,3 +363,4 @@ Release date: 2024-11-13
362363
[#579]: https://github.com/qutip/QuantumToolbox.jl/issues/579
363364
[#580]: https://github.com/qutip/QuantumToolbox.jl/issues/580
364365
[#581]: https://github.com/qutip/QuantumToolbox.jl/issues/581
366+
[#588]: https://github.com/qutip/QuantumToolbox.jl/issues/588

docs/src/users_guide/steadystate.md

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -39,10 +39,11 @@ To change a solver, use the keyword argument `solver`, for example:
3939

4040
Although it is not obvious, the [`SteadyStateDirectSolver()`](@ref SteadyStateDirectSolver) and [`SteadyStateEigenSolver()`](@ref SteadyStateEigenSolver) methods all use an LU decomposition internally and thus can have a large memory overhead. In contrast, for [`SteadyStateLinearSolver()`](@ref SteadyStateLinearSolver), iterative algorithms provided by [`LinearSolve.jl`](https://docs.sciml.ai/LinearSolve/stable/solvers/solvers/), such as `KrylovJL_GMRES()` and `KrylovJL_BICGSTAB()`, do not factor the matrix and thus take less memory than the LU methods and allow, in principle, for extremely large system sizes. The downside is that these methods can take much longer than the direct method as the condition number of the Liouvillian matrix is large, indicating that these iterative methods require a large number of iterations for convergence.
4141

42-
To overcome this, one can provide preconditioner that solves for an approximate inverse for the (modified) Liouvillian, thus better conditioning the problem, leading to faster convergence. The left and right preconditioner can be specified as the keyword argument `Pl` and `Pr`, respectively:
42+
To overcome this, preconditioners can be defined directly in the solver algorithm using LinearSolve's precs parameter. For example:
4343
```julia
44-
solver = SteadyStateLinearSolver(alg = KrylovJL_GMRES(), Pl = Pl, Pr = Pr)
44+
SteadyStateLinearSolver(alg = KrylovJL_GMRES(; precs = (A, p) -> (ilu(A, τ = 0.01), I)))
4545
```
46+
4647
The use of a preconditioner can actually make these iterative methods faster than the other solution methods. The problem with precondioning is that it is only well defined for Hermitian matrices. Since the Liouvillian is non-Hermitian, the ability to find a good preconditioner is not guaranteed. Moreover, if a preconditioner is found, it is not guaranteed to have a good condition number.
4748

4849
Furthermore, `QuantiumToolbox` can take advantage of the Intel (Pardiso) LU solver in the Intel Math Kernel library (Intel-MKL) by using [`LinearSolve.jl`](https://docs.sciml.ai/LinearSolve/stable/) and either [`Pardiso.jl`](https://github.com/JuliaSparse/Pardiso.jl) or [`MKL_jll.jl`](https://github.com/JuliaBinaryWrappers/MKL_jll.jl):

src/steadystate.jl

Lines changed: 10 additions & 25 deletions
Original file line numberDiff line numberDiff line change
@@ -23,23 +23,20 @@ A solver which solves [`steadystate`](@ref) by finding the zero (or lowest) eige
2323
struct SteadyStateEigenSolver <: SteadyStateSolver end
2424

2525
@doc raw"""
26-
SteadyStateLinearSolver(alg = KrylovJL_GMRES(), Pl = nothing, Pr = nothing)
26+
SteadyStateLinearSolver(
27+
alg = KrylovJL_GMRES(; precs = (A, p) -> A isa SparseMatrixCSC ? (ilu(A, τ = 0.01), I) : (I, I))
28+
)
2729
2830
A solver which solves [`steadystate`](@ref) by finding the inverse of Liouvillian [`SuperOperator`](@ref) using the `alg`orithms given in [`LinearSolve.jl`](https://docs.sciml.ai/LinearSolve/stable/).
2931
3032
# Arguments
3133
- `alg::SciMLLinearSolveAlgorithm=KrylovJL_GMRES()`: algorithms given in [`LinearSolve.jl`](https://docs.sciml.ai/LinearSolve/stable/)
32-
- `Pl::Union{Function,Nothing}=nothing`: left preconditioner, see documentation [Solving for Steady-State Solutions](@ref doc:Solving-for-Steady-State-Solutions) for more details.
33-
- `Pr::Union{Function,Nothing}=nothing`: right preconditioner, see documentation [Solving for Steady-State Solutions](@ref doc:Solving-for-Steady-State-Solutions) for more details.
34+
35+
# Note
36+
Refer to [`LinearSolve.jl`](https://docs.sciml.ai/LinearSolve/stable/) for more details about the available algorithms. For example, the preconditioners can be defined directly in the solver like: `SteadyStateLinearSolver(alg = KrylovJL_GMRES(; precs = (A, p) -> (I, Diagonal(A))))`.
3437
"""
35-
Base.@kwdef struct SteadyStateLinearSolver{
36-
MT<:Union{SciMLLinearSolveAlgorithm,Nothing},
37-
PlT<:Union{Function,Nothing},
38-
PrT<:Union{Function,Nothing},
39-
} <: SteadyStateSolver
40-
alg::MT = KrylovJL_GMRES()
41-
Pl::PlT = nothing
42-
Pr::PrT = nothing
38+
Base.@kwdef struct SteadyStateLinearSolver{MT<:Union{SciMLLinearSolveAlgorithm,Nothing}} <: SteadyStateSolver
39+
alg::MT = KrylovJL_GMRES(; precs = (A, p) -> A isa SparseMatrixCSC ? (ilu(A, τ = 0.01), I) : (I, I))
4340
end
4441

4542
@doc raw"""
@@ -159,14 +156,8 @@ function _steadystate(L::QuantumObject{SuperOperator}, solver::SteadyStateLinear
159156
L_tmp = L_tmp + Tn
160157

161158
(haskey(kwargs, :Pl) || haskey(kwargs, :Pr)) && error("The use of preconditioners must be defined in the solver.")
162-
if !isnothing(solver.Pl)
163-
kwargs = merge((; kwargs...), (Pl = solver.Pl(L_tmp),))
164-
elseif isa(L_tmp, SparseMatrixCSC)
165-
kwargs = merge((; kwargs...), (Pl = ilu(L_tmp, τ = 0.01),))
166-
end
167-
!isnothing(solver.Pr) && (kwargs = merge((; kwargs...), (Pr = solver.Pr(L_tmp),)))
168159

169-
prob = LinearProblem(L_tmp, v0)
160+
prob = LinearProblem{true}(L_tmp, v0)
170161
ρss_vec = solve(prob, solver.alg; kwargs...).u
171162

172163
ρss = reshape(ρss_vec, N, N)
@@ -407,16 +398,10 @@ function _steadystate_fourier(
407398
v0[n_max*N+1] = weight
408399

409400
(haskey(kwargs, :Pl) || haskey(kwargs, :Pr)) && error("The use of preconditioners must be defined in the solver.")
410-
if !isnothing(solver.Pl)
411-
kwargs = merge((; kwargs...), (Pl = solver.Pl(M),))
412-
elseif isa(M, SparseMatrixCSC)
413-
kwargs = merge((; kwargs...), (Pl = ilu(M, τ = 0.01),))
414-
end
415-
!isnothing(solver.Pr) && (kwargs = merge((; kwargs...), (Pr = solver.Pr(M),)))
416401
!haskey(kwargs, :abstol) && (kwargs = merge((; kwargs...), (abstol = tol,)))
417402
!haskey(kwargs, :reltol) && (kwargs = merge((; kwargs...), (reltol = tol,)))
418403

419-
prob = LinearProblem(M, v0)
404+
prob = LinearProblem{true}(M, v0)
420405
ρtot = solve(prob, solver.alg; kwargs...).u
421406

422407
offset1 = n_max * N

0 commit comments

Comments
 (0)