18.3 Matrix Factorizations

 
: R = chol (A)
: [R, p] = chol (A)
: [R, p, Q] = chol (A)
: [R, p, Q] = chol (A, "vector")
: [L, …] = chol (…, "lower")
: [R, …] = chol (…, "upper")

Compute the upper Cholesky factor, R, of the real symmetric or complex Hermitian positive definite matrix A.

The upper Cholesky factor R is computed by using the upper triangular part of matrix A and is defined by

R' * R = A.

Calling chol using the optional "upper" flag has the same behavior. In contrast, using the optional "lower" flag, chol returns the lower triangular factorization, computed by using the lower triangular part of matrix A, such that

L * L' = A.

Called with one output argument chol fails if matrix A is not positive definite. Note that if matrix A is not real symmetric or complex Hermitian then the lower triangular part is considered to be the (complex conjugate) transpose of the upper triangular part, or vice versa, given the "lower" flag.

Called with two or more output arguments p flags whether the matrix A was positive definite and chol does not fail. A zero value of p indicates that matrix A is positive definite and R gives the factorization. Otherwise, p will have a positive value.

If called with three output arguments matrix A must be sparse and a sparsity preserving row/column permutation is applied to matrix A prior to the factorization. That is R is the factorization of A(Q,Q) such that

R' * R = Q' * A * Q.

The sparsity preserving permutation is generally returned as a matrix. However, given the optional flag "vector", Q will be returned as a vector such that

R' * R = A(Q, Q).

In general the lower triangular factorization is significantly faster for sparse matrices.

See also: hess, lu, qr, qz, schur, svd, ichol, cholinv, chol2inv, cholupdate, cholinsert, choldelete, cholshift.

 
: Ainv = cholinv (A)

Compute the inverse of the symmetric positive definite matrix A using the Cholesky factorization.

See also: chol, chol2inv, inv.

 
: Ainv = chol2inv (R)

Invert a symmetric, positive definite square matrix from its Cholesky decomposition, R.

Note that R should be an upper-triangular matrix with positive diagonal elements. chol2inv (U) provides inv (R'*R) but is much faster than using inv.

See also: chol, cholinv, inv.

 
: [R1, info] = cholupdate (R, u, op)

Update or downdate a Cholesky factorization.

Given an upper triangular matrix R and a column vector u, attempt to determine another upper triangular matrix R1 such that

  • R1’*R1 = R’*R + u*u’ if op is "+"
  • R1’*R1 = R’*R - u*u’ if op is "-"

If op is "-", info is set to

  • 0 if the downdate was successful,
  • 1 if R’*R - u*u’ is not positive definite,
  • 2 if R is singular.

If info is not present, an error message is printed in cases 1 and 2.

See also: chol, cholinsert, choldelete, cholshift.

 
: R1 = cholinsert (R, j, u)
: [R1, info] = cholinsert (R, j, u)

Update a Cholesky factorization given a row or column to insert in the original factored matrix.

Given a Cholesky factorization of a real symmetric or complex Hermitian positive definite matrix A = R’*R, R upper triangular, return the Cholesky factorization of A1, where A1(p,p) = A, A1(:,j) = A1(j,:)’ = u and p = [1:j-1,j+1:n+1]. u(j) should be positive.

On return, info is set to

  • 0 if the insertion was successful,
  • 1 if A1 is not positive definite,
  • 2 if R is singular.

If info is not present, an error message is printed in cases 1 and 2.

See also: chol, cholupdate, choldelete, cholshift.

 
: R1 = choldelete (R, j)

Update a Cholesky factorization given a row or column to delete from the original factored matrix.

Given a Cholesky factorization of a real symmetric or complex Hermitian positive definite matrix A = R’*R, R upper triangular, return the Cholesky factorization of A(p,p), where p = [1:j-1,j+1:n+1].

See also: chol, cholupdate, cholinsert, cholshift.

 
: R1 = cholshift (R, i, j)

Update a Cholesky factorization given a range of columns to shift in the original factored matrix.

Given a Cholesky factorization of a real symmetric or complex Hermitian positive definite matrix A = R’*R, R upper triangular, return the Cholesky factorization of A(p,p), where p is the permutation
p = [1:i-1, shift(i:j, 1), j+1:n] if i < j
or
p = [1:j-1, shift(j:i,-1), i+1:n] if j < i.

See also: chol, cholupdate, cholinsert, choldelete.

 
: H = hess (A)
: [P, H] = hess (A)

Compute the Hessenberg decomposition of the matrix A.

The Hessenberg decomposition is P * H * P' = A where P is a square unitary matrix (P' * P = I, using complex-conjugate transposition) and H is upper Hessenberg (H(i, j) = 0 forall i > j+1).

