FastLapackInterface.BunchKaufmanWsType
BunchKaufmanWs

Workspace for LinearAlgebra.BunchKaufman factorization using the LAPACK.sytrf! or LAPACK.sytrf_rook! functions for symmetric matrices, and LAPACK.hetrf! or LAPACK.hetrf_rook! functions for hermitian matrices (e.g. with ComplexF64 or ComplexF32 elements).

Examples

julia> A = [1.2 7.8
            7.8 3.3]
2×2 Matrix{Float64}:
 1.2  7.8
 7.8  3.3

julia> ws = BunchKaufmanWs(A)
BunchKaufmanWs{Float64}
  work: 128-element Vector{Float64}
  ipiv: 2-element Vector{Int64}


julia> A, ipiv, info = LAPACK.sytrf!(ws, 'U', A)
([1.2 7.8; 7.8 3.3], [-1, -1], 0)

julia> t = LinearAlgebra.BunchKaufman(A, ipiv,'U', true, false, info)
BunchKaufman{Float64, Matrix{Float64}, Vector{Int64}}
D factor:
2×2 Tridiagonal{Float64, Vector{Float64}}:
 1.2  7.8
 7.8  3.3
U factor:
2×2 UnitUpperTriangular{Float64, Matrix{Float64}}:
 1.0  0.0
  ⋅   1.0
permutation:
2-element Vector{Int64}:
 1
 2

julia> A = [1.2 7.8
            7.8 3.3]
2×2 Matrix{Float64}:
 1.2  7.8
 7.8  3.3

julia> ws = BunchKaufmanWs(A)
BunchKaufmanWs{Float64}
  work: 128-element Vector{Float64}
  ipiv: 2-element Vector{Int64}


julia> A, ipiv, info = LAPACK.sytrf_rook!(ws, 'U', A)
([1.2 7.8; 7.8 3.3], [-1, -2], 0)

julia> t = LinearAlgebra.BunchKaufman(A, ipiv,'U', true, true, info)
BunchKaufman{Float64, Matrix{Float64}, Vector{Int64}}
D factor:
2×2 Tridiagonal{Float64, Vector{Float64}}:
 1.2  7.8
 7.8  3.3
U factor:
2×2 UnitUpperTriangular{Float64, Matrix{Float64}}:
 1.0  0.0
  ⋅   1.0
permutation:
2-element Vector{Int64}:
 1
 2
FastLapackInterface.CholeskyPivotedWsType
CholeskyPivotedWs

Workspace for LinearAlgebra.CholeskyPivoted factorization using the LAPACK.pstrf! function. The standard LinearAlgebra.Cholesky uses LAPACK.potrf! which is non-allocating and does not require a separate Workspace.

Examples

julia> A = [1.2 7.8
            7.8 3.3]
2×2 Matrix{Float64}:
 1.2  7.8
 7.8  3.3

julia> ws = CholeskyPivotedWs(A)
CholeskyPivotedWs{Float64}
  work: 4-element Vector{Float64}
  piv: 2-element Vector{Int64}


julia> AA, piv, rank, info = LAPACK.pstrf!(ws, 'U', A, 1e-6)
([1.816590212458495 4.293758683992806; 7.8 -17.236363636363635], [2, 1], 1, 1)

julia> CholeskyPivoted(AA, 'U', piv, rank, 1e-6, info)
CholeskyPivoted{Float64, Matrix{Float64}, Vector{Int64}}
U factor with rank 1:
2×2 UpperTriangular{Float64, Matrix{Float64}}:
 1.81659    4.29376
  ⋅       -17.2364
permutation:
2-element Vector{Int64}:
 2
 1
FastLapackInterface.EigenWsType
EigenWs

