Preconditioners overview

This page illustrates some of the method(s) in the Julia package Preconditioners.jl.

This page was generated from a single Julia file: 1-overview.jl.

In any such Julia documentation, you can access the source code using the "Edit on GitHub" link in the top right.

The corresponding notebook can be viewed in nbviewer here: 1-overview.ipynb, and opened in binder here: 1-overview.ipynb.

Setup

Packages needed here.

using Preconditioners: DiagonalPreconditioner, CholeskyPreconditioner
using InteractiveUtils: versioninfo
using SparseArrays: sprand
using LinearAlgebra: I, cond

Overview

Preconditioning is useful for accelerating solutions to systems of equations and optimization problems.

Examples

n = 1000
A = sprand(n, n, 0.01)
A = A + A' + 10I
1000×1000 SparseArrays.SparseMatrixCSC{Float64, Int64} with 20940 stored entries:
⡻⣮⣷⣷⣿⣷⣿⣿⣷⡿⣷⣿⣷⣿⣽⣿⣿⣿⣿⢿⣷⣿⣿⣷⣿⣿⡶⣿⢿⣾⣿⣷⣯⣞⣷⣿⣶⣿⣿⣯
⢽⣿⣿⣿⡿⣿⡿⣿⣿⣿⣿⣿⣾⡿⣽⣿⣿⣿⡿⣿⣿⢾⣿⣿⣿⣾⣿⣿⣾⣿⣿⣿⣿⣿⣿⡿⣿⣿⣿⡯
⢿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣟⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣺⣿⣿⣿⣿⣿⣿⣿⣿⣿⣫
⣿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⢿⣯⣿⣿⣿⣿⣧⣯⣿⣯⣾⣾⣿⣿⣿⣟⣟⣿⣿⣿⢿⡿⣿⣯
⣽⡿⣿⣿⣿⣿⣿⣿⣿⢟⣿⣿⣿⣹⣿⣾⣿⣿⣿⣿⣿⣿⣿⣷⣿⣿⣿⣿⣿⣿⣿⣿⣿⣻⣿⣿⣿⣿⣿⣾
⣽⣿⣿⣿⣿⣿⣿⣿⣿⣿⣵⣿⣿⣿⣿⣿⣯⣟⣿⣿⣿⣿⣷⣿⣿⣿⣿⣿⣟⡿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣽⣿⣾⡿⣿⣿⣿⣿⣟⣻⣿⣿⣻⣾⣷⣿⣿⣟⣿⣿⣯⣿⣿⣿⣿⣿⣿⡾⣿⣿⡿⣿⣿⣷⣻⣿⣽⣿⣿⣷
⣷⣿⣷⣿⣿⣿⣿⣿⣻⣿⣿⣿⣽⣿⣿⣿⣿⣿⣿⣿⣿⣿⣟⣿⣿⣻⣾⣯⣯⡟⣿⣿⣿⣶⣿⣿⣿⣟⣷⣿
⣿⣿⣿⣿⣿⢿⡿⣷⣿⣿⣯⢿⣿⢿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣷⣻⣿⣿⢶⣿⣿⣿⣿⢟⣿⣿⣿⣿⣿⡿⣷
⣿⣟⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣻⣾⣿⣿⣿⣿⣾⣟⣿⣿⣿⢿⣿⣿⣯⣿⡿⣿⣿⣿⣿⣿
⣽⣿⣻⣟⣿⣿⣿⣿⣿⣿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣯⢿⣽⣿⣿⣽⣿⣿⣿⣿⣿⣿⣿⣿⢽⣿
⢿⣿⣿⣿⣿⣿⡭⣿⢿⣿⣽⣿⣿⣿⣿⣽⢿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⢿⣿⣿⣿⣟⣿⢽⣟⣻⣿⡷⣯⣞
⣿⣿⣻⣿⣿⣿⡿⣿⣿⣿⣿⣿⣿⣿⣿⣻⣿⣾⣾⢿⣯⣟⣿⣿⣿⣿⣞⣿⣿⣿⣿⣟⡿⣿⣿⣟⣿⣿⣿⡯
⣼⣯⣿⣿⣿⣿⣺⣿⣿⣿⣿⣿⣻⡿⡾⣿⢻⣟⣿⣿⣷⣿⣿⣟⣾⣽⣿⣿⣿⣿⢿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣻⣷⣾⣿⣻⣻⣿⣿⣿⣿⣿⡽⣿⣿⣯⠿⣿⣿⣿⣟⣟⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⢿⣿⣿⣿⣿⣿⣿⢿⣿⣿⣿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⢿⣿⢿⣿⣷⣿⣿⢿⣷⣿⢿⣿⣯⣟⣿⣿⣯
⣫⢿⣿⣿⣿⣿⣿⣽⣿⣻⣿⣿⢿⣿⢻⣿⣿⣵⣯⣿⣿⣿⣟⣟⣿⣯⣿⣿⣿⣿⣿⣟⣿⣿⣿⣻⢿⣿⣿⣷
⣽⣿⣿⡿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣾⣿⣿⣿⣿⣿⣯⣿⣿⣿⣹⣿⢿⣿⣿⣿⣿⡿⣿⣿⣻⣿⣿⣯⣿⣿⣿
⣼⣿⣿⣿⣿⣿⣿⡷⣿⣿⣿⣿⣷⣿⣿⢿⣿⣿⣿⣿⣿⣿⢿⡿⣿⣿⣿⣿⣿⣿⣿⣽⣿⣷⣯⣿⣿⣿⣿⣷
⡿⣿⡿⡿⡿⣻⡿⣿⣻⣿⣿⣿⢿⣿⣽⣿⢿⣯⣿⣿⣷⣷⣫⢿⡿⡿⣿⣿⣿⣿⡿⣿⢿⣿⣿⣿⢿⣿⡿⢏

