pyamg.multilevel#

Generic AMG solver.

Functions

coarse_grid_solver(solver)

Return a coarse grid solver suitable for MultilevelSolver.

Classes

MultilevelSolver(levels[, coarse_solver])

Stores multigrid hierarchy and implements the multigrid cycle.

multilevel_solver(*args, **kwargs)

Deprecated level class.

class pyamg.multilevel.MultilevelSolver(levels, coarse_solver='pinv')[source]#

Stores multigrid hierarchy and implements the multigrid cycle.

The class constructs the cycling process and points to the methods for coarse grid solves. A MultilevelSolver object is typically returned from a particular AMG method (see ruge_stuben_solver or smoothed_aggregation_solver for example). A call to MultilevelSolver.solve() is a typical access point. The class also defines methods for constructing operator, cycle, and grid complexities.

Attributes:
levelslevel array

Array of level objects that contain A, R, and P.

coarse_solverstring

String passed to coarse_grid_solver indicating the solve type

Methods

aspreconditioner()

Create a preconditioner using this multigrid cycle

cycle_complexity()

A measure of the cost of a single multigrid cycle.

grid_complexity()

A measure of the rate of coarsening.

operator_complexity()

A measure of the size of the multigrid hierarchy.

solve()

Iteratively solves a linear system for the right hand side.

change_solve_matrix(A)

Change matrix solve/preconditioning matrix. This also changes the corresponding relaxation routines on the fine grid. This can be used, for example, to precondition a quadratic finite element discretization with AMG built from a linear discretization on quadratic quadrature points.

class Level[source]#

Stores one level of the multigrid hierarchy.

All level objects will have an ‘A’ attribute referencing the matrix of that level. All levels, except for the coarsest level, will also have ‘P’ and ‘R’ attributes referencing the prolongation and restriction operators that act between each level and the next coarser level.

Notes

The functionality of this class is a struct

Attributes:
Acsr_matrix

Problem matrix for Ax=b

Rcsr_matrix

Restriction matrix between levels (often R = P.T)

Pcsr_matrix

Prolongation or Interpolation matrix.

__init__()[source]#

Level construct (empty).

__init__(levels, coarse_solver='pinv')[source]#

Class constructor to initialize the cycle and ensure list of levels is complete.

Parameters:
levelslevel array

Array of level objects that contain A, R, and P.

coarse_solver: string, callable, tuple

The solver method is either (1) a string such as ‘splu’ or ‘pinv’ of a callable object which receives only parameters (A, b) and returns an (approximate or exact) solution to the linear system Ax = b, or (2) a callable object that takes parameters (A,b) and returns an (approximate or exact) solution to Ax = b, or (3) a tuple of the form (string|callable, args), where args is a dictionary of arguments to be passed to the function denoted by string or callable.

Sparse direct methods:

  • splu : sparse LU solver

Sparse iterative methods:

  • any method in scipy.sparse.linalg or pyamg.krylov (e.g. ‘cg’).

  • Methods in pyamg.krylov take precedence.

  • relaxation method, such as ‘gauss_seidel’ or ‘jacobi’,

Dense methods:

  • pinv : pseudoinverse (SVD)

  • lu : LU factorization

  • cholesky : Cholesky factorization

Notes

If not defined, the R attribute on each level is set to the transpose of P.

Examples

>>> # manual construction of a two-level AMG hierarchy
>>> from pyamg.gallery import poisson
>>> from pyamg.multilevel import MultilevelSolver
>>> from pyamg.strength import classical_strength_of_connection
>>> from pyamg.classical.interpolate import direct_interpolation
>>> from pyamg.classical.split import RS
>>> # compute necessary operators
>>> A = poisson((100, 100), format='csr')
>>> C = classical_strength_of_connection(A)
>>> splitting = RS(A)
>>> P = direct_interpolation(A, C, splitting)
>>> R = P.T
>>> # store first level data
>>> levels = []
>>> levels.append(MultilevelSolver.Level())
>>> levels.append(MultilevelSolver.Level())
>>> levels[0].A = A
>>> levels[0].C = C
>>> levels[0].splitting = splitting
>>> levels[0].P = P
>>> levels[0].R = R
>>> # store second level data
>>> levels[1].A = R @ A @ P                      # coarse-level matrix
>>> # create MultilevelSolver
>>> ml = MultilevelSolver(levels, coarse_solver='splu')
>>> print(ml)
MultilevelSolver
Number of Levels:     2
Operator Complexity:  1.891
Grid Complexity:      1.500
Coarse Solver:        'splu'
  level   unknowns     nonzeros
     0       10000        49600 [52.88%]
     1        5000        44202 [47.12%]
