mindspore.scipy.sparse.linalg.gmres

mindspore.scipy.sparse.linalg.gmres(A, b, x0=None, *, tol=1e-05, restart=20, maxiter=None, M=None, callback=None, restrt=None, atol=0.0, callback_type=None, solve_method='batched')[源代码]

Given given A and b, GMRES solves the linear system:

\[A x = b\]

A is specified as a function performing A(vi) -> vf = A @ vi, and in principle need not have any particular special properties, such as symmetry. However, convergence is often slow for nearly symmetric operators.

说明

  • gmres is not supported on Windows platform yet.

参数
  • A (Union[Tensor, function]) – 2D Tensor or function that calculates the linear map (matrix-vector product) \(Ax\) when called like \(A(x)\). As function, A must return Tensor with the same structure and shape as its input matrix.

  • b (Tensor) – Right hand side of the linear system representing a single vector. Can be stored as a Tensor.

  • x0 (Tensor, optional) – Starting guess for the solution. Must have the same structure as b. If this is unspecified, zeroes are used. Default: None.

  • tol (float, optional) – Tolerances for convergence, \(norm(residual) <= max(tol*norm(b), atol)\). We do not implement SciPy’s “legacy” behavior, so MindSpore’s tolerance will differ from SciPy unless you explicitly pass atol to SciPy’s gmres. Default: 1e-5.

  • restart (integer, optional) – Size of the Krylov subspace (“number of iterations”) built between restarts. GMRES works by approximating the true solution x as its projection into a Krylov space of this dimension - this parameter therefore bounds the maximum accuracy achievable from any guess solution. Larger values increase both number of iterations and iteration cost, but may be necessary for convergence. The algorithm terminates early if convergence is achieved before the full subspace is built. Default: 20.

  • maxiter (int) – Maximum number of times to rebuild the size-restart Krylov space starting from the solution found at the last iteration. If GMRES halts or is very slow, decreasing this parameter may help. Default: None.

  • M (Union[Tensor, function]) – Preconditioner for A. The preconditioner should approximate the inverse of A. Effective preconditioning dramatically improves the rate of convergence, which implies that fewer iterations are needed to reach a given error tolerance. Default: None.

  • callback (function) – User-supplied function to call after each iteration. It is called as callback(args), where args are selected by callback_type. Default: None.

  • restrt (int, optional) – Deprecated, use restart instead. Default: None.

  • atol (float, optional) – The same as tol. Default: 0.0.

  • callback_type (str, optional) –

    Callback function argument requested: Default: None.

    • x: current iterate (ndarray), called on every restart

    • pr_norm: relative (preconditioned) residual norm (float), called on every inner iteration

    • legacy (default): same as pr_norm, but also changes the meaning of ‘maxiter’ to count inner iterations instead of restart cycles.

  • solve_method (str) –

    There are two kinds of solve methods,’incremental’ or ‘batched’. Default: “batched”.

    • incremental: builds a QR decomposition for the Krylov subspace incrementally during the GMRES process using Givens rotations. This improves numerical stability and gives a free estimate of the residual norm that allows for early termination within a single “restart”.

    • batched: solve the least squares problem from scratch at the end of each GMRES iteration. It does not allow for early termination, but has much less overhead on GPUs.

返回

  • Tensor, the converged solution. Has the same structure as b.

  • Tensor, placeholder for convergence information: 0 : successful exit. >0 : convergence to tolerance not achieved, number of iterations. <0 : illegal input or breakdown.

Supported Platforms:

CPU GPU

样例

>>> import numpy as onp
>>> import mindspore.numpy as mnp
>>> from mindspore.common import Tensor
>>> from mindspore.scipy.sparse.linalg import gmres
>>> A = Tensor(mnp.array([[3, 2, 0], [1, -1, 0], [0, 5, 1]], dtype=mnp.float32))
>>> b = Tensor(mnp.array([2, 4, -1], dtype=mnp.float32))
>>> x, exitCode = gmres(A, b)
>>> print(exitCode)            # 0 indicates successful convergence
0
>>> print(onp.allclose(mnp.dot(A,x).asnumpy(), b.asnumpy()))
True