The Hessenberg decomposition is usually used as the first step in an eigenvalue computation, but has other applications as well (see Golub, Nash, and Van Loan, IEEE Transactions on Automatic Control, 1979).

See also: eig, chol, lu, qr, qz, schur, svd.

 
: [L, U] = lu (A)
: [L, U, P] = lu (A)
: [L, U, P, Q] = lu (S)
: [L, U, P, Q, R] = lu (S)
: […] = lu (S, thresh)
: y = lu (…)
: […] = lu (…, "vector")

Compute the LU decomposition of A.

If A is full then subroutines from LAPACK are used, and if A is sparse then UMFPACK is used.

The result is returned in a permuted form, according to the optional return value P. For example, given the matrix A = [1, 2; 3, 4],

[L, U, P] = lu (A)

returns

L =

  1.00000  0.00000
  0.33333  1.00000

U =

  3.00000  4.00000
  0.00000  0.66667

P =

  0  1
  1  0

The matrix is not required to be square.

When called with two or three output arguments and a sparse input matrix, lu does not attempt to perform sparsity preserving column permutations. Called with a fourth output argument, the sparsity preserving column transformation Q is returned, such that P * A * Q = L * U. This is the preferred way to call lu with sparse input matrices.

Called with a fifth output argument and a sparse input matrix, lu attempts to use a scaling factor R on the input matrix such that P * (R \ A) * Q = L * U. This typically leads to a sparser and more stable factorization.

An additional input argument thresh that defines the pivoting threshold can be given. thresh can be a scalar, in which case it defines the UMFPACK pivoting tolerance for both symmetric and unsymmetric cases. If thresh is a 2-element vector, then the first element defines the pivoting tolerance for the unsymmetric UMFPACK pivoting strategy and the second for the symmetric strategy. By default, the values defined by spparms are used ([0.1, 0.001]).

Given the string argument "vector", lu returns the values of P and Q as vector values, such that for full matrix, A(P,:) = L * U, and R(P,:) * A(:,Q) = L * U.

With two output arguments, returns the permuted forms of the upper and lower triangular matrices, such that A = L * U. With one output argument y, then the matrix returned by the LAPACK routines is returned. If the input matrix is sparse then the matrix L is embedded into U to give a return value similar to the full case. For both full and sparse matrices, lu loses the permutation information.

See also: luupdate, ilu, chol, hess, qr, qz, schur, svd.

 
: [L, U] = luupdate (L, U, x, y)
: [L, U, P] = luupdate (L, U, P, x, y)

Given an LU factorization of a real or complex matrix A = L*U, L lower unit trapezoidal and U upper trapezoidal, return the LU factorization of A + x*y.’, where x and y are column vectors (rank-1 update) or matrices with equal number of columns (rank-k update).

Optionally, row-pivoted updating can be used by supplying a row permutation (pivoting) matrix P; in that case, an updated permutation matrix is returned. Note that if L, U, P is a pivoted LU factorization as obtained by lu:

[L, U, P] = lu (A);

then a factorization of A+x*y.' can be obtained either as

