GridapSolvers.LinearSolvers
Krylov solvers
GridapSolvers.LinearSolvers.CGSolver
— Typestruct CGSolver <: LinearSolver
...
end
CGSolver(Pl;maxiter=1000,atol=1e-12,rtol=1.e-6,flexible=false,verbose=0,name="CG")
Left-Preconditioned Conjugate Gradient solver.
GridapSolvers.LinearSolvers.MINRESSolver
— Typestruct MINRESSolver <: LinearSolver
...
end
MINRESSolver(m;Pl=nothing,maxiter=100,atol=1e-12,rtol=1.e-6,verbose=false,name="MINRES")
MINRES solver, with optional left preconditioners Pl
. The preconditioner must be symmetric and positive definite.
GridapSolvers.LinearSolvers.GMRESSolver
— Typestruct GMRESSolver <: LinearSolver
...
end
GMRESSolver(m;Pr=nothing,Pl=nothing,restart=false,m_add=1,maxiter=100,atol=1e-12,rtol=1.e-6,verbose=false,name="GMRES")
GMRES solver, with optional right and left preconditioners Pr
and Pl
.
The solver starts by allocating a basis of size m
. Then:
- If
restart=true
, the basis size is fixed and restarted everym
iterations. - If
restart=false
, the basis size is allowed to increase. When full, the solver allocatesm_add
new basis vectors.
GridapSolvers.LinearSolvers.FGMRESSolver
— Typestruct FGMRESSolver <: LinearSolver
...
end
FGMRESSolver(m,Pr;Pl=nothing,restart=false,m_add=1,maxiter=100,atol=1e-12,rtol=1.e-6,verbose=false,name="FGMRES")
Flexible GMRES solver, with right-preconditioner Pr
and optional left-preconditioner Pl
.
The solver starts by allocating a basis of size m
. Then:
- If
restart=true
, the basis size is fixed and restarted everym
iterations. - If
restart=false
, the basis size is allowed to increase. When full, the solver allocatesm_add
new basis vectors at a time.
GridapSolvers.LinearSolvers.krylov_mul!
— FunctionComputes the Krylov matrix-vector product
y = Pl⁻¹⋅A⋅Pr⁻¹⋅x
by solving
Pr⋅wr = x wl = A⋅wr Pl⋅y = wl
GridapSolvers.LinearSolvers.krylov_residual!
— FunctionComputes the Krylov residual
r = Pl⁻¹(A⋅x - b)
by solving
w = A⋅x - b Pl⋅r = w
Richardson iterative solver
Given a linear system $Ax = b$ and a left preconditioner $Pl$, this iterative solver performs the following iteration until the solution converges.
\[ x_{k+1} = x_k + \omega Pl^{-1} (b - A x_k)\]
GridapSolvers.LinearSolvers.RichardsonLinearSolver
— Typestruct RichardsonLinearSolver <: LinearSolver
...
end
RichardsonLinearSolver(ω,maxiter;Pl=nothing,rtol=1e-10,atol=1e-6,verbose=true,name = "RichardsonLinearSolver")
Richardson Iteration, with an optional left preconditioners Pl
.
The relaxation parameter (ω) can either be of type Float64 or Vector{Float64}. This gives flexiblity in relaxation.
Smoothers
Given a linear system $Ax = b$, a smoother is an operator S
that takes an iterative solution $x_k$ and its residual $r_k = b - A x_k$, and modifies them in place
\[ S : (x_k,r_k) \rightarrow (x_{k+1},r_{k+1})\]
such that $|r_{k+1}| < |r_k|$.
GridapSolvers.LinearSolvers.RichardsonSmoother
— Typestruct RichardsonSmoother{A} <: LinearSolver
M :: A
niter :: Int64
ω :: Float64
end
Iterative Richardson smoother. Given a solution x
and a residual r
, performs niter
Richardson iterations with damping parameter ω
using the linear solver M
. A Richardson iteration is given by:
dx = ω * inv(M) * r
x = x + dx
r = r - A * dx
Updates both the solution x
and the residual r
in place.
GridapSolvers.LinearSolvers.RichardsonSmoother
— Methodfunction RichardsonSmoother(M::LinearSolver,niter::Int=1,ω::Float64=1.0)
Returns an instance of RichardsonSmoother
from its underlying properties.
Preconditioners
Given a linear system $Ax = b$, a preconditioner is an operator that takes an iterative residual $r_k$ and returns a correction $dx_k$.
GridapSolvers.LinearSolvers.JacobiLinearSolver
— Typestruct JacobiLinearSolver <: LinearSolver end
Given a matrix A
, the Jacobi or Diagonal preconditioner is defined as P = diag(A)
.
GridapSolvers.LinearSolvers.GMGLinearSolverFromMatrices
— Typestruct GMGLinearSolverFromMatrices <: LinearSolver
...
end
Geometric MultiGrid solver, from algebraic parts.
GridapSolvers.LinearSolvers.GMGLinearSolverFromWeakform
— Typestruct GMGLinearSolverFromWeakForm <: LinearSolver
...
end
Geometric MultiGrid solver, from FE parts.
GridapSolvers.LinearSolvers.GMGLinearSolver
— FunctionGMGLinearSolver(
mh::ModelHierarchy,
matrices::AbstractArray{<:AbstractMatrix},
prolongations,
restrictions;
pre_smoothers = Fill(RichardsonSmoother(JacobiLinearSolver(),10),num_levels(mh)-1),
post_smoothers = pre_smoothers,
coarsest_solver = LUSolver(),
mode::Symbol = :preconditioner,
maxiter = 100, atol = 1.0e-14, rtol = 1.0e-08, verbose = false,
)
Creates an instance of GMGLinearSolverFromMatrices
from the underlying model hierarchy, the system matrices at each level and the transfer operators and smoothers at each level except the coarsest.
The solver has two modes of operation, defined by the kwarg mode
:
:solver
: The GMG solver takes a rhsb
and returns a solutionx
.:preconditioner
: The GMG solver takes a residualr
and returns a correctiondx
.
GMGLinearSolver(
mh::ModelHierarchy,
trials::FESpaceHierarchy,
tests::FESpaceHierarchy,
biforms::AbstractArray{<:Function},
interp,
restrict;
pre_smoothers = Fill(RichardsonSmoother(JacobiLinearSolver(),10),num_levels(mh)-1),
post_smoothers = pre_smoothers,
coarsest_solver = Gridap.Algebra.LUSolver(),
mode::Symbol = :preconditioner,
is_nonlinear = false,
maxiter = 100, atol = 1.0e-14, rtol = 1.0e-08, verbose = false,
)
Creates an instance of GMGLinearSolverFromMatrices
from the underlying model hierarchy, the trial and test FEspace hierarchies, the weakform lhs at each level and the transfer operators and smoothers at each level except the coarsest.
The solver has two modes of operation, defined by the kwarg mode
:
:solver
: The GMG solver takes a rhsb
and returns a solutionx
.:preconditioner
: The GMG solver takes a residualr
and returns a correctiondx
.
Nullspaces
GridapSolvers.LinearSolvers.NullspaceSolver
— Typestruct NullspaceSolver{A,B}
solver :: A <: LinearSolver
kernel :: B <: Union{AbstractMatrix, Vector{<:AbstractVector}}
constrain_matrix :: Bool = true
end
Solver that computes the solution of a linear system A⋅x = b
with a kernel constraint, i.e the returned solution is orthogonal to the provided kernel.
We assume the kernel is provided as a matrix K
of dimensions (k,n)
with n
the dimension of the original system and k
the number of kernel vectors.
Two modes of operation are supported:
If
constrain_matrix
istrue
, the solver will explicitly constrain the system matrix. I.e we consider the augmented systemÂ⋅x̂ = b̂
whereAk = [A, K'; K, 0]
x̂ = [x; λ]
b̂ = [b; 0]
This is often the only option for direct solvers, which require the system matrix to be invertible. This is only supported in serial (its performance bottleneck in parallel).
If
constrain_matrix
isfalse
, the solver preserve the original system and simply project the initial guess and the solution onto the orthogonal complement of the kernel. This option is more suitable for iterative solvers, which usually do not require the system matrix to be invertible (e.g. GMRES, BiCGStab, etc).
Wrappers
PETSc
Building on top of GridapPETSc.jl, GridapSolvers provides specific solvers for some particularly complex PDEs:
GridapSolvers.LinearSolvers.ElasticitySolver
— Typestruct ElasticitySolver <: LinearSolver
...
end
GMRES + AMG solver, specifically designed for linear elasticity problems.
Follows PETSc's documentation for PCAMG and MatNullSpaceCreateRigidBody.
GridapSolvers.LinearSolvers.ElasticitySolver
— Methodfunction ElasticitySolver(space::FESpace; maxiter=500, atol=1.e-12, rtol=1.e-8)
Returns an instance of ElasticitySolver
from its underlying properties.
GridapSolvers.LinearSolvers.CachedPETScNS
— Typestruct CachedPETScNS <: NumericalSetup
...
end
Wrapper around a PETSc NumericalSetup, providing highly efficiend reusable caches:
When converting julia vectors/PVectors to PETSc vectors, we purposely create aliasing of the vector values. This means we can avoid copying data from one to another before solving, but we need to be careful about it.
This structure takes care of this, and makes sure you do not attempt to solve the system with julia vectors that are not the ones you used to create the solver cache.
GridapSolvers.LinearSolvers.CachedPETScNS
— Methodfunction CachedPETScNS(ns::PETScLinearSolverNS,x::AbstractVector,b::AbstractVector)
Create a new instance of CachedPETScNS
from its underlying properties. Once this structure is created, you can only solve the system with the same vectors you used to create it.
GridapSolvers.LinearSolvers.get_dof_coordinates
— Functionget_dof_coordinates(space::FESpace)
Given a lagrangian FESpace, returns the physical coordinates of the DoFs, as required by some PETSc solvers. See PETSc documentation.
IterativeSolvers.jl
GridapSolvers provides wrappers for some iterative solvers from the package IterativeSolvers.jl:
GridapSolvers.LinearSolvers.IterativeLinearSolver
— Typestruct IterativeLinearSolver <: LinearSolver
...
end
Wrappers for IterativeSolvers.jl krylov-like iterative solvers.
All wrappers take the same kwargs as the corresponding solver in IterativeSolvers.jl.
The following solvers are available:
GridapSolvers.LinearSolvers.IS_ConjugateGradientSolver
— FunctionIS_ConjugateGradientSolver(;kwargs...)
Wrapper for the Conjugate Gradient solver.
GridapSolvers.LinearSolvers.IS_GMRESSolver
— FunctionIS_GMRESSolver(;kwargs...)
Wrapper for the GMRES solver.
GridapSolvers.LinearSolvers.IS_MINRESSolver
— FunctionIS_MINRESSolver(;kwargs...)
Wrapper for the MINRES solver.
GridapSolvers.LinearSolvers.IS_SSORSolver
— FunctionIS_SSORSolver(ω;kwargs...)
Wrapper for the SSOR solver.