__repr__()[source]#

Print basic statistics about the multigrid hierarchy.

__solve(lvl, x, b, cycle, cycles_per_level=1)#

Multigrid cycling.

Parameters:
lvlint

Solve problem on level lvl

xnumpy array

Initial guess x and return correction

bnumpy array

Right-hand side for Ax=b

cycle{‘V’,’W’,’F’,’AMLI’}

Recursively called cycling function. The Defines the cycling used: cycle = ‘V’, V-cycle cycle = ‘W’, W-cycle cycle = ‘F’, F-cycle cycle = ‘AMLI’, AMLI-cycle

cycles_per_levelint, default 1

Number of V-cycles on each level of an F-cycle

aspreconditioner(cycle='V')[source]#

Create a preconditioner using this multigrid cycle.

Parameters:
cycle{‘V’,’W’,’F’,’AMLI’}

Type of multigrid cycle to perform in each iteration.

Returns:
precondLinearOperator

Preconditioner suitable for the iterative solvers in defined in the scipy.sparse.linalg module (e.g. cg, gmres) and any other solver that uses the LinearOperator interface. Refer to the LinearOperator documentation in scipy.sparse.linalg

See also

MultilevelSolver.solve, scipy.sparse.linalg.LinearOperator

Examples

>>> from pyamg.aggregation import smoothed_aggregation_solver
>>> from pyamg.gallery import poisson
>>> from scipy.sparse.linalg import cg
>>> import scipy as sp
>>> A = poisson((100, 100), format='csr')          # matrix
>>> b = np.random.rand(A.shape[0])                 # random RHS
>>> ml = smoothed_aggregation_solver(A)            # AMG solver
>>> M = ml.aspreconditioner(cycle='V')             # preconditioner
>>> x, info = cg(A, b, tol=1e-8, maxiter=30, M=M)  # solve with CG
change_solve_matrix(A)[source]#

Change matrix solve/preconditioning matrix.

Parameters:
Acsr_matrix

Target solution matrix

Notes

This also changes the corresponding relaxation routines on the fine grid. This can be used, for example, to precondition a quadratic finite element discretization with linears.

cycle_complexity(cycle='V')[source]#

Cycle complexity of V, W, AMLI, and F(1,1) cycle with simple relaxation.

Cycle complexity is an approximate measure of the number of floating point operations (FLOPs) required to perform a single multigrid cycle relative to the cost a single smoothing operation.

Parameters:
cycle{‘V’,’W’,’F’,’AMLI’}

Type of multigrid cycle to perform in each iteration.

Returns:
ccfloat

Defined as F_sum / F_0, where F_sum is the total number of nonzeros in the matrix on all levels encountered during a cycle and F_0 is the number of nonzeros in the matrix on the finest level.

Notes

This is only a rough estimate of the true cycle complexity. The estimate assumes that the cost of pre and post-smoothing are (each) equal to the number of nonzeros in the matrix on that level. This assumption holds for smoothers like Jacobi and Gauss-Seidel. However, the true cycle complexity of cycle using more expensive methods, like block Gauss-Seidel will be underestimated.

Additionally, if the cycle used in practice isn’t a (1,1)-cycle, then this cost estimate will be off.

grid_complexity()[source]#

Grid complexity of this multigrid hierarchy.

Defined as:

Number of unknowns on all levels / Number of unknowns on the finest level

class level[source]#

Deprecated level class.

__init__()[source]#

