Factorizations & Preconditioners
This functionality is deprecated as of version 1.6 an will be removed with version 2.0. Use the integration with LinearSolve.jl.
Factorizations
In this package, preconditioners and LU factorizations are both seen as complete or approximate factorizations. Correspondingly we provide a common API for them.
ExtendableSparse.AbstractLUFactorization
ExtendableSparse.factorize! — Method
factorize!(factorization, matrix)Update or create factorization, possibly reusing information from the current state. This method is aware of pattern changes.
ExtendableSparse.issolver — Method
issolver(factorization)Determine if factorization is a solver or not
ExtendableSparse.update! — Method
update!(factorization)Update factorization after matrix update.
ExtendableSparse.AbstractFactorization — Type
abstract type AbstractFactorizationAbstract type for a factorization with ExtandableSparseMatrix.
This type is meant to be a "type flexible" (with respect to the matrix element type) and lazily construcdet (can be constructed without knowing the matrix, and updated later) LU factorization or preconditioner. It wraps different concrete, type fixed factorizations which shall provide the usual ldiv! methods.
Any such preconditioner/factorization MyFact should have the following fields
A::ExtendableSparseMatrix
factorization
phash::UInt64and provide methods
MyFact(;kwargs...)
update!(precon::MyFact)The idea is that, depending if the matrix pattern has changed, different steps are needed to update the preconditioner.
LU Factorizations
Handling of the LU factorizations is meant to support a workflow where sequences of problems are solved based on the same matrix, where one possibly wants to reuse existing symbolic factorization data.
The support comes in two flavors.
- Using
factorize!which can work as a drop-in replacement forlu!:
using ExtendableSparse, LinearAlgebra
A = fdrand(20, 20, 1; matrixtype = ExtendableSparseMatrix)
n = size(A, 1)
b = rand(n)
factorization = SparspakLU()
factorize!(factorization, A)
nm1 = norm(factorization \ b)
# mock update from Newton etc.
for i = 4:(n - 3)
A[i, i + 3] -= 1.0e-4
end
factorize!(factorization, A)
nm2 = norm(factorization \ b)
nm1, nm2(2102.4842258563503, 2145.4580760918147)- Using
update!, where the matrix only needs to be given at construction time:
using ExtendableSparse, LinearAlgebra
A = fdrand(20, 20, 1; matrixtype = ExtendableSparseMatrix)
n = size(A, 1)
b = rand(n)
factorization = CholeskyFactorization(A)
nm1 = norm(factorization \ b)
# mock update from Newton etc.
for i = 4:(n - 3)
A[i, i + 3] -= 1.0e-4
A[i - 3, i] -= 1.0e-4
end
update!(factorization)
nm2 = norm(factorization \ b)
nm1, nm2(1844.0479123559924, 1988.812734244643)ExtendableSparse.LUFactorization — Type
LUFactorization()
LUFactorization(matrix)Default LU Factorization. Maps to Sparspak.jl for non-GPL builds, otherwise to UMFPACK.
ExtendableSparse.SparspakLU — Type
SparspakLU()
SparspakLU(matrix) LU factorization based on Sparspak.jl (P.Krysl's Julia re-implementation of Sparspak by George & Liu)
ExtendableSparse.AbstractLUFactorization — Type
abstract type AbstractLUFactorization <: AbstractFactorizationAbstract subtype for (full) LU factorizations
ExtendableSparse.CholeskyFactorization — Type
CholeskyFactorization(;valuetype=Float64, indextype=Int64) CholeskyFactorization(matrix)
Default Cholesky factorization via cholmod.
Base.:\ — Function
A\ for ExtendableSparse. It calls the LU factorization form Sparspak.jl, unless GPL components are allowed in the Julia sysimage and the floating point type of the matrix is Float64 or Complex64. In that case, Julias standard `` is called, which is realized via UMFPACK.
\(symm_ext, b)
\ for Symmetric{ExtendableSparse}
\(symm_ext, b)
\ for Hermitian{ExtendableSparse}
lufact
hsSolve LU factorization problem.
Preconditioners
The API is similar to that for LU factorizations.
The support comes in two flavors.
- Using
factorize!:
using ExtendableSparse, LinearAlgebra
using IterativeSolvers, IncompleteLU
A = fdrand(20, 20, 1; matrixtype = ExtendableSparseMatrix)
n = size(A, 1)
b = rand(n)
preconditioner = ILUTPreconditioner(; droptol = 1.0e-2)
factorize!(preconditioner, A)
# mock update from Newton etc.
nm1 = norm(bicgstabl(A, b, 1; Pl = preconditioner))
for i = 4:(n - 3)
A[i, i + 3] -= 1.0e-4
end
factorize!(preconditioner, A)
nm2 = norm(bicgstabl(A, b, 1; Pl = preconditioner))
nm1, nm2(2031.839083792524, 2073.1495516885534)- Using
update!:
using ExtendableSparse, LinearAlgebra
using IterativeSolvers
A = fdrand(20, 20, 1; matrixtype = ExtendableSparseMatrix)
n = size(A, 1)
b = rand(n)
preconditioner = ILU0Preconditioner(A)
nm1 = norm(cg(A, b; Pl = preconditioner))
# mock update from Newton etc.
for i = 4:(n - 3)
A[i, i + 3] -= 1.0e-4
A[i - 3, i] -= 1.0e-4
end
update!(preconditioner)
nm2 = norm(cg(A, b; Pl = preconditioner))
nm1, nm2(1987.0905724726947, 2067.565125992075)ExtendableSparse.AbstractPreconditioner — Type
abstract type AbstractPreconditioner <: AbstractFactorizationAbstract subtype for preconditioners
ExtendableSparse.ILUZeroPreconditioner — Type
ILUZeroPreconditioner()
ILUZeroPreconditioner(matrix)Incomplete LU preconditioner with zero fill-in using ILUZero.jl. This preconditioner also calculates and stores updates to the off-diagonal entries and thus delivers better convergence than the ILU0Preconditioner.
ExtendableSparse.PointBlockILUZeroPreconditioner — Type
PointBlockILUZeroPreconditioner(;blocksize)
PointBlockILUZeroPreconditioner(matrix;blocksize)Incomplete LU preconditioner with zero fill-in using ILUZero.jl. This preconditioner also calculates and stores updates to the off-diagonal entries and thus delivers better convergence than the ILU0Preconditioner.
ExtendableSparse.BlockPreconBuilder — Type
BlockPreconBuilder(;precs=UMFPACKPreconBuilder(),
partitioning = A -> [1:size(A,1)]Return callable object constructing a left block Jacobi preconditioner from partition of unknowns.
partitioning(A)shall return a vector of AbstractVectors describing the indices of the partitions of the matrix. For a matrix of sizen x n, e.g. partitioning could be[ 1:n÷2, (n÷2+1):n]or [ 1:2:n, 2:2:n].precs(A,p)shall return a left precondioner for a matrix block.
ExtendableSparse.BlockPreconditioner — Type
BlockPreconditioner(;partitioning, factorization=LUFactorization)Create a block preconditioner from partition of unknowns given by partitioning, a vector of AbstractVectors describing the indices of the partitions of the matrix. For a matrix of size n x n, e.g. partitioning could be [ 1:n÷2, (n÷2+1):n] or [ 1:2:n, 2:2:n]. Factorization is a callable (Function or struct) which allows to create a factorization (with ldiv! methods) from a submatrix of A.
ExtendableSparse.allow_views — Method
allow_views(::preconditioner_type)Factorizations on matrix partitions within a block preconditioner may or may not work with array views. E.g. the umfpack factorization cannot work with views, while ILUZeroPreconditioner can. Implementing a method for allow_views returning false resp. true allows to dispatch to the proper case.
ExtendableSparse.AMGCL_AMGPreconditioner — Function
AMGCL_AMGPreconditioner(;kwargs...)
AMGCL_AMGPreconditioner(matrix;kwargs...)Create the AMGCL_AMGPreconditioner wrapping AMG preconditioner from AMGCLWrap.jl For kwargs see there.
ExtendableSparse.AMGCL_RLXPreconditioner — Function
AMGCL_RLXPreconditioner(;kwargs...)
AMGCL_RLXPreconditioner(matrix;kwargs...)Create the AMGCL_RLXPreconditioner wrapping RLX preconditioner from AMGCLWrap.jl
ExtendableSparse.SA_AMGPreconditioner — Function
SA_AMGPreconditioner(;kwargs...)
SA_AMGPreconditioner(matrix;kwargs...)Create the SA_AMGPreconditioner wrapping the smoothed aggregation AMG preconditioner from AlgebraicMultigrid.jl For kwargs see there.
ExtendableSparse.RS_AMGPreconditioner — Function
RS_AMGPreconditioner(;kwargs...)
RS_AMGPreconditioner(matrix;kwargs...)Create the RS_AMGPreconditioner wrapping the Ruge-Stüben AMG preconditioner from AlgebraicMultigrid.jl For kwargs see there.
ExtendableSparse.ILUTPreconditioner — Function
ILUTPreconditioner(;droptol=1.0e-3)
ILUTPreconditioner(matrix; droptol=1.0e-3)Create the ILUTPreconditioner wrapping the one from IncompleteLU.jl For using this, you need to issue using IncompleteLU.
Experimental/deprecated
ExtendableSparse.PardisoLU — Function
PardisoLU(;iparm::Vector,
dparm::Vector,
mtype::Int)
PardisoLU(matrix; iparm,dparm,mtype)LU factorization based on pardiso. For using this, you need to issue using Pardiso and have the pardiso library from pardiso-project.org installed.
The optional keyword arguments mtype, iparm and dparm are Pardiso internal parameters.
Forsetting them, one can also access the PardisoSolver e.g. like
using Pardiso
plu=PardisoLU()
Pardiso.set_iparm!(plu.ps,5,13.0)ExtendableSparse.MKLPardisoLU — Function
MKLPardisoLU(;iparm::Vector, mtype::Int)
MKLPardisoLU(matrix; iparm, mtype)LU factorization based on pardiso. For using this, you need to issue using Pardiso. This version uses the early 2000's fork in Intel's MKL library.
The optional keyword arguments mtype and iparm are Pardiso internal parameters.
For setting them you can also access the PardisoSolver e.g. like
using Pardiso
plu=MKLPardisoLU()
Pardiso.set_iparm!(plu.ps,5,13.0)ExtendableSparse.AMGPreconditioner — Function
AMGPreconditioner(;max_levels=10, max_coarse=10)
AMGPreconditioner(matrix;max_levels=10, max_coarse=10)Create the AMGPreconditioner wrapping the Ruge-Stüben AMG preconditioner from AlgebraicMultigrid.jl
Deprecated in favor of RS_AMGPreconditioner
ExtendableSparse.JacobiPreconditioner — Type
JacobiPreconditioner()
JacobiPreconditioner(matrix)Jacobi preconditioner.
ExtendableSparse.ParallelJacobiPreconditioner — Type
ParallelJacobiPreconditioner()
ParallelJacobiPreconditioner(matrix)ParallelJacobi preconditioner.
ExtendableSparse.ILU0Preconditioner — Type
ILU0Preconditioner()
ILU0Preconditioner(matrix)Incomplete LU preconditioner with zero fill-in, without modification of off-diagonal entries, so it delivers slower convergende than ILUZeroPreconditioner.
ExtendableSparse.ParallelILU0Preconditioner — Type
ParallelILU0Preconditioner()
ParallelILU0Preconditioner(matrix)Parallel ILU preconditioner with zero fill-in.
Iteration schemes
ExtendableSparse.simple! — Function
simple!(u,A,b;
abstol::Real = zero(real(eltype(b))),
reltol::Real = sqrt(eps(real(eltype(b)))),
log=false,
maxiter=100,
P=nothing
) -> solution, [history]Simple iteration scheme $u_{i+1}= u_i - P^{-1} (A u_i -b)$ with similar API as the methods in IterativeSolvers.jl.