Examine conditioning:

cA = cond(Matrix(A), 2)
4.391057834122149

Diagonal preconditioner

Dp = DiagonalPreconditioner(A)
DiagonalPreconditioner{Float64, Vector{Float64}}([10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0  …  10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0])

Apply preconditioner

DA = Dp \ A
1000×1000 SparseArrays.SparseMatrixCSC{Float64, Int64} with 20940 stored entries:
⡻⣮⣷⣷⣿⣷⣿⣿⣷⡿⣷⣿⣷⣿⣽⣿⣿⣿⣿⢿⣷⣿⣿⣷⣿⣿⡶⣿⢿⣾⣿⣷⣯⣞⣷⣿⣶⣿⣿⣯
⢽⣿⣿⣿⡿⣿⡿⣿⣿⣿⣿⣿⣾⡿⣽⣿⣿⣿⡿⣿⣿⢾⣿⣿⣿⣾⣿⣿⣾⣿⣿⣿⣿⣿⣿⡿⣿⣿⣿⡯
⢿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣟⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣺⣿⣿⣿⣿⣿⣿⣿⣿⣿⣫
⣿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⢿⣯⣿⣿⣿⣿⣧⣯⣿⣯⣾⣾⣿⣿⣿⣟⣟⣿⣿⣿⢿⡿⣿⣯
⣽⡿⣿⣿⣿⣿⣿⣿⣿⢟⣿⣿⣿⣹⣿⣾⣿⣿⣿⣿⣿⣿⣿⣷⣿⣿⣿⣿⣿⣿⣿⣿⣿⣻⣿⣿⣿⣿⣿⣾
⣽⣿⣿⣿⣿⣿⣿⣿⣿⣿⣵⣿⣿⣿⣿⣿⣯⣟⣿⣿⣿⣿⣷⣿⣿⣿⣿⣿⣟⡿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣽⣿⣾⡿⣿⣿⣿⣿⣟⣻⣿⣿⣻⣾⣷⣿⣿⣟⣿⣿⣯⣿⣿⣿⣿⣿⣿⡾⣿⣿⡿⣿⣿⣷⣻⣿⣽⣿⣿⣷
⣷⣿⣷⣿⣿⣿⣿⣿⣻⣿⣿⣿⣽⣿⣿⣿⣿⣿⣿⣿⣿⣿⣟⣿⣿⣻⣾⣯⣯⡟⣿⣿⣿⣶⣿⣿⣿⣟⣷⣿
⣿⣿⣿⣿⣿⢿⡿⣷⣿⣿⣯⢿⣿⢿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣷⣻⣿⣿⢶⣿⣿⣿⣿⢟⣿⣿⣿⣿⣿⡿⣷
⣿⣟⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣻⣾⣿⣿⣿⣿⣾⣟⣿⣿⣿⢿⣿⣿⣯⣿⡿⣿⣿⣿⣿⣿
⣽⣿⣻⣟⣿⣿⣿⣿⣿⣿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣯⢿⣽⣿⣿⣽⣿⣿⣿⣿⣿⣿⣿⣿⢽⣿
⢿⣿⣿⣿⣿⣿⡭⣿⢿⣿⣽⣿⣿⣿⣿⣽⢿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⢿⣿⣿⣿⣟⣿⢽⣟⣻⣿⡷⣯⣞
⣿⣿⣻⣿⣿⣿⡿⣿⣿⣿⣿⣿⣿⣿⣿⣻⣿⣾⣾⢿⣯⣟⣿⣿⣿⣿⣞⣿⣿⣿⣿⣟⡿⣿⣿⣟⣿⣿⣿⡯
⣼⣯⣿⣿⣿⣿⣺⣿⣿⣿⣿⣿⣻⡿⡾⣿⢻⣟⣿⣿⣷⣿⣿⣟⣾⣽⣿⣿⣿⣿⢿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣻⣷⣾⣿⣻⣻⣿⣿⣿⣿⣿⡽⣿⣿⣯⠿⣿⣿⣿⣟⣟⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⢿⣿⣿⣿⣿⣿⣿⢿⣿⣿⣿⣿⣿⣯⣿⣿⣿⣿⣿⣿⣿⣿⣿⢿⣿⢿⣿⣷⣿⣿⢿⣷⣿⢿⣿⣯⣟⣿⣿⣯
⣫⢿⣿⣿⣿⣿⣿⣽⣿⣻⣿⣿⢿⣿⢻⣿⣿⣵⣯⣿⣿⣿⣟⣟⣿⣯⣿⣿⣿⣿⣿⣟⣿⣿⣿⣻⢿⣿⣿⣷
⣽⣿⣿⡿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣾⣿⣿⣿⣿⣿⣯⣿⣿⣿⣹⣿⢿⣿⣿⣿⣿⡿⣿⣿⣻⣿⣿⣯⣿⣿⣿
⣼⣿⣿⣿⣿⣿⣿⡷⣿⣿⣿⣿⣷⣿⣿⢿⣿⣿⣿⣿⣿⣿⢿⡿⣿⣿⣿⣿⣿⣿⣿⣽⣿⣷⣯⣿⣿⣿⣿⣷
⡿⣿⡿⡿⡿⣻⡿⣿⣻⣿⣿⣿⢿⣿⣽⣿⢿⣯⣿⣿⣷⣷⣫⢿⡿⡿⣿⣿⣿⣿⡿⣿⢿⣿⣿⣿⢿⣿⡿⢏