[L1, U1] = lu (L, U, P*x, y)

or

[L1, U1, P1] = lu (L, U, P, x, y)

The first form uses the unpivoted algorithm, which is faster, but less stable. The second form uses a slower pivoted algorithm, which is more stable.

The matrix case is done as a sequence of rank-1 updates; thus, for large enough k, it will be both faster and more accurate to recompute the factorization from scratch.

See also: lu, cholupdate, qrupdate.

 
: [Q, R] = qr (A)
: [Q, R, P] = qr (A)
: X = qr (A) # non-sparse A
: R = qr (A) # sparse A
: X = qr (A, B) # sparse A
: [C, R] = qr (A, B)
: […] = qr (…, 0)
: […] = qr (…, "econ")
: […] = qr (…, "vector")
: […] = qr (…, "matrix")

Compute the QR factorization of A, using standard LAPACK subroutines.

The QR factorization is

Q * R = A

where Q is an orthogonal matrix and R is upper triangular.

For example, given the matrix A = [1, 2; 3, 4],

[Q, R] = qr (A)

returns

Q =

  -0.31623  -0.94868
  -0.94868   0.31623

R =

  -3.16228  -4.42719
   0.00000  -0.63246

which multiplied together return the original matrix

Q * R
  ⇒
     1.0000   2.0000
     3.0000   4.0000

If just a single return value is requested then it is either R, if A is sparse, or X, such that R = triu (X) if A is full. (Note: unlike most commands, the single return value is not the first return value when multiple values are requested.)

If a third output P is requested, then qr calculates the permuted QR factorization

Q * R = A * P

where Q is an orthogonal matrix, R is upper triangular, and P is a permutation matrix.

If A is dense, the permuted QR factorization has the additional property that the diagonal entries of R are ordered by decreasing magnitude. In other words, abs (diag (R)) will be ordered from largest to smallest.

If A is sparse, P is a fill-reducing ordering of the columns of A. In that case, the diagonal entries of R are not ordered by decreasing magnitude.

For example, given the matrix A = [1, 2; 3, 4],

[Q, R, P] = qr (A)

returns

Q =

  -0.44721  -0.89443
  -0.89443   0.44721

R =

  -4.47214  -3.13050
   0.00000   0.44721

P =

   0  1
   1  0

