LAPACK

This section details the LAPACK functions that are supported for use with various Workspaces. Each function has a resize keyword argument that is true by default, allowing for automatic resizing of the workspaces to accomodate larger or smaller Matrices or different features than they were originally constructed for. This is provided as a convenience but involves an efficiency cost. When working with matrices of different sizes, the best strategy is to successively apply a function to all matrices of the same size and to minimize triggering the resizing mechanism.

Unified Interface

After having created the Workspace that corresponds to the targeted factorization or decomposition, one of the following two aliases can be used to dispatch the call to the correct LAPACK function.

QR

LinearAlgebra.LAPACK.geqrf!Method
geqrf!(ws, A; resize=true) -> (A, ws.τ)

Compute the QR factorization of A, A = QR, using previously allocated QRWs workspace ws. ws.τ contains scalars which parameterize the elementary reflectors of the factorization. ws.τ must have length greater than or equal to the smallest dimension of A. If this is not the case, and resize==true the workspace will be automatically resized to the appropriate size.

A and ws.τ modified in-place.

source
LinearAlgebra.LAPACK.geqrt!Method
geqrt!(ws, A; resize=true) -> (A, ws.T)

Compute the blocked QR factorization of A, A = QR, using a preallocated QRWYWs workspace ws. ws.T contains upper triangular block reflectors which parameterize the elementary reflectors of the factorization. The first dimension of ws.T sets the block size and it must satisfy 1 <= size(ws.T, 1) <= min(size(A)...). The second dimension of T must equal the smallest dimension of A, i.e. size(ws.T, 2) == size(A, 2). If this is not the case and resize==true, the workspace will automatically be resized to the appropriate dimensions.

A and ws.T are modified in-place.

source
LinearAlgebra.LAPACK.geqp3!Method
geqp3!(ws, A; resize=true) -> (A, ws.τ, ws.jpvt)

Compute the pivoted QR factorization of A, AP = QR using BLAS level 3, using the preallocated QRPivotedWs workspace ws. P is a pivoting matrix, represented by ws.jpvt. ws.τ stores the elementary reflectors. ws.jpvt must have length greater than or equal to n if A is an (m x n) matrix and ws.τ must have length greater than or equal to the smallest dimension of A. If this is not the case and resize == true the workspace will be appropriately resized.

A, ws.jpvt, and ws.τ are modified in-place.

source
LinearAlgebra.LAPACK.ormqr!Method
ormqr!(ws, side, trans, A, C) -> C

Computes Q * C (trans = N), transpose(Q) * C (trans = T), adjoint(Q) * C (trans = C) for side = L or the equivalent right-sided multiplication for side = R using Q from a QR factorization of A. Uses preallocated workspace ws::QROrmWs and the factors are assumed to be stored in ws.τ. C is overwritten.

source

Schur

LinearAlgebra.LAPACK.gees!Method
gees!(ws, jobvs, A; select=nothing, resize=true) -> (A, vs, ws.eigen_values)

Computes the eigenvalues (jobvs = N) or the eigenvalues and Schur vectors (jobvs = V) of matrix A, using the preallocated SchurWs worspace ws. If ws is not of the appropriate size and resize==true it will be resized for A. A is overwritten by its Schur form, and ws.eigen_values is overwritten with the eigenvalues.

It is possible to specify select, a function used to sort the eigenvalues during the decomponsition. The function should have the signature f(wr::T, wi::T) -> Bool, where wr and wi are the real and imaginary parts of the eigenvalue, and T == eltype(A).

Returns A, vs containing the Schur vectors, and ws.eigen_values.

source
LinearAlgebra.LAPACK.gges!Method
gges!(ws, jobvsl, jobvsr, A, B; select=nothing, resize=true) -> (A, B, ws.α, ws.β, ws.vsl, ws.vsr)

Computes the generalized eigenvalues, generalized Schur form, left Schur vectors (jobsvl = V), or right Schur vectors (jobvsr = V) of A and B, using preallocated GeneralizedSchurWs workspace ws. If ws is not of the right size, and resize==true it will be resized appropriately.

It is possible to specify select, a function used to sort the eigenvalues to the top left corner of the Schur form. The function should have the signature f(αr::T, αi::T, β::T) -> Bool where T == eltype(A). An eigenvalue (αr[j]+αi[j])/β[j] is selected if f(αr[j],αi[j],β[j]) is true, i.e. if either one of a complex conjugate pair of eigenvalues is selected, then both complex eigenvalues are selected. The generalized eigenvalues components are returned in ws.α and ws.β where ws.α is a complex vector and ẁs.β, a real vector. The generalized eigenvalues (ws.α./ws.β) are returned in ws.eigen_values, a complex vector. The left Schur vectors are returned in ws.vsl and the right Schur vectors are returned in ws.vsr.

source

LU

LinearAlgebra.LAPACK.getrf!Method
getrf!(ws, A; resize=true) -> (A, ws.ipiv, info)

