Go to the first, previous, next, last section, table of contents.
A general rectangular @math{M}-by-@math{N} matrix @math{A} has a
@math{QR} decomposition into the product of an orthogonal
@math{M}-by-@math{M} square matrix @math{Q} (where @math{Q^T Q = I}) and
an @math{M}-by-@math{N} right-triangular matrix @math{R},
This decomposition can be used to convert the linear system @math{A x =
b} into the triangular system @math{R x = Q^T b}, which can be solved by
back-substitution. Another use of the @math{QR} decomposition is to
compute an orthonormal basis for a set of vectors. The first @math{N}
columns of @math{Q} form an orthonormal basis for the range of @math{A},
@math{ran(A)}, when @math{A} has full column rank.
- Function: int gsl_linalg_QR_decomp (gsl_matrix * A, gsl_vector * tau)
-
This function factorizes the @math{M}-by-@math{N} matrix A into
the @math{QR} decomposition @math{A = Q R}. On output the diagonal and
upper triangular part of the input matrix contain the matrix
@math{R}. The vector tau and the columns of the lower triangular
part of the matrix A contain the Householder coefficients and
Householder vectors which encode the orthogonal matrix Q. The
vector tau must be of length @math{k=\min(M,N)}. The matrix
@math{Q} is related to these components by, @math{Q = Q_k ... Q_2 Q_1}
where @math{Q_i = I - \tau_i v_i v_i^T} and @math{v_i} is the
Householder vector @math{v_i =
(0,...,1,A(i+1,i),A(i+2,i),...,A(m,i))}. This is the same storage scheme
as used by LAPACK.
The algorithm used to perform the decomposition is Householder QR (Golub
& Van Loan, Matrix Computations, Algorithm 5.2.1).
- Function: int gsl_linalg_QR_solve (const gsl_matrix * QR, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x)
-
This function solves the system @math{A x = b} using the @math{QR}
decomposition of @math{A} into (QR, tau) given by
gsl_linalg_QR_decomp
.
- Function: int gsl_linalg_QR_svx (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * x)
-
This function solves the system @math{A x = b} in-place using the
@math{QR} decomposition of @math{A} into (QR,tau) given by
gsl_linalg_QR_decomp
. On input x should contain the
right-hand side @math{b}, which is replaced by the solution on output.
- Function: int gsl_linalg_QR_lssolve (const gsl_matrix * QR, const gsl_vector * tau, const gsl_vector * b, gsl_vector * x, gsl_vector * residual)
-
This function finds the least squares solution to the overdetermined
system @math{A x = b} where the matrix A has more rows than
columns. The least squares solution minimizes the Euclidean norm of the
residual, @math{||Ax - b||}.The routine uses the @math{QR} decomposition
of @math{A} into (QR, tau) given by
gsl_linalg_QR_decomp
. The solution is returned in x. The
residual is computed as a by-product and stored in residual.
- Function: int gsl_linalg_QR_QTvec (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * v)
-
This function applies the matrix @math{Q^T} encoded in the decomposition
(QR,tau) to the vector v, storing the result @math{Q^T
v} in v. The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix @math{Q^T}.
- Function: int gsl_linalg_QR_Qvec (const gsl_matrix * QR, const gsl_vector * tau, gsl_vector * v)
-
This function applies the matrix @math{Q} encoded in the decomposition
(QR,tau) to the vector v, storing the result @math{Q
v} in v. The matrix multiplication is carried out directly using
the encoding of the Householder vectors without needing to form the full
matrix @math{Q}.
- Function: int gsl_linalg_QR_Rsolve (const gsl_matrix * QR, const gsl_vector * b, gsl_vector * x)
-
This function solves the triangular system @math{R x = b} for
x. It may be useful if the product @math{b' = Q^T b} has already
been computed using
gsl_linalg_QR_QTvec
.
- Function: int gsl_linalg_QR_Rsvx (const gsl_matrix * QR, gsl_vector * x)
-
This function solves the triangular system @math{R x = b} for x
in-place. On input x should contain the right-hand side @math{b}
and is replaced by the solution on output. This function may be useful if
the product @math{b' = Q^T b} has already been computed using
gsl_linalg_QR_QTvec
.
- Function: int gsl_linalg_QR_unpack (const gsl_matrix * QR, const gsl_vector * tau, gsl_matrix * Q, gsl_matrix * R)
-
This function unpacks the encoded @math{QR} decomposition
(QR,tau) into the matrices Q and R, where
Q is @math{M}-by-@math{M} and R is @math{M}-by-@math{N}.
- Function: int gsl_linalg_QR_QRsolve (gsl_matrix * Q, gsl_matrix * R, const gsl_vector * b, gsl_vector * x)
-
This function solves the system @math{R x = Q^T b} for x. It can
be used when the @math{QR} decomposition of a matrix is available in
unpacked form as (Q,R).
- Function: int gsl_linalg_QR_update (gsl_matrix * Q, gsl_matrix * R, gsl_vector * w, const gsl_vector * v)
-
This function performs a rank-1 update @math{w v^T} of the @math{QR}
decomposition (Q, R). The update is given by @math{Q'R' = Q
R + w v^T} where the output matrices @math{Q'} and @math{R'} are also
orthogonal and right triangular. Note that w is destroyed by the
update.
- Function: int gsl_linalg_R_solve (const gsl_matrix * R, const gsl_vector * b, gsl_vector * x)
-
This function solves the triangular system @math{R x = b} for the
@math{N}-by-@math{N} matrix R.
- Function: int gsl_linalg_R_svx (const gsl_matrix * R, gsl_vector * x)
-
This function solves the triangular system @math{R x = b} in-place. On
input x should contain the right-hand side @math{b}, which is
replaced by the solution on output.
Go to the first, previous, next, last section, table of contents.