Raise deprecation warning on use, not import.

operator_complexity()[source]#

Operator complexity of this multigrid hierarchy.

Defined as:

Number of nonzeros in the matrix on all levels / Number of nonzeros in the matrix on the finest level

psolve(b)[source]#

Legacy solve interface.

solve(b, x0=None, tol=1e-05, maxiter=100, cycle='V', accel=None, callback=None, residuals=None, cycles_per_level=1, return_info=False)[source]#

Execute multigrid cycling.

Parameters:
barray

Right hand side.

x0array

Initial guess.

tolfloat

Stopping criteria: relative residual r[k]/||b|| tolerance. If accel is used, the stopping criteria is set by the Krylov method.

maxiterint

Stopping criteria: maximum number of allowable iterations.

cycle{‘V’,’W’,’F’,’AMLI’}

Type of multigrid cycle to perform in each iteration.

accelstring, function

Defines acceleration method. Can be a string such as ‘cg’ or ‘gmres’ which is the name of an iterative solver in pyamg.krylov (preferred) or scipy.sparse.linalg. If accel is not a string, it will be treated like a function with the same interface provided by the iterative solvers in SciPy.

callbackfunction

User-defined function called after each iteration. It is called as callback(xk) where xk is the k-th iterate vector.

residualslist

List to contain residual norms at each iteration. The residuals will be the residuals from the Krylov iteration – see the accel method to see verify whether this ||r|| or ||Mr|| (as in the case of GMRES).

cycles_per_level: int, default 1

Number of V-cycles on each level of an F-cycle

return_infobool

If true, will return (x, info) If false, will return x (default)

Returns:
xarray

Approximate solution to Ax=b after k iterations

infostring

Halting status

0

successful exit

>0

convergence to tolerance not achieved, return iteration count instead.

See also

aspreconditioner

Examples

>>> from numpy import ones
>>> from pyamg import ruge_stuben_solver
>>> from pyamg.gallery import poisson
>>> A = poisson((100, 100), format='csr')
>>> b = A * ones(A.shape[0])
>>> ml = ruge_stuben_solver(A, max_coarse=10)
>>> residuals = []
>>> x = ml.solve(b, tol=1e-12, residuals=residuals) # standalone solver
pyamg.multilevel.coarse_grid_solver(solver)[source]#

Return a coarse grid solver suitable for MultilevelSolver.

Parameters:
solverstring, callable, tuple

The solver method is either (1) a string such as ‘splu’ or ‘pinv’ of a callable object which receives only parameters (A, b) and returns an (approximate or exact) solution to the linear system Ax = b, or (2) a callable object that takes parameters (A,b) and returns an (approximate or exact) solution to Ax = b, or (3) a tuple of the form (string|callable, args), where args is a dictionary of arguments to be passed to the function denoted by string or callable.

The set of valid string arguments is:
  • Sparse direct methods:
    • splu : sparse LU solver

  • Sparse iterative methods:
    • the name of any method in scipy.sparse.linalg or pyamg.krylov (e.g. ‘cg’). Methods in pyamg.krylov take precedence.

    • relaxation method, such as ‘gauss_seidel’ or ‘jacobi’, present in pyamg.relaxation

  • Dense methods:
    • pinv : pseudoinverse (SVD)

    • lu : LU factorization

    • cholesky : Cholesky factorization

Returns:
ptrGenericSolver

A class for use as a standalone or coarse grids solver

Examples

>>> import numpy as np
>>> from scipy.sparse import spdiags
>>> from pyamg.gallery import poisson
>>> from pyamg import coarse_grid_solver
>>> A = poisson((10, 10), format='csr')
>>> b = A * np.ones(A.shape[0])
>>> cgs = coarse_grid_solver('lu')
>>> x = cgs(A, b)
class pyamg.multilevel.multilevel_solver(*args, **kwargs)[source]#

Deprecated level class.

Deprecated since version 4.2.3: Use MultilevelSolver instead.

__init__(*args, **kwargs)[source]#

Raise deprecation warning on use, not import.