Compute the pivoted LU factorization of A, A = LU, using the preallocated LUWs workspace ws. If the workspace is too small and resize==true it will be resized appropriately for A.

Returns A, modified in-place, ws.ipiv, the pivoting information, and the ws.info code which indicates success (info = 0), a singular value in U (info = i, in which case U[i,i] is singular), or an error code (info < 0).

source
LinearAlgebra.LAPACK.getrs!Method
getrs!(ws, trans, A, B)

Solves the linear equation A * X = B, transpose(A) * X = B, or adjoint(A) * X = B for square A. Modifies the matrix/vector B in place with the solution. A is the LU factorization from getrf! with the pivoting information stored in ws.ipiv. trans may be one of N (no modification), T (transpose), or C (conjugate transpose).

source

Eigen

LinearAlgebra.LAPACK.geevx!Method
geevx!(ws, balanc, jobvl, jobvr, sense, A; resize=true) -> (A, ws.W, [ws.rwork,] ws.VL, ws.VR, ilo, ihi, ws.scale, abnrm, ws.rconde, ws.rcondv)

Finds the eigensystem of A with matrix balancing using a preallocated EigenWs. If jobvl = N, the left eigenvectors of A aren't computed. If jobvr = N, the right eigenvectors of A aren't computed. If jobvl = V or jobvr = V, the corresponding eigenvectors are computed. If balanc = N, no balancing is performed. If balanc = P, A is permuted but not scaled. If balanc = S, A is scaled but not permuted. If balanc = B, A is permuted and scaled. If sense = N, no reciprocal condition numbers are computed. If sense = E, reciprocal condition numbers are computed for the eigenvalues only. If sense = V, reciprocal condition numbers are computed for the right eigenvectors only. If sense = B, reciprocal condition numbers are computed for the right eigenvectors and the eigenvectors. If sense = E,B, the right and left eigenvectors must be computed. ws.rwork is only returned in the Real case. If ws does not have the appropriate size for A and the work to be done, if resize=true, it will be automatically resized accordingly.

source
LinearAlgebra.LAPACK.syevr!Method
syevr!(ws, jobz, range, uplo, A, vl, vu, il, iu, abstol; resize=true) -> (ws.W, ws.Z)

Finds the eigenvalues (jobz = N) or eigenvalues and eigenvectors (jobz = V) of a symmetric matrix A using a preallocated HermitianEigenWs. If the workspace is not appropriate for A and resize==true it will be automatically resized. If uplo = U, the upper triangle of A is used. If uplo = L, the lower triangle of A is used. If range = A, all the eigenvalues are found. If range = V, the eigenvalues in the half-open interval (vl, vu] are found. If range = I, the eigenvalues with indices between il and iu are found. abstol can be set as a tolerance for convergence.

The eigenvalues are returned as ws.W and the eigenvectors in ws.Z.

source
LinearAlgebra.LAPACK.ggev!Method
ggev!(ws, jobvl, jobvr, A, B; resize=true) -> (ws.αr, [ws.αi,], ws.β, ws.vl, ws.vr)

Finds the generalized eigendecomposition of A and B usin a preallocated GeneralizedEigenWs. If the workspace is not appropriately sized and resize == true, it will automatically be resized. If jobvl = N, the left eigenvectors aren't computed. If jobvr = N, the right eigenvectors aren't computed. If jobvl = V or jobvr = V, the corresponding eigenvectors are computed. ws.αi is only returned in the Real case.

source

BunchKaufman

LinearAlgebra.LAPACK.sytrf!Method
sytrf!(ws, uplo, A; resize=true) -> (A, ws.ipiv, info)

Computes the Bunch-Kaufman factorization of a symmetric matrix A, using previously allocated workspace ws. If the workspace was too small and resize==true it will automatically resized. If uplo = U, the upper half of A is stored. If uplo = L, the lower half is stored.

Returns A, overwritten by the factorization, a pivot vector ws.ipiv, and the error code info which is a non-negative integer. If info is positive the matrix is singular and the diagonal part of the factorization is exactly zero at position info.

source

Cholesky

LinearAlgebra.LAPACK.pstrf!Method
pstrf!(ws, uplo, A, tol; resize=true) -> (A, ws.piv, rank, info)

Computes the (upper if uplo = U, lower if uplo = L) pivoted Cholesky decomposition of positive-definite matrix A with a user-set tolerance tol, using a preallocated CholeskyPivotedWs. If the workspace was too small and resize==true it will be automatically resized. A is overwritten by its Cholesky decomposition.

Returns A, the pivots piv, the rank of A, and an info code. If info = 0, the factorization succeeded. If info = i > 0, then A is indefinite or rank-deficient.

source

LSE

LinearAlgebra.LAPACK.gglse!Method
gglse!(ws, A, c, B, d) -> (ws.X,res)

Solves the equation A * x = c where x is subject to the equality constraint B * x = d. Uses the formula ||c - A*x||^2 = 0 to solve. Uses preallocated LSEWs to store X and work buffers. Returns ws.X and the residual sum-of-squares.

source