Preconditioners
AMGCLWrap.AMGPrecon
— TypeAMGPrecon(sparsematrix::AbstractSparseMatrix;
blocksize=1,
param=nothing)
Create algebraic multigrid preconditioner with ldiv!
and \
methods solving the preconditioning system.
Parameters:
sparsematrix
:SparseArrays.AbstractSparseMatrixCSC
orSparseMatricesCSR.SparseMatrixCSR
.blocksize
: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix ofblocksize x blocksize
static matrices. Block sizes 1...8 are instantiated.param
: Any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string byJSON3.write
. Ifparams
is an emtpy string ornothing
a default value is used.
AMGPrecon(sparsematrix::AbstractSparseMatrix; blocksize::Int=1, param=nothing, verbose::Bool=false, coarsening::Union{AbstractCoarsening, NamedTuple}=SmoothedAggregationCoarsening(), relax::Union{AbstractRelaxation, NamedTuple}=SPAI0Relaxation()) Create algebraic multigrid preconditioner with ldiv!
and \
methods solving the preconditioning system.
Parameters:
sparsematrix
:SparseArrays.AbstractSparseMatrixCSC
orSparseMatricesCSR.SparseMatrixCSR
.blocksize
: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix ofblocksize x blocksize
static matrices. Block sizes 1...8 are instantiated.verbose
: if true, print generated JSON string passed to amgcl.param:
Ignored ifnothing
(default). Otherwise, any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string byJSON3.write
.coarsening
: A coarsening strategyrelax
: A relaxation method
AMGCLWrap.AMGPreconBuilder
— Typestruct AMGPreconBuilder
Preconditioner strategy (e.g. for the new precs
kwarg in LinearSolve) for interfacing AMGPrecon
.
Fields (for documentation, see AMGPrecon
):
param::Any
verbose::Bool
blocksize::Int64
coarsening::Union{AMGCLWrap.AbstractCoarsening, NamedTuple}
relax::Union{AMGCLWrap.AbstractRelaxation, NamedTuple}
AMGCLWrap.RLXPrecon
— TypeRLXPrecon(sparsematrix::AbstractSparseMatrix;
blocksize=1,
param=nothing)
Create single level relaxation preconditioner with ldiv!
and \
methods solving the preconditioning system.
Parameters:
sparsematrix
:SparseArrays.AbstractSparseMatrixCSC
orSparseMatricesCSR.SparseMatrixCSR
.blocksize
: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix ofblocksize x blocksize
static matrices. Block sizes 1...8 are instantiated.param
: Any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string byJSON3.write
. Ifparams
is an emtpy string ornothing
a default value is used.
RLXPrecon(sparsematrix::AbstractSparseMatrix;
blocksize::Int=1,
param=nothing,
verbose::Bool=false,
precond::Union{AbstractRelaxation, NamedTuple}=ILU0Relaxation())
Create single level relaxation preconditioner with ldiv!
and \
methods solving the preconditioning system.
Parameters:
sparsematrix
:SparseArrays.AbstractSparseMatrixCSC
orSparseMatricesCSR.SparseMatrixCSR
.blocksize
: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix ofblocksize x blocksize
static matrices. Block sizes 1...8 are instantiated.verbose
: if true, print generated JSON string passed to amgcl.param:
Ignored ifnothing
(default). Otherwise, any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string byJSON3.write
.precond
: A preconditioning method
AMGCLWrap.RLXPreconBuilder
— Typestruct RLXPreconBuilder
Preconditioner strategy (e.g. for the new precs
kwarg in LinearSolve) for interfacing RLXPrecon
.
Fields (for documentation, see RLXPrecon
):
param::Any
verbose::Bool
blocksize::Int64
precond::Union{AMGCLWrap.AbstractRelaxation, NamedTuple}