Skip to content

openseespy_solvers.cupy

CUDA sparse linear and eigen solver constructors implemented with cupyx.scipy.sparse.linalg.

This module follows cupyx.scipy.sparse.linalg naming where possible. Constructors return OpenSeesPy-compatible solver objects; OpenSeesPy supplies A, b, K, and M at solve time. Importing this module requires cupy.

Solving Linear Problems

Direct methods:

Constructor cupy analogue Description
spsolve cupyx.scipy.sparse.linalg.splu Sparse direct solver on CUDA

Iterative methods:

Constructor cupy analogue Description
cg cupyx.scipy.sparse.linalg.cg Conjugate Gradient on CUDA
gmres cupyx.scipy.sparse.linalg.gmres GMRES on CUDA

Eigenvalue Problems

Constructor cupy/scipy analogue Description
eigsh cupyx.scipy.sparse.linalg.eigsh / scipy.sparse.linalg.eigsh CUDA-assisted generalized symmetric eigen solve
lobpcg cupyx.scipy.sparse.linalg.lobpcg LOBPCG on CUDA

Installation

Install the CUDA extra matching your driver:

python -m pip install "openseespy-solvers[cuda13]"   # or [cuda12]

OpenSeesPy Usage

from openseespy_solvers.cupy import eigsh

solver = eigsh(tol=1e-8)
lam = ops.eigen("PythonSparse", num_modes, solver.to_openseespy())

Notes

OpenSeesPy assembles matrices on the CPU. cupy solvers copy matrix data to the GPU before solving, so GPU speedups depend on problem size, sparsity, and transfer overhead.

eigsh solves K x = lambda M x with shift-invert support:

  • mass_mode="general" uses scipy.sparse.linalg.eigsh with GPU inner solves when shift-invert is active.
  • mass_mode="diagonal" and mass_mode="lumped" use GPU shift-invert with diagonal mass.

Raw cupyx.scipy.sparse.linalg.eigsh does not accept a mass matrix; use lobpcg when you need a different iterative eigen strategy.

Function Reference

Sparse linear algebra solvers for OpenSeesPy (cupyx.scipy.sparse.linalg backend).

This module provides solver objects that wrap :mod:cupyx.scipy.sparse.linalg for OpenSeesPy's PythonSparse commands. Factory signatures match the cupyx.scipy.sparse.linalg API except that A and b (or K and M) are supplied by OpenSees at solve time.

Importing this module requires cupy. Install a CUDA-matched extra, for example python -m pip install "openseespy-solvers[cuda13]" or use [cuda12].

Notes

Generalized :func:cupyx.scipy.sparse.linalg.eigsh is not available. :func:eigsh supports:

  • mass_mode='general' (default): matches OpenSees eigsh (shift-invert with sigma=0 for smallest modes; plain ARPACK for largest). GPU inner solves when shift-invert is used.
  • mass_mode='diagonal' / 'lumped': shift-invert on GPU (smallest modes, or when sigma is set); use general for largest modes without a shift.

Use :func:lobpcg when these modes are not appropriate.

Submodules

precond GPU preconditioner factories for use with M=.

See Also

cupyx.scipy.sparse.linalg openseespy_solvers.scipy CPU backend with the same interface.

spsolve

spsolve(*, permc_spec: str = 'COLAMD', scheme: str | None = None, writable: str | list[str] = 'none', debug: bool = False, dtype: Any = float64) -> _SpSolve
Source code in src/openseespy_solvers/cupy/__init__.py
def spsolve(
    *,
    permc_spec: str = "COLAMD",
    scheme: str | None = None,
    writable: str | list[str] = "none",
    debug: bool = False,
    dtype: Any = np.float64,
) -> _SpSolve:
    r"""Configure a sparse direct solver for OpenSees ``PythonSparse`` (GPU).

    Uses :func:`cupyx.scipy.sparse.linalg.splu` internally. The LU factorization
    is reused when OpenSees reports ``matrix_status='UNCHANGED'``.

    Parameters
    ----------
    permc_spec : str, optional
        Column permutation applied during LU factorization. Passed to
        :func:`cupyx.scipy.sparse.linalg.splu`. Default is ``'COLAMD'``.
    """ + _OPENSEES_LINEAR + _LINEAR_RETURNS + _LINEAR_NOTES + """
    See Also
    --------
    cupyx.scipy.sparse.linalg.splu
    cupyx.scipy.sparse.linalg.spsolve
    openseespy_solvers.scipy.spsolve

    Examples
    --------
    >>> from openseespy_solvers.cupy import spsolve
    >>> solver = spsolve()
    >>> solver.backend
    'cupy'
    """
    return _SpSolve(
        permc_spec=permc_spec,
        scheme=scheme or "CSR",
        writable=writable,
        debug=debug,
        dtype=dtype,
    )

cg