If the input matrix A is sparse, the sparse QR factorization is computed by using SPQR or CXSPARSE (e.g., if SPQR is not available). Because the matrix Q is, in general, a full matrix, it is recommended to request only one return value R. In that case, the computation avoids the construction of Q and returns a sparse R such that R = chol (A' * A).

If A is dense, an additional matrix B is supplied and two return values are requested, then qr returns C, where C = Q' * B. This allows the least squares approximation of A \ B to be calculated as

[C, R] = qr (A, B)
X = R \ C

If A is a sparse MxN matrix and an additional matrix B is supplied, one or two return values are possible. If one return value X is requested and M < N, then X is the minimum 2-norm solution of A \ B. If M >= N, X is the least squares approximation of A \ B. If two return values are requested, C and R have the same meaning as in the dense case (C is dense and R is sparse). The version with one return parameter should be preferred because it uses less memory and can handle rank-deficient matrices better.

If the final argument is the string "vector" then P is a permutation vector (of the columns of A) instead of a permutation matrix. In this case, the defining relationship is:

Q * R = A(:, P)

The default, however, is to return a permutation matrix and this may be explicitly specified by using a final argument of "matrix".

If the final argument is the scalar 0 or the string "econ", an economy factorization is returned. If the original matrix A has size MxN and M > N, then the economy factorization will calculate just N rows in R and N columns in Q and omit the zeros in R. If M ≤ N, there is no difference between the economy and standard factorizations. When calculating an economy factorization and A is dense, the output P is always a vector rather than a matrix. If A is sparse, output P is a sparse permutation matrix.

Background: The QR factorization has applications in the solution of least squares problems

min norm (A*x - b)

for overdetermined systems of equations (i.e., A is a tall, thin matrix).

The permuted QR factorization [Q, R, P] = qr (A) allows the construction of an orthogonal basis of span (A).

See also: chol, hess, lu, qz, schur, svd, qrupdate, qrinsert, qrdelete, qrshift.

 
: [Q1, R1] = qrupdate (Q, R, u, v)

Update a QR factorization given update vectors or matrices.

Given a QR factorization of a real or complex matrix A = Q*R, Q unitary and R upper trapezoidal, return the QR factorization of A + u*v, where u and v are column vectors (rank-1 update) or matrices with equal number of columns (rank-k update). Notice that the latter case is done as a sequence of rank-1 updates; thus, for k large enough, it will be both faster and more accurate to recompute the factorization from scratch.

The QR factorization supplied may be either full (Q is square) or economized (R is square).

See also: qr, qrinsert, qrdelete, qrshift.

 
: [Q1, R1] = qrinsert (Q, R, j, x, orient)

Update a QR factorization given a row or column to insert in the original factored matrix.

Given a QR factorization of a real or complex matrix A = Q*R, Q unitary and R upper trapezoidal, return the QR factorization of [A(:,1:j-1) x A(:,j:n)], where u is a column vector to be inserted into A (if orient is "col"), or the QR factorization of [A(1:j-1,:);x;A(:,j:n)], where x is a row vector to be inserted into A (if orient is "row").

The default value of orient is "col". If orient is "col", u may be a matrix and j an index vector resulting in the QR factorization of a matrix B such that B(:,j) gives u and B(:,j) = [] gives A. Notice that the latter case is done as a sequence of k insertions; thus, for k large enough, it will be both faster and more accurate to recompute the factorization from scratch.

If orient is "col", the QR factorization supplied may be either full (Q is square) or economized (R is square).

If orient is "row", full factorization is needed.

See also: qr, qrupdate, qrdelete, qrshift.

 
: [Q1, R1] = qrdelete (Q, R, j, orient)

Update a QR factorization given a row or column to delete from the original factored matrix.

Given a QR factorization of a real or complex matrix A = Q*R, Q unitary and R upper trapezoidal, return the QR factorization of [A(:,1:j-1), U, A(:,j:n)], where u is a column vector to be inserted into A (if orient is "col"), or the QR factorization of [A(1:j-1,:);X;A(:,j:n)], where x is a row orient is "row"). The default value of orient is "col".

If orient is "col", j may be an index vector resulting in the QR factorization of a matrix B such that A(:,j) = [] gives B. Notice that the latter case is done as a sequence of k deletions; thus, for k large enough, it will be both faster and more accurate to recompute the factorization from scratch.

If orient is "col", the QR factorization supplied may be either full (Q is square) or economized (R is square).

If orient is "row", full factorization is needed.

See also: qr, qrupdate, qrinsert, qrshift.

 
: [Q1, R1] = qrshift (Q, R, i, j)

Update a QR factorization given a range of columns to shift in the original factored matrix.

Given a QR factorization of a real or complex matrix A = Q*R, Q unitary and R upper trapezoidal, return the QR factorization of A(:,p), where p is the permutation
p = [1:i-1, shift(i:j, 1), j+1:n] if i < j
or
p = [1:j-1, shift(j:i,-1), i+1:n] if j < i.

See also: qr, qrupdate, qrinsert, qrdelete.

 
: [AA, BB, Q, Z, V, W] = qz (A, B)
: [AA, BB, Q, Z, V, W] = qz (A, B, opt)

Compute the QZ decomposition of a generalized eigenvalue problem.

The generalized eigenvalue problem is defined as

A x = lambda B x

