#include <stdio.h>#include <stdlib.h>#include <string.h>#include <math.h>#include "include/hholder.h"#include "include/qr.h"Go to the source code of this file.
Functions | |
| int | qr (double **A, int m, int n, double *B, double *X, double *rnorm, double *tau, double *res, double **wws, double *ws) |
| int | qr_decomp (double **a, int M, int N, double *tau, double **cchain, double *chain) |
| int | qr_solve (double **QR, int M, int N, double *tau, double *b, double *x, double *residual, double *resNorm, double **cchain, double *chain) |
| int | qr_QTvec (double **QR, double *tau, double *v, double *help) |
| int | qr_Qvec (double **QR, double *tau, double *v, double *help) |
| int | qr_invRvec (double **A, double *X) |
| int | qr_weight (int N, int M, double **A, double *b, double *weight) |
Variables | |
| int | qr_M |
| int | qr_N |
| int | qr_MNmin |
| int qr | ( | double ** | A, |
| int | m, | ||
| int | n, | ||
| double * | B, | ||
| double * | X, | ||
| double * | rnorm, | ||
| double * | tau, | ||
| double * | res, | ||
| double ** | wws, | ||
| double * | ws | ||
| ) |
| A | On entry, a[m][n] contains the m by n matrix A. On exit, a[][] contains the QR factorization |
| m | Dimensions of matrix A are a[m][n] |
| n | Dimensions of matrix A are a[m][n] |
| B | B[] is an m-vector containing the right-hand side vector b. |
| X | On exit, x[] will contain the solution vector x. |
| rnorm | On exit, rnorm contains the Euclidean norm of the residual vector |
| tau | On exit, tau[] will contain the householder coefficients. |
| res | An m-array of working space, res[]. On output contains residual b - Ax |
| wws | m*n array of working space |
| ws | 2m-array of working space |
Definition at line 52 of file qr.c.
References qr_decomp(), and qr_solve().
| int qr_decomp | ( | double ** | a, |
| int | M, | ||
| int | N, | ||
| double * | tau, | ||
| double ** | cchain, | ||
| double * | chain | ||
| ) |
Factorise a general M x N matrix A into
A = Q R
where Q is orthogonal (M x M) and R is upper triangular (M x N).
Q is stored as a packed set of Householder vectors in the strict lower triangular part of the input matrix A and a set of coefficients in vector tau.
R is stored in the diagonal and upper triangle of the input matrix.
The full matrix for Q can be obtained as the product
Q = Q_1 Q_2 .. Q_k and it's transform as the product Q^T = Q_k .. Q_2 Q_1
where k = min(M,N) and
Q_i = (I - tau_i * h_i * h_i^T)
and where h_i is a Householder vector
h_i = [1, A(i+1,i), A(i+2,i), ... , A(M,i)]
This storage scheme is the same as in LAPACK.
NOTICE! The calling program must take care that pointer tau is of size n.
| a | contains coefficient matrix A (m*n) as input and factorisation QR as output |
| M | nr of rows in matrix A |
| N | nr of columns in matrix A |
| tau | n-vector for householder coefficients |
| cchain | m*n array of working space |
| chain | 2m array of working space |
Definition at line 171 of file qr.c.
References householder_hm(), householder_transform(), M, qr_M, qr_MNmin, and qr_N.
Referenced by qr().
| int qr_invRvec | ( | double ** | A, |
| double * | X | ||
| ) |
Form the product R^-1 v. (R is saved in the upper triangel of QR matrix.)
Definition at line 409 of file qr.c.
References qr_N.
Referenced by qr_solve().
| int qr_QTvec | ( | double ** | QR, |
| double * | tau, | ||
| double * | v, | ||
| double * | help | ||
| ) |
Form the product Q^T v from householder vectors saved in the lower triangel of QR matrix and householder coefficients saved in vector tau.
Definition at line 345 of file qr.c.
References householder_hv(), qr_M, and qr_MNmin.
Referenced by qr_solve().
| int qr_Qvec | ( | double ** | QR, |
| double * | tau, | ||
| double * | v, | ||
| double * | help | ||
| ) |
Form the product Q v from householder vectors saved in the lower triangel of QR matrix and householder coefficients saved in vector tau.
Definition at line 379 of file qr.c.
References householder_hv(), qr_M, and qr_MNmin.
Referenced by qr_solve().
| int qr_solve | ( | double ** | QR, |
| int | M, | ||
| int | N, | ||
| double * | tau, | ||
| double * | b, | ||
| double * | x, | ||
| double * | residual, | ||
| double * | resNorm, | ||
| double ** | cchain, | ||
| double * | chain | ||
| ) |
Find the least squares solution to the overdetermined system
A x = b
for m >= n using the QR factorisation A = Q R. qr_decomp() must be used prior to this function in order to form the QR factorisation of A. Solution is formed in the following order: QR x = b => R x = Q^T b => x = R^-1 (Q^T b)
NOTICE! The calling program must take care that pointers b, x and residual are of the right size.
| QR | m*n matrix containing householder vectors of A |
| M | nr of rows in matrix A |
| N | nr of columns in matrix A |
| tau | vector containing householder coefficients tau |
| b | Contains m-vector b of A x = b |
| x | solution vector x of length n |
| residual | residual vector of length m |
| resNorm | norm^2 of the residual vector |
| cchain | m*n array of the working space |
| chain | 2m lenght array of the working space |
Definition at line 258 of file qr.c.
References M, qr_invRvec(), qr_M, qr_MNmin, qr_N, qr_QTvec(), and qr_Qvec().
Referenced by qr().
| int qr_weight | ( | int | N, |
| int | M, | ||
| double ** | A, | ||
| double * | b, | ||
| double * | weight | ||
| ) |
| int qr_M |
Global variables for this routine
Definition at line 34 of file qr.c.
Referenced by qr_decomp(), qr_QTvec(), qr_Qvec(), and qr_solve().
| int qr_MNmin |
Definition at line 36 of file qr.c.
Referenced by qr_decomp(), qr_QTvec(), qr_Qvec(), and qr_solve().
| int qr_N |
Definition at line 35 of file qr.c.
Referenced by qr_decomp(), qr_invRvec(), and qr_solve().
1.8.0