Examine effect on condition number

cd = cond(Matrix(DA), 2)
4.3890070706966755

Incomplete Cholesky preconditioner with cut-off level 2

Pc = CholeskyPreconditioner(A, 2)
CholeskyPreconditioner{LimitedLDLFactorizations.LimitedLDLFactorization{Float64, Int64, Vector{Int64}, Vector{Int64}}}(LimitedLDLFactorizations.LimitedLDLFactorization{Float64, Int64, Vector{Int64}, Vector{Int64}}(true, 1000, [1, 23, 38, 59, 74, 98, 112, 132, 149, 159  …  11613, 11615, 11617, 11619, 11620, 11622, 11625, 11627, 11628, 11628], [351, 367, 368, 378, 387, 424, 428, 464, 472, 503  …  999, 984, 999, 990, 992, 998, 988, 995, 999, 1000], [351, 367, 368, 378, 387, 424, 428, 464, 472, 503  …  999, 997, 997, 999, 998, 999, 1000, 999, 1000, 1000], [0.06819749513056723, 0.0061231551073375416, 0.0903303344367049, 0.09067505874536265, 0.011265044459775808, 0.012165409189835684, 0.02113227916412107, 0.07066077220517626, 0.09862783401855937, 0.017905755410669254  …  0.017815643245973352, 0.012435792815526932, 0.03261912473080676, 0.03675850570674623, 0.08176601271926723, 0.05237384192027496, 0.08595581413395696, 0.06117564324247227, 0.08616078193686531, 0.020266686009621587], [0.06819749513056723, 0.0061231551073375416, 0.0903303344367049, 0.09067505874536265, 0.011265044459775808, 0.012165409189835684, 0.02113227916412107, 0.07066077220517626, 0.09862783401855937, 0.017905755410669254  …  -0.001592952935775128, 0.00022752748772937036, -0.005261783503069264, 0.0002993675477421, -2.9771995748751693e-5, 0.09402129934175961, 0.00045210993385300005, 2.8656161381966756e-6, 0.022550916519346223, -0.002289025522383134], 1000, [10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0  …  10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0], [10.0, 10.0, 9.977291517253173, 10.0, 9.954829994182312, 10.0, 9.938426287938743, 10.0, 10.0, 10.0  …  9.423477150523308, 8.916095260406568, 8.9138115442482, 9.545424619912444, 9.434407221329323, 9.400730954112701, 9.49399816859915, 9.273962304218477, 9.521858544695831, 9.367138991469652], [61, 721, 95, 389, 153, 615, 190, 542, 611, 313  …  986, 987, 988, 992, 993, 994, 996, 999, 1000, 425], 0.0, 10.0, 0.0, 2, [313, 167, 314, 135, 315, 316, 317, 228, 318, 319  …  257, 994, 995, 996, 299, 997, 290, 225, 998, 999], [10.275132519181179, 10.25944778275103, 10.361414804522685, 10.270537616677949, 10.331016145807128, 10.192835097427322, 10.244330511671262, 10.185349031615031, 10.169861450234123, 10.146750668819461  …  10.285412871190017, 10.518698072782989, 10.525284592085116, 10.216383767299916, 10.27022605996291, 10.284717923064381, 10.27961334981372, 10.341524265824619, 10.223104517409404, 10.29694713173352], [0.3119652957993958, 0.3122036724070051, 0.3106636709657647, 0.3120350725667515, 0.3111203938086087, 0.31322217498432603, 0.3124339425869731, 0.313337260442811, 0.31357575878253313, 0.31393266412196696  …  0.31180935092484224, 0.3083322877713262, 0.30823579848587723, 0.31286097962145215, 0.3120398054657325, 0.31181988537339084, 0.3118972963538068, 0.31096228745671306, 0.3127581240676659, 0.3116346633220077], [-2.22809558096e-313, -2.65249473925e-313, 0.047418359130034886, 0.0, 0.06701162382604312, 3.97413086e-316, 0.07827144948956215, 0.0, 0.0, 0.0  …  -0.00458798391725156, 0.0009247418077317936, 0.0008211072060041868, -0.007260107855622218, -0.0013151642934147347, -0.0027938726703073636, -0.00526308976940703, -2.968274495615187e-5, 2.882165342333294e-6, -0.0022808030970572586], [1000, 1000, 1000, 996, 999, 996, 998, 994, 992, 997  …  9969, 9969, 9969, 9969, 9969, 9969, 9970, 9971, 9971, 9971], [21, 36, 57, 72, 96, 110, 130, 147, 157, 171  …  11611, 11613, 11615, 11617, 0, 11620, 11623, 11625, 0, 0], [784, 879, 820, 756, 764, 7, 909, 553, 579, 659  …  942, 979, 992, 991, 978, 832, 826, 997, 998, 0], [1, 2, 3, 4, 5, 6, 7, 8, 9, 10  …  991, 992, 993, 994, 995, 996, 997, 998, 999, 1000], Int64[], false), 2)

