Solvers
Standalone
AMGCLWrap.AMGSolver
— TypeAMGSolver(sparsematrix::AbstractSparseMatrix;
blocksize=1,
param=nothing)
Create Algebraic multigrid preconditioned Krylov subspace solver with ldiv!
and \
methods solving the matrix 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.
AMGSolver(sparsematrix::AbstractSparseMatrix;
blocksize::Int=1,
param=nothing,
verbose::Bool=false,
coarsening::Union{AbstractCoarsening, NamedTuple}=SmoothedAggregationCoarsening(),
relax::Union{AbstractRelaxation, NamedTuple}=SPAI0Relaxation(),
solver::Union{AbstractSolver, NamedTuple}=BICGStabSolver(;verbose))
Create Algebraic multigrid preconditioned Krylov subspace solver with ldiv!
and \
methods solving the matrix 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
: One of the Coarsening strategiesrelax
: One of the Relaxation strategiessolver
: One of the Iterative solver strategies
AMGCLWrap.RLXSolver
— TypeRLXSolver(sparsematrix::AbstractSparseMatrix;
blocksize=1,
param=nothing)
Create single level relaxation preconditioned Krylov subspace solver with ldiv!
and \
methods solving the matrix 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.
RLXSolver(sparsematrix::AbstractSparseMatrix;
blocksize::Int=1,
param=nothing,
verbose::Bool=false,
precond::Union{AbstractRelaxation, NamedTuple}=ILU0Relaxation(),
solver::Union{AbstractSolver, NamedTuple}=BICGStabSolver(;verbose))
Create single level relaxation preconditioned Krylov subspace solver with ldiv!
and \
methods solving the matrix 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
: One of the Relaxation strategies seen as preconditionersolver
: One of the Iterative solver strategies
LinearSolve algorithms
AMGCLWrap.AMGSolverAlgorithm
— FunctionAMGSolverAlgorithm(;blocksize::Int=1,
param=nothing,
verbose::Bool=false,
coarsening::Union{AbstractCoarsening, NamedTuple}=SmoothedAggregationCoarsening(),
relax::Union{AbstractRelaxation, NamedTuple}=SPAI0Relaxation(),
solver::Union{AbstractSolver, NamedTuple}=BICGStabSolver(;verbose))
Algebraic multigrid preconditioned Krylov subspace solver algorithm for LinearSolve.jk Parameters:
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
: One of the Coarsening strategiesrelax
: One of the Relaxation strategiessolver
: One of the Iterative solver strategies
AMGCLWrap.RLXSolverAlgorithm
— FunctionRLXSolverAlgorithm(;blocksize::Int=1,
param=nothing,
verbose::Bool=false,
precond::Union{AbstractRelaxation, NamedTuple}=ILU0Relaxation(),
solver::Union{AbstractSolver, NamedTuple}=BICGStabSolver(;verbose))
Algebraic multigrid preconditioned Krylov subspace solver algorithm for LinearSolve.jk Parameters:
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
: One of the Relaxation strategies seen as preconditionersolver
: One of the Iterative solver strategies