Workspace for LinearAlgebra.Eigen factorization using the LAPACK.geevx! function.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = EigenWs(A, rvecs=true)
EigenWs{Float64, Matrix{Float64}, Float64}
  work: 260-element Vector{Float64}
  rwork: 2-element Vector{Float64}
  VL: 0×2 Matrix{Float64}
  VR: 2×2 Matrix{Float64}
  W: 2-element Vector{Float64}
  scale: 2-element Vector{Float64}
  iwork: 0-element Vector{Int64}
  rconde: 0-element Vector{Float64}
  rcondv: 0-element Vector{Float64}


julia> t = LAPACK.geevx!(ws, 'N', 'N', 'V', 'N', A);

julia> LinearAlgebra.Eigen(t[2], t[5])
Eigen{Float64, Float64, Matrix{Float64}, Vector{Float64}}
values:
2-element Vector{Float64}:
 -1.6695025194532018
  6.169502519453203
vectors:
2×2 Matrix{Float64}:
 -0.625424  -0.420019
  0.780285  -0.907515
FastLapackInterface.GeneralizedEigenWsType
GeneralizedEigenWs

Workspace that can be used for LinearAlgebra.GeneralizedEigen factorization using LAPACK.ggev!. Supports Real and Complex matrices.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> B = [8.2 1.7
            5.9 2.1]
2×2 Matrix{Float64}:
 8.2  1.7
 5.9  2.1

julia> ws = GeneralizedEigenWs(A, rvecs=true)
GeneralizedEigenWs{Float64, Matrix{Float64}, Float64}
  work: 78-element Vector{Float64}
  vl: 0×2 Matrix{Float64}
  vr: 2×2 Matrix{Float64}
  αr: 2-element Vector{Float64}
  αi: 2-element Vector{Float64}
  β: 2-element Vector{Float64}


julia> αr, αi, β, _, vr = LAPACK.ggev!(ws, 'N', 'V', A, B);

julia> LinearAlgebra.GeneralizedEigen(αr ./ β, vr)
GeneralizedEigen{Float64, Float64, Matrix{Float64}, Vector{Float64}}
values:
2-element Vector{Float64}:
 -0.8754932558185097
  1.6362721153456299
vectors:
2×2 Matrix{Float64}:
 -0.452121  -0.0394242
  1.0        1.0
FastLapackInterface.GeneralizedSchurWsType
GeneralizedSchurWs

Workspace to be used with the LinearAlgebra.GeneralizedSchur representation of the Generalized Schur decomposition which uses the LAPACK.gges! function.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> B = [8.2 0.3
            1.7 4.3]
2×2 Matrix{Float64}:
 8.2  0.3
 1.7  4.3

julia> ws = GeneralizedSchurWs(A)
GeneralizedSchurWs{Float64}
  work: 90-element Vector{Float64}
  αr: 2-element Vector{Float64}
  αi: 2-element Vector{Float64}
  β: 2-element Vector{Float64}
  vsl: 2×2 Matrix{Float64}
  vsr: 2×2 Matrix{Float64}
  sdim: Base.RefValue{Int64}
  bwork: 2-element Vector{Int64}
  eigen_values: 2-element Vector{ComplexF64}
  
julia> t = GeneralizedSchur(LAPACK.gges!(ws, 'V','V', A, B)...)
GeneralizedSchur{Float64, Matrix{Float64}, Vector{ComplexF64}, Vector{Float64}}
S factor:
2×2 Matrix{Float64}:
 -1.43796  1.63843
  0.0      7.16295
T factor:
2×2 Matrix{Float64}:
 5.06887  -4.00221
 0.0       6.85558
Q factor:
2×2 Matrix{Float64}:
 -0.857329  0.514769
  0.514769  0.857329
Z factor:
2×2 Matrix{Float64}:
 -0.560266  0.828313
  0.828313  0.560266
α:
2-element Vector{ComplexF64}:
 -1.4379554610733563 + 0.0im
   7.162947865097022 + 0.0im
β:
2-element Vector{Float64}:
 5.068865029631368
 6.855578082442485
FastLapackInterface.HermitianEigenWsType
HermitianEigenWs

