AMGCLWrap.AMGCLWrap
— ModuleAMGCLWrap
AMGCLWrap
Julia wrapper around a subset of AMGCL via the C wrapper AMGCL_C.
If you find this package useful, please give credits to the author of AMGCL.
AMGCLWrap.AMGPrecon
— MethodAMGPrecon(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
.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.AMGPrecon
— MethodAMGPrecon(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.
AMGCLWrap.AMGSolver
— MethodAMGSolver(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
.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.AMGSolver
— MethodAMGSolver(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.
AMGCLWrap.AbstractAMGCLParams
— Typeabstract type AbstractAMGCLParams
Abstract parameter type.
AMGCLWrap.AbstractCoarsening
— Typeabstract type AbstractCoarsening <: AMGCLWrap.AbstractAMGCLParams
Abstract coarsening parameter type.
AMGCLWrap.AbstractRelaxation
— Typeabstract type AbstractRelaxation <: AMGCLWrap.AbstractAMGCLParams
Abstract relaxation parameter type.
AMGCLWrap.AbstractSolver
— Typeabstract type AbstractSolver <: AMGCLWrap.AbstractAMGCLParams
Abstract solver parameter type.
AMGCLWrap.BICGStabSolver
— Typestruct BICGStabSolver <: AMGCLWrap.AbstractSolver
BICGStab solver
type::String
pside::String
tol::Float64
abstol::Float64
maxiter::Int64
verbose::Bool
AMGCLWrap.CGSolver
— Typestruct CGSolver <: AMGCLWrap.AbstractSolver
type::String
tol::Float64
abstol::Float64
maxiter::Int64
verbose::Bool
AMGCLWrap.GMRESSolver
— Typestruct GMRESSolver <: AMGCLWrap.AbstractSolver
type::String
M::Int64
pside::String
tol::Float64
abstol::Float64
maxiter::Int64
verbose::Bool
AMGCLWrap.ILU0Relaxation
— Typestruct ILU0Relaxation <: AMGCLWrap.AbstractRelaxation
type::String
AMGCLWrap.RLXPrecon
— MethodRLXPrecon(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.
AMGCLWrap.RLXSolver
— MethodRLXSolver(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
.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
AMGCLWrap.RLXSolver
— MethodRLXSolver(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.
AMGCLWrap.RugeStubenCoarsening
— Typestruct RugeStubenCoarsening <: AMGCLWrap.AbstractCoarsening
type::String
eps_strong::Float64
do_trunc::Bool
eps_trunc::Float64
AMGCLWrap.SPAI0Relaxation
— Typestruct SPAI0Relaxation <: AMGCLWrap.AbstractRelaxation
type::String
AMGCLWrap.SmoothedAggregationCoarsening
— Typestruct SmoothedAggregationCoarsening <: AMGCLWrap.AbstractCoarsening
type::String
relax::Float64
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:
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
Only available for Julia version >=1.9
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:
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
Only available for Julia version >=1.9
AMGCLWrap.blocksize_instantiated
— Methodblocksize_instantiated(blocksize)::Bool
Check if given blocksize has been instantiated.
AMGCLWrap.error_state
— Methoderror_state(operator)::Int
Return error state of operator. Ok if 0.