There are two calling forms of the function:

  1. [AA, BB, Q, Z, V, W, lambda] = qz (A, B)

    Compute the complex QZ decomposition, generalized eigenvectors, and generalized eigenvalues.

    
    AA = Q * A * Z, BB = Q * B * Z
    A * V * diag (diag (BB)) = B * V * diag (diag (AA))
    diag (diag (BB)) * W' * A = diag (diag (AA)) * W' * B
    
    

    with AA and BB upper triangular, and Q and Z unitary. The matrices V and W respectively contain the right and left generalized eigenvectors.

  2. [AA, BB, Z {, lambda}] = qz (A, B, opt)

    The opt argument must be equal to either "real" or "complex". If it is equal to "complex", then this calling form is equivalent to the first one with only two input arguments.

    If opt is equal to "real", then the real QZ decomposition is computed. In particular, AA is only guaranteed to be quasi-upper triangular with 1-by-1 and 2-by-2 blocks on the diagonal, and Q and Z are orthogonal. The identities mentioned above for right and left generalized eigenvectors are only verified if AA is upper triangular (i.e., when all the generalized eigenvalues are real, in which case the real and complex QZ coincide).

Note: qz performs permutation balancing, but not scaling (see balance), which may be lead to less accurate results than eig. The order of output arguments was selected for compatibility with MATLAB.

See also: eig, gsvd, balance, chol, hess, lu, qr, qzhess, schur.

 
: [aa, bb, q, z] = qzhess (A, B)

Compute the Hessenberg-triangular decomposition of the matrix pencil (A, B), returning aa = q * A * z, bb = q * B * z, with q and z orthogonal.

For example:

[aa, bb, q, z] = qzhess ([1, 2; 3, 4], [5, 6; 7, 8])
  ⇒ aa =
      -3.02244  -4.41741
       0.92998   0.69749
  ⇒ bb =
      -8.60233  -9.99730
       0.00000  -0.23250
  ⇒ q =
      -0.58124  -0.81373
      -0.81373   0.58124
  ⇒ z =
     Diagonal Matrix
       1   0
       0   1

The Hessenberg-triangular decomposition is the first step in Moler and Stewart’s QZ decomposition algorithm.

Algorithm taken from Golub and Van Loan, Matrix Computations, 2nd edition.

See also: lu, chol, hess, qr, qz, schur, svd.

 
: S = schur (A)
: S = schur (A, "real")
: S = schur (A, "complex")
: S = schur (A, opt)
: [U, S] = schur (…)

Compute the Schur decomposition of A.

The Schur decomposition of a square matrix A is defined as

S = U' * A * U