cg(x0: Any = None, *, rtol: float = 1e-05, atol: float = 0.0, maxiter: int | None = None, M: Any = None, callback: Any = None, scheme: str | None = None, writable: str | list[str] = 'none', debug: bool = False, dtype: Any = float64) -> _CG
Source code in src/openseespy_solvers/cupy/__init__.py
def cg(
    x0: Any = None,
    *,
    rtol: float = 1e-5,
    atol: float = 0.0,
    maxiter: int | None = None,
    M: Any = None,
    callback: Any = None,
    scheme: str | None = None,
    writable: str | list[str] = "none",
    debug: bool = False,
    dtype: Any = np.float64,
) -> _CG:
    r"""Configure a Conjugate Gradient solver for OpenSees ``PythonSparse`` (GPU).

    Uses :func:`cupyx.scipy.sparse.linalg.cg` internally to solve ``Ax = b``.

    Parameters
    ----------
    x0 : cupy.ndarray, optional
        Starting guess for the solution.
    rtol, atol : float, optional
        Relative and absolute tolerances for convergence. Defaults are
        ``rtol=1e-5`` and ``atol=0.0``.
    maxiter : int, optional
        Maximum number of iterations.
    M : {sparse matrix, LinearOperator, callable}, optional
        Preconditioner for ``A``.
    callback : callable, optional
        User-supplied function called after each iteration.
    """ + _OPENSEES_LINEAR + _LINEAR_RETURNS + _LINEAR_NOTES + """
    See Also
    --------
    cupyx.scipy.sparse.linalg.cg
    scipy.sparse.linalg.cg
    openseespy_solvers.scipy.cg

    Examples
    --------
    >>> from openseespy_solvers.cupy import cg
    >>> solver = cg(rtol=1e-8)
    >>> solver._on_device
    True
    """
    return _CG(
        x0,
        rtol=rtol,
        atol=atol,
        maxiter=maxiter,
        M=M,
        callback=callback,
        scheme=scheme or "CSR",
        writable=writable,
        debug=debug,
        dtype=dtype,
    )

gmres

gmres(x0: Any = None, *, rtol: float = 1e-05, atol: float = 0.0, restart: int = 20, maxiter: int | None = None, M: Any = None, callback: Any = None, scheme: str | None = None, writable: str | list[str] = 'none', debug: bool = False, dtype: Any = float64) -> _GMRES
Source code in src/openseespy_solvers/cupy/__init__.py
def gmres(
    x0: Any = None,
    *,
    rtol: float = 1e-5,
    atol: float = 0.0,
    restart: int = 20,
    maxiter: int | None = None,
    M: Any = None,
    callback: Any = None,
    scheme: str | None = None,
    writable: str | list[str] = "none",
    debug: bool = False,
    dtype: Any = np.float64,
) -> _GMRES:
    r"""Configure a GMRES solver for OpenSees ``PythonSparse`` (GPU).

    Uses :func:`cupyx.scipy.sparse.linalg.gmres` internally.

    Parameters
    ----------
    x0 : cupy.ndarray, optional
        Starting guess for the solution.
    rtol, atol : float, optional
        Relative and absolute tolerances. Defaults are ``rtol=1e-5`` and
        ``atol=0.0``.
    restart : int, optional
        Number of iterations between restarts. Default is ``20``.
    maxiter : int, optional
        Maximum number of iterations.
    M : {sparse matrix, LinearOperator, callable}, optional
        Preconditioner for ``A``.
    callback : callable, optional
        User-supplied function called after each iteration.
    """ + _OPENSEES_LINEAR + _LINEAR_RETURNS + _LINEAR_NOTES + """
    See Also
    --------
    cupyx.scipy.sparse.linalg.gmres
    openseespy_solvers.scipy.gmres
    """
    return _GMRES(
        x0,
        rtol=rtol,
        atol=atol,
        restart=restart,
        maxiter=maxiter,
        M=M,
        callback=callback,
        scheme=scheme or "CSR",
        writable=writable,
        debug=debug,
        dtype=dtype,
    )

eigsh