Workspace to be used with Hermitian diagonalization using the LAPACK.syevr! function. Supports both Real and Complex Hermitian matrices.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = HermitianEigenWs(A, vecs=true)
HermitianEigenWs{Float64, Matrix{Float64}, Float64}
  work: 66-element Vector{Float64}
  rwork: 0-element Vector{Float64}
  iwork: 20-element Vector{Int64}
  w: 2-element Vector{Float64}
  Z: 2×2 Matrix{Float64}
  isuppz: 4-element Vector{Int64}


julia> LinearAlgebra.Eigen(LAPACK.syevr!(ws, 'V', 'A', 'U', A, 0.0, 0.0, 0, 0, 1e-6)...)
Eigen{Float64, Float64, Matrix{Float64}, Vector{Float64}}
values:
2-element Vector{Float64}:
 -0.2783393759541063
  4.778339375954106
vectors:
2×2 Matrix{Float64}:
 -0.841217  0.540698
  0.540698  0.841217
FastLapackInterface.LSEWsType
LSEWs

Workspace for the least squares solving function LAPACK.geqrf!.

Examples

julia> A = [1.2 2.3 6.2
            6.2 3.3 8.8
            9.1 2.1 5.5]
3×3 Matrix{Float64}:
 1.2  2.3  6.2
 6.2  3.3  8.8
 9.1  2.1  5.5

julia> B = [2.7 3.1 7.7
            4.1 8.1 1.8]
2×3 Matrix{Float64}:
 2.7  3.1  7.7
 4.1  8.1  1.8

julia> c = [0.2, 7.2, 2.9]
3-element Vector{Float64}:
 0.2
 7.2
 2.9

julia> d = [3.9, 2.1]
2-element Vector{Float64}:
 3.9
 2.1

julia> ws = LSEWs(A, B)
LSEWs{Float64}
  work: 101-element Vector{Float64}
  X: 3-element Vector{Float64}

julia> LAPACK.gglse!(ws, A, c, B, d)
([0.19723156207005318, 0.0683561362406917, 0.40981438442398854], 13.750943845251626)
FastLapackInterface.LUWsType
LUWs

Workspace to be used with the LinearAlgebra.LU representation of the LU factorization which uses the LAPACK.getrf! function.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = LUWs(A)
LUWs
  ipiv: 2-element Vector{Int64}

julia> t = LU(LAPACK.getrf!(ws, A)...)
LU{Float64, Matrix{Float64}, Vector{Int64}}
L factor:
2×2 Matrix{Float64}:
 1.0       0.0
 0.193548  1.0
U factor:
2×2 Matrix{Float64}:
 6.2  3.3
 0.0  1.66129
FastLapackInterface.QROrgWsType
QROrgWs

Workspace to be used with the LinearAlgebra.LAPACK.orgqr! function. It requires the workspace of a QR or a QRPivoted previous factorization

Examples

```jldoctest julia> A = [1.2 2.3 6.2 3.3] 2×2 Matrix{Float64}: 1.2 2.3 6.2 3.3

julia> ws = QRPivotedWs(A) QRPivotedWs{Float64, Float64} work: 100-element Vector{Float64} rwork: 0-element Vector{Float64} τ: 2-element Vector{Float64} jpvt: 2-element Vector{Int64}

julia> C=[1 0.5; 2 1] 2×2 Matrix{Float64}: 1.0 0.5 2.0 1.0

julia> orgws = QROrgWs(ws, 'L', 'N', A, C) QROrgWs{Float64} work: 4224-element Vector{Float64} τ: 2-element Vector{Float64}

FastLapackInterface.QROrmWsType
QROrmWs

Workspace to be used with the LinearAlgebra.LAPACK.ormqr! function. It requires the workspace of a QR or a QRPivoted previous factorization

Examples

```jldoctest julia> A = [1.2 2.3 6.2 3.3] 2×2 Matrix{Float64}: 1.2 2.3 6.2 3.3

