Solvers

Standalone

AMGCLWrap.AMGSolverType
AMGSolver(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 or SparseMatricesCSR.SparseMatrixCSR.
  • blocksize: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix of blocksize 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 by JSON3.write. If params is an emtpy string or nothing a default value is used.
source
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 or SparseMatricesCSR.SparseMatrixCSR.

  • blocksize: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix of blocksize x blocksize static matrices. Block sizes 1...8 are instantiated.

  • verbose: if true, print generated JSON string passed to amgcl.

  • param: Ignored if nothing (default). Otherwise, any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string by JSON3.write.

  • coarsening: One of the Coarsening strategies

  • relax: One of the Relaxation strategies

  • solver: One of the Iterative solver strategies

source
AMGCLWrap.RLXSolverType
RLXSolver(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 or SparseMatricesCSR.SparseMatrixCSR.
  • blocksize: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix of blocksize 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 by JSON3.write. If params is an emtpy string or nothing a default value is used.
source
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 or SparseMatricesCSR.SparseMatrixCSR.

  • blocksize: If blocksize >1, group unknowns into blocks of given size and cast the matrix internally to a sparse matrix of blocksize x blocksize static matrices. Block sizes 1...8 are instantiated.

  • verbose: if true, print generated JSON string passed to amgcl.

  • param: Ignored if nothing (default). Otherwise, any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string by JSON3.write.

  • precond: One of the Relaxation strategies seen as preconditioner

  • solver: One of the Iterative solver strategies

source

LinearSolve algorithms

AMGCLWrap.AMGSolverAlgorithmFunction
AMGSolverAlgorithm(;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 of blocksize x blocksize static matrices. Block sizes 1...8 are instantiated.

  • verbose: if true, print generated JSON string passed to amgcl.

  • param: Ignored if nothing (default). Otherwise, any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string by JSON3.write.

  • coarsening: One of the Coarsening strategies

  • relax: One of the Relaxation strategies

  • solver: One of the Iterative solver strategies

source
AMGCLWrap.RLXSolverAlgorithmFunction
RLXSolverAlgorithm(;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 of blocksize x blocksize static matrices. Block sizes 1...8 are instantiated.

  • verbose: if true, print generated JSON string passed to amgcl.

  • param: Ignored if nothing (default). Otherwise, any object (e.g. Tuple, Dict or JSON string) which can be turned into a JSON string by JSON3.write.

  • precond: One of the Relaxation strategies seen as preconditioner

  • solver: One of the Iterative solver strategies

source