eigsh(*, sigma: float | None = None, which: str = 'LM', mass_mode: str = 'general', linear_solver: LinearSolver | None = None, mode: str = 'normal', v0: Any = None, ncv: int | None = None, maxiter: int | None = None, tol: float = 0.0, scheme: str | None = None, debug: bool = False, dtype: Any = float64) -> _Eigsh
Source code in src/openseespy_solvers/cupy/__init__.py
def eigsh(
    *,
    sigma: float | None = None,
    which: str = "LM",
    mass_mode: str = "general",
    linear_solver: LinearSolver | None = None,
    mode: str = "normal",
    v0: Any = None,
    ncv: int | None = None,
    maxiter: int | None = None,
    tol: float = 0.0,
    scheme: str | None = None,
    debug: bool = False,
    dtype: Any = np.float64,
) -> _Eigsh:
    r"""Configure a generalized eigen solver for OpenSees ``PythonSparse`` (GPU).

    Matches OpenSees ``eigsh``: ``sigma=0`` and shift-invert when
    ``find_smallest=True``; no shift when ``find_smallest=False`` unless you set
    ``sigma``. ARPACK always uses ``which='LM'``.

    - ``mass_mode='general'`` (default): `scipy.sparse.linalg.eigsh`; GPU inner solves only when
      shift-invert is active.
    - ``mass_mode='diagonal'`` / ``'lumped'``: GPU shift-invert (smallest modes).

    Parameters
    ----------
    sigma : float, optional
        Shift for shift-invert. Default ``0.0`` when ``find_smallest=True``.
    which : {'LM', 'SM', 'LR', 'SR', 'LI', 'SI', 'LA', 'SA', 'BE'}, optional
        Operator eigenvalues for ARPACK. Default ``'LM'`` (shift-invert).
    mass_mode : {'diagonal', 'lumped', 'general'}, optional
        Eigen mode strategy. Default is ``'general'``.
    linear_solver : LinearSolver, optional
        Inner solver for shift-invert (``general`` and diagonal/lumped paths).
    mode : {'normal', 'buckling', 'cayley'}, optional
        `scipy.sparse.linalg.eigsh` mode (``mass_mode='general'`` only). Default ``'normal'``.
    v0 : cupy.ndarray, optional
        Starting vector for ARPACK.
    ncv : int, optional
        Number of Lanczos vectors.
    maxiter : int, optional
        Maximum number of iterations.
    tol : float, optional
        Convergence tolerance. Default is ``0.0``.
    """ + _OPENSEES_EIGEN + _EIGEN_RETURNS + _EIGEN_NOTES + """
    See Also
    --------
    lobpcg
    openseespy_solvers.scipy.eigsh
    openseespy_solvers.nvmath.direct_solver

    Notes
    -----
    `cupyx.scipy.sparse.linalg.eigsh` does not expose generalized eigen solves, so
    ``diagonal`` / ``lumped`` run cupy Lanczos on
    ``OP = (K - sigma diag(m))^{-1} diag(m)`` while ``general``
    uses `scipy.sparse.linalg.eigsh` with the same shift-invert idea and a full mass matrix.
    """
    return _Eigsh(
        sigma=sigma,
        which=which,
        mass_mode=mass_mode,
        linear_solver=linear_solver,
        mode=mode,
        v0=v0,
        ncv=ncv,
        maxiter=maxiter,
        tol=tol,
        scheme=scheme or "CSR",
        debug=debug,
        dtype=dtype,
    )

lobpcg

lobpcg(*, X: Any = None, M: Any = None, tol: float | None = None, maxiter: int = 20, largest: bool = False, rng: Any = None, scheme: str | None = None, debug: bool = False, dtype: Any = float64) -> _Lobpcg
Source code in src/openseespy_solvers/cupy/__init__.py
def lobpcg(
    *,
    X: Any = None,
    M: Any = None,
    tol: float | None = None,
    maxiter: int = 20,
    largest: bool = False,
    rng: Any = None,
    scheme: str | None = None,
    debug: bool = False,
    dtype: Any = np.float64,
) -> _Lobpcg:
    r"""Configure a LOBPCG eigenvalue solver for OpenSees ``PythonSparse`` (GPU).

    Uses :func:`cupyx.scipy.sparse.linalg.lobpcg` internally to solve
    ``K x = \lambda M x``.

    Parameters
    ----------
    X : cupy.ndarray, optional
        Initial approximation to the eigenvectors.
    M : {sparse matrix, LinearOperator, callable}, optional
        Preconditioner for ``K`` (approximate ``K^{-1}``). A callable is invoked
        as ``M(K)``; see :mod:`openseespy_solvers.cupy.precond`.
    tol : float, optional
        Convergence tolerance.
    maxiter : int, optional
        Maximum number of iterations. Default is ``20``.
    largest : bool, optional
        If ``True``, compute largest eigenvalues. Default is ``False``.
    rng : {None, int, `numpy.random.Generator`}, optional
        Random number generator used when ``X`` is ``None``.
    """ + _OPENSEES_EIGEN + _EIGEN_RETURNS + _EIGEN_NOTES + """
    See Also
    --------
    cupyx.scipy.sparse.linalg.lobpcg
    openseespy_solvers.scipy.lobpcg

    Notes
    -----
    ``M`` in LOBPCG approximates ``K^{-1}`` (not the mass matrix ``B``). Use the same
    callable preconditioner pattern as :func:`cg` / :func:`gmres``, e.g.
    ``M=precond.jacobi`` (diagonal) or ``M=precond.direct(direct_solver())`` (sparse
    direct ``K^{-1}``).
    Prefer :func:`eigsh` for shift-invert modal analysis when that fits the mass model.
    """
    return _Lobpcg(
        X=X,
        M=M,
        tol=tol,
        maxiter=maxiter,
        largest=largest,
        rng=rng,
        scheme=scheme or "CSR",
        debug=debug,
        dtype=dtype,
    )