Here is a quick way to handle a matrix argument. This could be done more efficiently eventually; see Issue #33.

import Base.\
\(C::CholeskyPreconditioner, A::AbstractMatrix) = reduce(hcat, [C \ c for c in eachcol(A)])
\ (generic function with 152 methods)

apply preconditioner

CA = Pc \ A
1000×1000 SparseArrays.SparseMatrixCSC{Float64, Int64} with 1000000 stored entries:
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿
⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿⣿

examine effect on condition number

cc = cond(Matrix(CA), 2)
1.423797827023208

Reproducibility

This page was generated with the following version of Julia:

io = IOBuffer(); versioninfo(io); split(String(take!(io)), '\n')
10-element Vector{SubString{String}}:
 "Julia Version 1.8.5"
 "Commit 17cfb8e65ea (2023-01-08 06:45 UTC)"
 "Platform Info:"
 "  OS: Linux (x86_64-linux-gnu)"
 "  CPU: 2 × Intel(R) Xeon(R) Platinum 8370C CPU @ 2.80GHz"
 "  WORD_SIZE: 64"
 "  LIBM: libopenlibm"
 "  LLVM: libLLVM-13.0.1 (ORCJIT, icelake-server)"
 "  Threads: 1 on 2 virtual cores"
 ""

And with the following package versions

import Pkg; Pkg.status()
Status `~/work/Preconditioners.jl/Preconditioners.jl/docs/Project.toml`
  [e30172f5] Documenter v0.27.24
  [98b081ad] Literate v2.14.0
  [af69fa37] Preconditioners v0.6.0 `~/work/Preconditioners.jl/Preconditioners.jl`
  [2f01184e] SparseArrays

This page was generated using Literate.jl.