julia> ws = QRPivotedWs(A) QRPivotedWs{Float64, Float64} work: 100-element Vector{Float64} rwork: 0-element Vector{Float64} τ: 2-element Vector{Float64} jpvt: 2-element Vector{Int64}

julia> C=[1 0.5; 2 1] 2×2 Matrix{Float64}: 1.0 0.5 2.0 1.0

julia> ormws = QROrmWs(ws, 'L', 'N', A, C) QROrmWs{Float64} work: 4224-element Vector{Float64} τ: 2-element Vector{Float64}

FastLapackInterface.QRPivotedWsType
QRPivotedWs

Workspace to be used with the LinearAlgebra.QRPivoted representation of the QR factorization which uses the LAPACK.geqp3! function.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = QRPivotedWs(A)
QRPivotedWs{Float64, Float64}
  work: 100-element Vector{Float64}
  rwork: 0-element Vector{Float64}
  τ: 2-element Vector{Float64}
  jpvt: 2-element Vector{Int64}

julia> t = QRPivoted(LAPACK.geqp3!(ws, A)...)
QRPivoted{Float64, Matrix{Float64}, Vector{Float64}, Vector{Int64}}
Q factor: 2×2 LinearAlgebra.QRPackedQ{Float64, Matrix{Float64}, Vector{Float64}}
R factor:
2×2 Matrix{Float64}:
 -6.31506  -3.67692
  0.0      -1.63102
permutation:
2-element Vector{Int64}:
 1
 2

julia> Matrix(t)
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3
FastLapackInterface.QRWYWsType
QRWYWs

Workspace to be used with the LinearAlgebra.QRCompactWY representation of the blocked QR factorization which uses the LAPACK.geqrt! function. By default the blocksize for the algorithm is taken as min(36, min(size(template))), this can be overridden by using the blocksize keyword of the constructor.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = QRWYWs(A)
QRWYWs{Float64, Matrix{Float64}}
  work: 4-element Vector{Float64}
  T: 2×2 Matrix{Float64}

julia> t = LinearAlgebra.QRCompactWY(LAPACK.geqrt!(ws, A)...)
LinearAlgebra.QRCompactWY{Float64, Matrix{Float64}, Matrix{Float64}}
Q factor: 2×2 LinearAlgebra.QRCompactWYQ{Float64, Matrix{Float64}, Matrix{Float64}}
R factor:
2×2 Matrix{Float64}:
 -6.31506  -3.67692
  0.0      -1.63102

julia> Matrix(t)
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3
FastLapackInterface.QRWsType
QRWs

Workspace for standard LinearAlgebra.QR factorization using the LAPACK.geqrf! function.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = QRWs(A)
QRWs{Float64}
  work: 64-element Vector{Float64}
  τ: 2-element Vector{Float64}

julia> t = QR(LAPACK.geqrf!(ws, A)...)
QR{Float64, Matrix{Float64}, Vector{Float64}}
Q factor: 2×2 LinearAlgebra.QRPackedQ{Float64, Matrix{Float64}, Vector{Float64}}
R factor:
2×2 Matrix{Float64}:
 -6.31506  -3.67692
  0.0      -1.63102

julia> Matrix(t)
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3
FastLapackInterface.SCHURORDERType

enum SCHURORDER lp rp id ed

Keywords

- `lp`: Left plane (real(eigenvalue) < criterium)
- `rp`: Right plane (real(eigenvalue) >= criterium)
- `id`: Interior of disk (abs(eigenvalue)^2 < criterium) 
- `ed`: Exterior of disk (abs(eigenvalue)^2 >= criterium)

Note

- the left half-plane is obtained with criterium = 0
- the unit disk is obtained with criterium = 1
- because of numerical error in computing repeated eigenvalues, you need to adapt
  criterium depending whether you want to include or not 0 is the left half-plane or
  1 in the unit disk
