KSP
The KSP (Krylov Subspace Methods) module provides iterative linear solvers for systems of the form Ax = b. PETSc offers a wide variety of Krylov methods and preconditioners.
Overview
KSP provides:
- Krylov methods: GMRES, CG, BiCGStab, and many more
- Preconditioners: Jacobi, ILU, multigrid, direct solvers, etc.
- Runtime configuration: Choose methods via command-line options
- Convergence monitoring: Built-in residual tracking
Creating a KSP Solver
From a Matrix
# Basic creation with default options
ksp = KSP(A)
# With preconditioner matrix P (for different preconditioning)
ksp = KSP(A, P)
# With options
ksp = KSP(A;
ksp_type = "gmres",
pc_type = "ilu",
ksp_rtol = 1e-8
)From a DM
# Create KSP associated with a DM (for multigrid, etc.)
ksp = KSP(dm;
ksp_type = "cg",
pc_type = "mg"
)From a Sparse Matrix
# Directly from Julia SparseMatrixCSC
using SparseArrays
S = sprand(100, 100, 0.1) + 10I
ksp = KSP(petsclib, MPI.COMM_SELF, S)Solving
# Solve Ax = b
solve!(x, ksp, b)
# Or allocate solution vector
x = solve(ksp, b)Common Solver/Preconditioner Options
Krylov Methods (ksp_type)
cg- Conjugate Gradient (symmetric positive definite)gmres- Generalized Minimum Residualbicgstab- BiConjugate Gradient Stabilizedrichardson- Richardson iterationpreonly- Apply preconditioner only (for direct solvers)
Preconditioners (pc_type)
jacobi- Diagonal scalingilu- Incomplete LU factorizationlu- Direct LU factorizationmg- Geometric multigridgamg- Algebraic multigridnone- No preconditioning
Convergence Options
ksp_rtol- Relative tolerance (default: 1e-5)ksp_atol- Absolute toleranceksp_max_it- Maximum iterationsksp_monitor- Print residual each iteration
Example: Multigrid Solver
ksp = KSP(dm;
ksp_type = "cg",
pc_type = "mg",
pc_mg_levels = 4,
pc_mg_galerkin = true,
mg_levels_ksp_type = "richardson",
mg_levels_pc_type = "jacobi",
mg_coarse_pc_type = "lu"
)Functions
PETSc.KSP — Method
KSP(dm::AbstractPetscDM; prefix="", options...)Create a KSP associated with the dm with optional prefix and options.
The communicator is obtained from dm. The KSP can be used with geometric multigrid when the DM provides grid hierarchy information.
Arguments
dm::AbstractPetscDM: The DM object to associate with the KSPprefix::String: Optional prefix for command-line optionsoptions...: Additional PETSc options as keyword arguments
External Links
- PETSc Manual:
KSP/KSPCreate
- PETSc Manual:
KSP/KSPSetDM
- PETSc Manual:
KSP/KSPSetFromOptions
PETSc.KSP — Method
KSP(comm::MPI.Comm, A::AbstractPetscMat, P::AbstractPetscMat{PetscLib} = A; prefix="", options...)Create a KSP using the matrix A and preconditioner construction matrix P with optional prefix and options.
The communicator is obtained from A and if it has size 1 then the garbage collector is set, otherwise the user is responsible for calling destroy.
External Links
- PETSc Manual:
KSP/KSPCreate
- PETSc Manual:
KSP/KSPSetOperators
- PETSc Manual:
KSP/KSPSetFromOptions
PETSc.getDMDA — Method
getDMDA(ksp::AbstractKSP)Get dmda for ksp
The returned dmda is owned by the ksp
External Links
- PETSc Manual:
KSP/KSPGetDM
PETSc.get_solution — Method
sol = get_solution(ksp::AbstractPetscKSP)Returns the soluteion vector associated with the KSP object.
PETSc.setcomputeoperators! — Method
setcomputeoperators!(ksp::PetscKSP, ops!::Function)
setcomputeoperators!(ops!::Function, ksp::PetscKSP)Define ops! to be the compute operators function for the ksp. A call to ops!(A, P, new_ksp) should set the elements of the PETSc matrix linear operator A and preconditioning matrix P based on the new_ksp.
External Links
- PETSc Manual:
KSP/KSPSetComputeOperators
PETSc.setcomputerhs! — Method
setcomputerhs!(ksp::AbstractKSP, rhs!::Function)
setcomputerhs!(rhs!::Function, ksp::AbstractKSP)Define rhs! to be the right-hand side function of the ksp. A call to rhs!(b, new_ksp) should set the elements of the PETSc vector b based on the new_ksp.
External Links
- PETSc Manual:
KSP/KSPSetComputeRHS