where U is a unitary matrix (U'* U is identity) and S is upper triangular. The eigenvalues of A (and S) are the diagonal elements of S. If the matrix A is real, then the real Schur decomposition is computed, in which the matrix U is orthogonal and S is block upper triangular with blocks of size at most 2 x 2 along the diagonal.

The default for real matrices is a real Schur decomposition. A complex decomposition may be forced by passing the flag "complex".

The eigenvalues are optionally ordered along the diagonal according to the value of opt:

opt = "a"

Move eigenvalues with negative real parts to the leading block of S. Mnemonic: "a" for Algebraic Riccati Equations, where this ordering is useful.

opt = "d"

Move eigenvalues with magnitude less than one to the leading block of S. Mnemonic: "d" for Discrete Algebraic Riccati Equations, where this ordering is useful.

opt = "u"

Unordered. No particular ordering of eigenvalues (default).

The leading k columns of U always span the A-invariant subspace corresponding to the k leading eigenvalues of S.

See also: rsf2csf, ordschur, ordeig, lu, chol, hess, qr, qz, svd, eig.

 
: [U, T] = rsf2csf (UR, TR)

Convert a real, upper quasi-triangular Schur form TR to a complex, upper triangular Schur form T.

Note that the following relations hold:

UR * TR * UR' = U * T * U' and U' * U is the identity matrix I.

Note also that U and T are not unique.

See also: schur.

 
: [UR, SR] = ordschur (U, S, select)

Reorder the real Schur factorization (U,S) obtained with the schur function, so that selected eigenvalues appear in the upper left diagonal blocks of the quasi triangular Schur matrix.

The logical vector select specifies the selected eigenvalues as they appear along S’s diagonal.

For example, given the matrix A = [1, 2; 3, 4], and its Schur decomposition

[U, S] = schur (A)

which returns

U =

  -0.82456  -0.56577
   0.56577  -0.82456

S =

  -0.37228  -1.00000
   0.00000   5.37228

It is possible to reorder the decomposition so that the positive eigenvalue is in the upper left corner, by doing:

[U, S] = ordschur (U, S, [0,1])

See also: schur, ordeig, ordqz.

 
: [AR, BR, QR, ZR] = ordqz (AA, BB, Q, Z, keyword)
: [AR, BR, QR, ZR] = ordqz (AA, BB, Q, Z, select)

Reorder the QZ decomposition of a generalized eigenvalue problem.

The generalized eigenvalue problem is defined as

A x = lambda B x

Its generalized Schur decomposition is computed using the qz algorithm:

[AA, BB, Q, Z] = qz (A, B)

where AA, BB, Q, and Z fulfill


AA = Q * A * Z, BB = Q * B * Z

The ordqz function computes a unitary transformation QR and ZR such that the order of the eigenvalue on the diagonal of AA and BB is changed. The resulting reordered matrices AR and BR fulfill:


AR = QR * A * ZR, BR = QR * B * ZR

The function can either be called with the keyword argument which selects the eigenvalues in the top left block of AR and BR in the following way:

"S", "udi"

small: leading block has all |lambda| < 1

"B", "udo"

big: leading block has all |lambda| ≥ 1

"-", "lhp"

negative real part: leading block has all eigenvalues in the open left half-plane

"+", "rhp"

non-negative real part: leading block has all eigenvalues in the closed right half-plane

If a logical vector select is given instead of a keyword the ordqz function reorders all eigenvalues k to the left block for which select(k) is true.

Note: The keywords are compatible with the ones from qr.

See also: eig, ordeig, qz, schur, ordschur.

 
: lambda = ordeig (A)
: lambda = ordeig (A, B)

Return the eigenvalues of quasi-triangular matrices in their order of appearance in the matrix A.

The quasi-triangular matrix A is usually the result of a Schur factorization. If called with a second input B then the generalized eigenvalues of the pair A, B are returned in the order of appearance of the matrix A-lambda*B. The pair A, B is usually the result of a QZ decomposition.

See also: ordschur, ordqz, eig, schur, qz.

 
: angle = subspace (A, B)

Determine the largest principal angle between two subspaces spanned by the columns of matrices A and B.

 
: s = svd (A)
: [U, S, V] = svd (A)
: [U, S, V] = svd (A, "econ")
: [U, S, V] = svd (A, 0)

Compute the singular value decomposition of A.

The singular value decomposition is defined by the relation

A = U*S*V'

The function svd normally returns only the vector of singular values. When called with three return values, it computes U, S, and V. For example,

svd (hilb (3))

returns

ans =

  1.4083189
  0.1223271
  0.0026873

and

[u, s, v] = svd (hilb (3))

returns

u =

  -0.82704   0.54745   0.12766
  -0.45986  -0.52829  -0.71375
  -0.32330  -0.64901   0.68867

s =

  1.40832  0.00000  0.00000
  0.00000  0.12233  0.00000
  0.00000  0.00000  0.00269

v =

  -0.82704   0.54745   0.12766
  -0.45986  -0.52829  -0.71375
  -0.32330  -0.64901   0.68867

When given a second argument that is not 0, svd returns an economy-sized decomposition, eliminating the unnecessary rows or columns of U or V.

If the second argument is exactly 0, then the choice of decomposition is based on the matrix A. If A has more rows than columns then an economy-sized decomposition is returned, otherwise a regular decomposition is calculated.

Algorithm Notes: When calculating the full decomposition (left and right singular matrices in addition to singular values) there is a choice of two routines in LAPACK. The default routine used by Octave is gesvd. The alternative is gesdd which is 5X faster, but may use more memory and may be inaccurate for some input matrices. There is a third routine gejsv, suitable for better accuracy at extreme scale. See the documentation for svd_driver for more information on choosing a driver.

See also: svd_driver, svds, eig, lu, chol, hess, qr, qz.

 
: val = svd_driver ()
: old_val = svd_driver (new_val)
: old_val = svd_driver (new_val, "local")

Query or set the underlying LAPACK driver used by svd.

Currently recognized values are "gesdd", "gesvd", and "gejsv". The default is "gesvd".

When called from inside a function with the "local" option, the variable is changed locally for the function and any subroutines it calls. The original variable value is restored when exiting the function.

Algorithm Notes: The LAPACK library routines gesvd and gesdd are different only when calculating the full singular value decomposition (left and right singular matrices as well as singular values). When calculating just the singular values the following discussion is not relevant.

The newer gesdd routine is based on a Divide-and-Conquer algorithm that is 5X faster than the alternative gesvd, which is based on QR factorization. However, the new algorithm can use significantly more memory. For an MxN input matrix the memory usage is of order O(min(M,N) ^ 2), whereas the alternative is of order O(max(M,N)).

The routine gejsv uses a preconditioned Jacobi SVD algorithm. Unlike gesvd and gesdd, in gejsv, there is no bidiagonalization step that could contaminate accuracy in some extreme cases. Also, gejsv is known to be optimally accurate in some sense. However, the speed is slower (single threaded at its core) and uses more memory (O(min(M,N) ^ 2 + M + N)).

Beyond speed and memory issues, there have been instances where some input matrices were not accurately decomposed by gesdd. See currently active bug https://savannah.gnu.org/bugs/?55564. Until these accuracy issues are resolved in a new version of the LAPACK library, the default driver in Octave has been set to "gesvd".

See also: svd.

 
: [housv, beta, zer] = housh (x, j, z)

Compute Householder reflection vector housv to reflect x to be the j-th column of identity, i.e.,

(I - beta*housv*housv')x =  norm (x)*e(j) if x(j) < 0,
(I - beta*housv*housv')x = -norm (x)*e(j) if x(j) >= 0

Inputs

x

vector

j

index into vector

z

threshold for zero (usually should be the number 0)

Outputs (see Golub and Van Loan):

beta

If beta = 0, then no reflection need be applied (zer set to 0)

housv

householder vector

 
: [u, h, nu] = krylov (A, V, k, eps1, pflg)

Construct an orthogonal basis u of a block Krylov subspace.

The block Krylov subspace has the following form:

[v a*v a^2*v ... a^(k+1)*v]

The construction is made with Householder reflections to guard against loss of orthogonality.

If V is a vector, then h contains the Hessenberg matrix such that a*u == u*h+rk*ek', in which rk = a*u(:,k)-u*h(:,k), and ek' is the vector [0, 0, …, 1] of length k. Otherwise, h is meaningless.

If V is a vector and k is greater than length (A) - 1, then h contains the Hessenberg matrix such that a*u == u*h.

The value of nu is the dimension of the span of the Krylov subspace (based on eps1).

If b is a vector and k is greater than m-1, then h contains the Hessenberg decomposition of A.

The optional parameter eps1 is the threshold for zero. The default value is 1e-12.

If the optional parameter pflg is nonzero, row pivoting is used to improve numerical behavior. The default value is 0.

Reference: A. Hodel, P. Misra, Partial Pivoting in the Computation of Krylov Subspaces of Large Sparse Systems, Proceedings of the 42nd IEEE Conference on Decision and Control, December 2003.