- criterium is passed as optional parameter to `gees` and `gges` functions
FastLapackInterface.SchurWsType
SchurWs

Workspace to be used with the LinearAlgebra.Schur representation of the Schur decomposition which uses the LAPACK.gees! function.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = SchurWs(A)
SchurWs{Float64}
  work: 68-element Vector{Float64}
  wr: 2-element Vector{Float64}
  wi: 2-element Vector{Float64}
  vs: 2×2 Matrix{Float64}
  sdim: Base.RefValue{Int64}
  bwork: 2-element Vector{Int64}
  eigen_values: 2-element Vector{ComplexF64}

julia> t = Schur(LAPACK.gees!(ws, 'V', A)...)
Schur{Float64, Matrix{Float64}, Vector{Float64}}
T factor:
2×2 Matrix{Float64}:
 -1.6695  -3.9
  0.0      6.1695
Z factor:
2×2 Matrix{Float64}:
 -0.625424  -0.780285
  0.780285  -0.625424
eigenvalues:
2-element Vector{Float64}:
 -1.6695025194532018
  6.169502519453203

julia> Matrix(t)
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3
FastLapackInterface.WorkspaceMethod
Workspace(lapack_function, A)

Will create the correct Workspace for the target lapack_function and matrix A.

Examples

julia> A = [1.2 2.3
            6.2 3.3]
2×2 Matrix{Float64}:
 1.2  2.3
 6.2  3.3

julia> ws = Workspace(LAPACK.geqrt!, A)
QRWYWs{Float64, Matrix{Float64}}
  work: 4-element Vector{Float64}
  T: 2×2 Matrix{Float64}


julia> LinearAlgebra.QRCompactWY(factorize!(ws, A)...)
LinearAlgebra.QRCompactWY{Float64, Matrix{Float64}, Matrix{Float64}}
Q factor: 2×2 LinearAlgebra.QRCompactWYQ{Float64, Matrix{Float64}, Matrix{Float64}}
R factor:
2×2 Matrix{Float64}:
 -6.31506  -3.67692
  0.0      -1.63102
Base.resize!Method
resize!(ws, A; kwargs...)

Resizes the ws to be appropriate for use with matrix A. The kwargs can be used to communicate which features should be supported by the Workspace, such as left and right eigenvectors while using EigenWs. This function is mainly used for automatic resizing inside LAPACK functions.

LinearAlgebra.LAPACK.gees!Method
gees!(ws, jobvs, A; select=nothing, criterium=0.0, 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 select the eigenvalues appearing in the top left corner of the Schur form:

  • by setting the select option to one of
    • lp: Left plane (real(eigenvalue) < criterium)
    • rp: Right plane (real(eigenvalue) >= criterium)
    • id: Interior of disk (abs(eigenvalue)^2 < criterium)
    • ed: Exterior of disk (abs(eigenvalue)^2 >= criterium)
    and setting criterium.
  • by setting select equal to a function used to sort the eigenvalues during the decomponsition. In this case, the criterium keyword isn't used.

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.

See also FastLapackInterface.SCHURODER

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.

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.

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.

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.

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).

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).

LinearAlgebra.LAPACK.gges!Method
gges!(ws, jobvsl, jobvsr, A, B; select=nothing, criterium = 0, 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 select the eigenvalues appearing in the top left corner of the Schur form:

  • by setting the select option to one of
    • lp: Left plane (real(eigenvalue) < criterium)
    • rp: Right plane (real(eigenvalue) >= criterium)
    • id: Interior of disk (abs(eigenvalue)^2 < criterium)
    • ed: Exterior of disk (abs(eigenvalue)^2 >= criterium)
    and setting criterium.
  • by setting select equal to a function used to sort the eigenvalues during the decomposition. In this case, the criterium keyword isn't used.

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.

See also FastLapackInterface.SCHURODER

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.

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.

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.

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.

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.

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.