linalg 1.8.2
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
|
#include <stdbool.h>
#include <complex.h>
Go to the source code of this file.
Macros | |
#define | LA_NO_OPERATION 0 |
#define | LA_TRANSPOSE 1 |
#define | LA_HERMITIAN_TRANSPOSE 2 |
#define | LA_NO_ERROR 0 |
#define | LA_INVALID_INPUT_ERROR 101 |
#define | LA_ARRAY_SIZE_ERROR 102 |
#define | LA_SINGULAR_MATRIX_ERROR 103 |
#define | LA_MATRIX_FORMAT_ERROR 104 |
#define | LA_OUT_OF_MEMORY_ERROR 105 |
#define | LA_CONVERGENCE_ERROR 106 |
#define | LA_INVALID_OPERATION_ERROR 107 |
Functions | |
int | la_rank1_update (int m, int n, double alpha, const double *x, const double *y, double *a, int lda) |
int | la_rank1_update_cmplx (int m, int n, double complex alpha, const double complex *x, const double complex *y, double complex *a, int lda) |
int | la_trace (int m, int n, const double *a, int lda, double *rst) |
int | la_trace_cmplx (int m, int n, const double complex *a, int lda, double complex *rst) |
int | la_mtx_mult (bool transa, bool transb, int m, int n, int k, double alpha, const double *a, int lda, const double *b, int ldb, double beta, double *c, int ldc) |
int | la_mtx_mult_cmplx (int opa, int opb, int m, int n, int k, double complex alpha, const double complex *a, int lda, const double complex *b, int ldb, double complex beta, double complex *c, int ldc) |
int | la_diag_mtx_mult (bool lside, bool transb, int m, int n, int k, double alpha, const double *a, const double *b, int ldb, double beta, double *c, int ldc) |
int | la_diag_mtx_mult_cmplx (bool lside, int opb, int m, int n, int k, double complex alpha, const double complex *a, const double complex *b, int ldb, double complex beta, double complex *c, int ldc) |
int | la_diag_mtx_mult_mixed (bool lside, int opb, int m, int n, int k, double complex alpha, const double *a, const double complex *b, int ldb, double complex beta, double complex *c, int ldc) |
int | la_rank (int m, int n, double *a, int lda, int *rnk) |
int | la_rank_cmplx (int m, int n, double complex *a, int lda, int *rnk) |
int | la_det (int n, double *a, int lda, double *d) |
int | la_det_cmplx (int n, double complex *a, int lda, double complex *d) |
int | la_tri_mtx_mult (bool upper, double alpha, int n, const double *a, int lda, double beta, double *b, int ldb) |
int | la_tri_mtx_mult_cmplx (bool upper, double complex alpha, int n, const double complex *a, int lda, double complex beta, double complex *b, int ldb) |
int | la_lu_factor (int m, int n, double *a, int lda, int *ipvt) |
int | la_lu_factor_cmplx (int m, int n, double complex *a, int lda, int *ipvt) |
int | la_form_lu (int n, double *a, int lda, int *ipvt, double *u, int ldu, double *p, int ldp) |
int | la_form_lu_cmplx (int n, double complex *a, int lda, int *ipvt, double complex *u, int ldu, double *p, int ldp) |
int | la_qr_factor (int m, int n, double *a, int lda, double *tau) |
int | la_qr_factor_cmplx (int m, int n, double complex *a, int lda, double complex *tau) |
int | la_qr_factor_pvt (int m, int n, double *a, int lda, double *tau, int *jpvt) |
int | la_qr_factor_cmplx_pvt (int m, int n, double complex *a, int lda, double complex *tau, int *jpvt) |
int | la_form_qr (bool fullq, int m, int n, double *r, int ldr, const double *tau, double *q, int ldq) |
int | la_form_qr_cmplx (bool fullq, int m, int n, double complex *r, int ldr, const double complex *tau, double complex *q, int ldq) |
int | la_form_qr_pvt (bool fullq, int m, int n, double *r, int ldr, const double *tau, const int *pvt, double *q, int ldq, double *p, int ldp) |
int | la_form_qr_cmplx_pvt (bool fullq, int m, int n, double complex *r, int ldr, const double complex *tau, const int *pvt, double complex *q, int ldq, double complex *p, int ldp) |
int | la_mult_qr (bool lside, bool trans, int m, int n, int k, double *a, int lda, const double *tau, double *c, int ldc) |
int | la_mult_qr_cmplx (bool lside, bool trans, int m, int n, int k, double complex *a, int lda, const double complex *tau, double complex *c, int ldc) |
int | la_qr_rank1_update (int m, int n, double *q, int ldq, double *r, int ldr, double *u, double *v) |
int | la_cholesky_factor (bool upper, int n, double *a, int lda) |
int | la_cholesky_factor_cmplx (bool upper, int n, double complex *a, int lda) |
int | la_cholesky_rank1_update (int n, double *r, int ldr, double *u) |
int | la_cholesky_rank1_downdate (int n, double *r, int ldr, double *u) |
int | la_cholesky_rank1_downdate_cmplx (int n, double complex *r, int ldr, double complex *u) |
int | la_svd (int m, int n, double *a, int lda, double *s, double *u, int ldu, double *vt, int ldv) |
int | la_svd_cmplx (int m, int n, double complex *a, int lda, double *s, double complex *u, int ldu, double complex *vt, int ldv) |
int | la_solve_tri_mtx (bool lside, bool upper, bool trans, bool nounit, int m, int n, double alpha, const double *a, int lda, double *b, int ldb) |
int | la_solve_tri_mtx_cmplx (bool lside, bool upper, bool trans, bool nounit, int m, int n, double complex alpha, const double complex *a, int lda, double complex *b, int ldb) |
int | la_solve_lu (int m, int n, const double *a, int lda, const int *ipvt, double *b, int ldb) |
int | la_solve_lu_cmplx (int m, int n, const double complex *a, int lda, const int *ipvt, double complex *b, int ldb) |
int | la_solve_qr (int m, int n, int k, double *a, int lda, const double *tau, double *b, int ldb) |
int | la_solve_qr_cmplx (int m, int n, int k, double complex *a, int lda, const double complex *tau, double complex *b, int ldb) |
int | la_solve_qr_pvt (int m, int n, int k, double *a, int lda, const double *tau, const int *jpvt, double *b, int ldb) |
int | la_solve_qr_cmplx_pvt (int m, int n, int k, double complex *a, int lda, const double complex *tau, const int *jpvt, double complex *b, int ldb) |
int | la_solve_cholesky (bool upper, int m, int n, const double *a, int lda, double *b, int ldb) |
int | la_solve_cholesky_cmplx (bool upper, int m, int n, const double complex *a, int lda, double complex *b, int ldb) |
int | la_solve_least_squares (int m, int n, int k, double *a, int lda, double *b, int ldb) |
int | la_solve_least_squares_cmplx (int m, int n, int k, double complex *a, int lda, double complex *b, int ldb) |
int | la_inverse (int n, double *a, int lda) |
int | la_inverse_cmplx (int n, double complex *a, int lda) |
int | la_pinverse (int m, int n, double *a, int lda, double *ainv, int ldai) |
int | la_pinverse_cmplx (int m, int n, double complex *a, int lda, double complex *ainv, int ldai) |
int | la_eigen_symm (bool vecs, int n, double *a, int lda, double *vals) |
int | la_eigen_asymm (bool vecs, int n, double *a, int lda, double complex *vals, double complex *v, int ldv) |
int | la_eigen_gen (bool vecs, int n, double *a, int lda, double *b, int ldb, double complex *alpha, double *beta, double complex *v, int ldv) |
int | la_eigen_cmplx (bool vecs, int n, double complex *a, int lda, double complex *vals, double complex *v, int ldv) |
int | la_sort_eigen (bool ascend, int n, double *vals, double *vecs, int ldv) |
int | la_sort_eigen_cmplx (bool ascend, int n, double complex *vals, double complex *vecs, int ldv) |
int | la_lq_factor (int m, int n, double *a, int lda, double *tau) |
int | la_form_lq (int m, int n, double *l, int ldl, const double *tau, double *q, int ldq) |
int | la_mult_lq (bool lside, bool trans, int m, int n, int k, const double *a, int lda, const double *tau, double *c, int ldc) |
int | la_solve_lq (int m, int n, int k, const double *a, int lda, const double *tau, double *b, int ldb) |
int | la_qr_rank1_update_cmplx (int m, int n, double complex *q, int ldq, double complex *r, int ldr, double complex *u, double complex *v) |
int | la_cholesky_rank1_update_cmplx (int n, double complex *r, int ldr, double complex *u) |
int | la_lq_factor_cmplx (int m, int n, double complex *a, int lda, double complex *tau) |
int | la_form_lq_cmplx (int m, int n, double complex *l, int ldl, const double complex *tau, double complex *q, int ldq) |
int | la_mult_lq_cmplx (bool lside, bool trans, int m, int n, int k, const double complex *a, int lda, const double complex *tau, double complex *c, int ldc) |
int | la_solve_lq_cmplx (int m, int n, int k, const double complex *a, int lda, const double complex *tau, double complex *b, int ldb) |
int | la_band_mtx_vec_mult (bool trans, int m, int n, int kl, int ku, double alpha, const double *a, int lda, const double *x, double beta, double *y) |
int | la_band_mtx_vec_mult_cmplx (int trans, int m, int n, int kl, int ku, double complex alpha, const double complex *a, int lda, const double complex *x, double complex beta, double complex *y) |
int | la_band_to_full_mtx (int m, int n, int kl, int ku, const double *b, int ldb, double *f, int ldf) |
int | la_band_to_full_mtx_cmplx (int m, int n, int kl, int ku, const double complex *b, int ldb, double complex *f, int ldf) |
int | la_band_diag_mtx_mult (bool left, int m, int n, int kl, int ku, double alpha, double *a, int lda, const double *b) |
int | la_band_diag_mtx_mult_cmplx (bool left, int m, int n, int kl, int ku, double complex alpha, double complex *a, int lda, const double complex *b) |
int la_band_diag_mtx_mult | ( | bool | left, |
int | m, | ||
int | n, | ||
int | kl, | ||
int | ku, | ||
double | alpha, | ||
double * | a, | ||
int | lda, | ||
const double * | b ) |
Multiplies a banded matrix, A, with a diagonal matrix, B, such that A = alpha * A * B, or A = alpha * B * A.
left | Set to true to compute A = alpha * A * B; else, set to false to compute A = alpha * B * A. |
m | The number of rows in matrix A. |
kl | The number of subdiagonals. Must be at least 0. |
ku | The number of superdiagonals. Must be at least 0. |
alpha | The scalar multiplier. |
a | The M-by-N matrix A storing the banded matrix in a compressed form supplied column by column. The following code segment transfers between a full matrix to the banded matrix storage scheme. do j = 1, n
k = ku + 1 - j
do i = max(1, j - ku), min(m, j + kl)
a(k + i, j) = matrix(i, j)
end do
end do
|
b | An array containing the diagonal elements of matrix B. |
b
and f
are not compatible in size.a
and b
are not compatible in terms of internal dimensions.ku
or kl
are not zero or greater. int la_band_diag_mtx_mult_cmplx | ( | bool | left, |
int | m, | ||
int | n, | ||
int | kl, | ||
int | ku, | ||
double complex | alpha, | ||
double complex * | a, | ||
int | lda, | ||
const double complex * | b ) |
Multiplies a banded matrix, A, with a diagonal matrix, B, such that A = alpha * A * B, or A = alpha * B * A.
left | Set to true to compute A = alpha * A * B; else, set to false to compute A = alpha * B * A. |
m | The number of rows in matrix A. |
kl | The number of subdiagonals. Must be at least 0. |
ku | The number of superdiagonals. Must be at least 0. |
alpha | The scalar multiplier. |
a | The M-by-N matrix A storing the banded matrix in a compressed form supplied column by column. The following code segment transfers between a full matrix to the banded matrix storage scheme. do j = 1, n
k = ku + 1 - j
do i = max(1, j - ku), min(m, j + kl)
a(k + i, j) = matrix(i, j)
end do
end do
|
b | An array containing the diagonal elements of matrix B. |
b
and f
are not compatible in size.a
and b
are not compatible in terms of internal dimensions.ku
or kl
are not zero or greater. int la_band_mtx_vec_mult | ( | bool | trans, |
int | m, | ||
int | n, | ||
int | kl, | ||
int | ku, | ||
double | alpha, | ||
const double * | a, | ||
int | lda, | ||
const double * | x, | ||
double | beta, | ||
double * | y ) |
Multiplies a banded matrix, A, by a vector x such that alpha * op(A) * x + beta * y = y.
trans | Set to true for op(A) == A^T; else, false for op(A) == A. |
kl | The number of subdiagonals. Must be at least 0. |
ku | The number of superdiagonals. Must be at least 0. |
alpha | A scalar multiplier. |
a | The M-by-N matrix A storing the banded matrix in a compressed form supplied column by column. |
x | If trans is true, this is an M-element vector; else, if trans is false, this is an N-element vector. |
beta | A scalar multiplier. |
y | On input, the vector Y. On output, the resulting vector. if trans is true, this vector is an N-element vector; else, it is an M-element vector. |
ku
or kl
are not zero or greater. int la_band_mtx_vec_mult_cmplx | ( | int | trans, |
int | m, | ||
int | n, | ||
int | kl, | ||
int | ku, | ||
double complex | alpha, | ||
const double complex * | a, | ||
int | lda, | ||
const double complex * | x, | ||
double complex | beta, | ||
double complex * | y ) |
Multiplies a banded matrix, A, by a vector x such that alpha * op(A) * x + beta * y = y.
trans | Set to LA_TRANSPOSE if op(A) = A^T, set to LA_HERMITIAN_TRANSPOSE if op(A) = A^H, otherwise set to LA_NO_OPERATION if op(A) = A. |
kl | The number of subdiagonals. Must be at least 0. |
ku | The number of superdiagonals. Must be at least 0. |
alpha | A scalar multiplier. |
a | The M-by-N matrix A storing the banded matrix in a compressed form supplied column by column. |
x | If trans is true, this is an M-element vector; else, if trans is false, this is an N-element vector. |
beta | A scalar multiplier. |
y | On input, the vector Y. On output, the resulting vector. if trans is true, this vector is an N-element vector; else, it is an M-element vector. |
ku
or kl
are not zero or greater. int la_band_to_full_mtx | ( | int | m, |
int | n, | ||
int | kl, | ||
int | ku, | ||
const double * | b, | ||
int | ldb, | ||
double * | f, | ||
int | ldf ) |
Converts a banded matrix stored in dense form to a full matrix.
kl | The number of subdiagonals. Must be at least 0. |
ku | The number of superdiagonals. Must be at least 0. |
b | The banded matrix to convert, stored in dense form. See band_mtx_vec_mult for details on this storage method. |
f | The M-by-N element full matrix. |
b
and f
are not compatible in size.ku
or kl
are not zero or greater. int la_band_to_full_mtx_cmplx | ( | int | m, |
int | n, | ||
int | kl, | ||
int | ku, | ||
const double complex * | b, | ||
int | ldb, | ||
double complex * | f, | ||
int | ldf ) |
Converts a banded matrix stored in dense form to a full matrix.
kl | The number of subdiagonals. Must be at least 0. |
ku | The number of superdiagonals. Must be at least 0. |
b | The banded matrix to convert, stored in dense form. See band_mtx_vec_mult for details on this storage method. |
f | The M-by-N element full matrix. |
b
and f
are not compatible in size.ku
or kl
are not zero or greater. int la_cholesky_factor | ( | bool | upper, |
int | n, | ||
double * | a, | ||
int | lda ) |
Computes the Cholesky factorization of a symmetric, positive definite matrix.
upper | Set to true to compute the upper triangular factoriztion \( A = U^T U \); else, set to false to compute the lower triangular factorzation \( A = L L^T \). |
n | The dimension of matrix A. |
a | On input, the N-by-N matrix to factor. On output, the factored matrix is returned in either the upper or lower triangular portion of the matrix, dependent upon the value of upper . |
lda | The leading dimension of matrix A. |
lda
is not correct.a
is not positive definite. int la_cholesky_factor_cmplx | ( | bool | upper, |
int | n, | ||
double complex * | a, | ||
int | lda ) |
Computes the Cholesky factorization of a symmetric, positive definite matrix.
upper | Set to true to compute the upper triangular factoriztion \( A = U^T U \); else, set to false to compute the lower triangular factorzation \( A = L L^T \). |
n | The dimension of matrix A. |
a | On input, the N-by-N matrix to factor. On output, the factored matrix is returned in either the upper or lower triangular portion of the matrix, dependent upon the value of upper . |
lda | The leading dimension of matrix A. |
lda
is not correct.a
is not positive definite. int la_cholesky_rank1_downdate | ( | int | n, |
double * | r, | ||
int | ldr, | ||
double * | u ) |
Computes the rank 1 downdate to a Cholesky factored matrix (upper triangular).
n | The dimension of the matrix. |
r | On input, the N-by-N upper triangular matrix R. On output, the updated matrix R1. |
ldr | The leading dimension of matrix R. |
u | On input, the N-element update vector U. On output, the rotation sines used to transform R to R1. |
ldr
is not correct.int la_cholesky_rank1_downdate_cmplx | ( | int | n, |
double complex * | r, | ||
int | ldr, | ||
double complex * | u ) |
Computes the rank 1 downdate to a Cholesky factored matrix (upper triangular).
n | The dimension of the matrix. |
r | On input, the N-by-N upper triangular matrix R. On output, the updated matrix R1. |
ldr | The leading dimension of matrix R. |
u | On input, the N-element update vector U. On output, the rotation sines used to transform R to R1. |
ldr
is not correct.int la_cholesky_rank1_update | ( | int | n, |
double * | r, | ||
int | ldr, | ||
double * | u ) |
Computes the rank 1 update to a Cholesky factored matrix (upper triangular).
n | The dimension of the matrix. |
r | On input, the N-by-N upper triangular matrix R. On output, the updated matrix R1. |
ldr | The leading dimension of matrix R. |
u | On input, the N-element update vector U. On output, the rotation sines used to transform R to R1. |
ldr
is not correct.int la_cholesky_rank1_update_cmplx | ( | int | n, |
double complex * | r, | ||
int | ldr, | ||
double complex * | u ) |
Computes the rank 1 update to a Cholesky factored matrix (upper triangular).
n | The dimension of the matrix. |
r | On input, the N-by-N upper triangular matrix R. On output, the updated matrix R1. |
ldr | The leading dimension of matrix R. |
u | On input, the N-element update vector U. On output, the rotation sines used to transform R to R1. |
ldr
is not correct.int la_det | ( | int | n, |
double * | a, | ||
int | lda, | ||
double * | d ) |
Computes the determinant of a square matrix.
n | The dimension of the matrix. |
a | The N-by-N matrix. The matrix is overwritten on output. |
lda | The leading dimension of the matrix. |
d | The determinant of a . |
int la_det_cmplx | ( | int | n, |
double complex * | a, | ||
int | lda, | ||
double complex * | d ) |
Computes the determinant of a square matrix.
n | The dimension of the matrix. |
a | The N-by-N matrix. The matrix is overwritten on output. |
lda | The leading dimension of the matrix. |
d | The determinant of a . |
int la_diag_mtx_mult | ( | bool | lside, |
bool | transb, | ||
int | m, | ||
int | n, | ||
int | k, | ||
double | alpha, | ||
const double * | a, | ||
const double * | b, | ||
int | ldb, | ||
double | beta, | ||
double * | c, | ||
int | ldc ) |
Computes the matrix operation: \( C = \alpha A op(B) + \beta C \), or \( C = \alpha op(B) A + \beta C \).
lside | Set to true to apply matrix A from the left; else, set to false to apply matrix A from the left. |
trans | Set to true if op(B) == B**T; else, set to false if op(B) == B. |
m | The number of rows in the matrix C. |
n | The number of columns in the matrix C. |
k | The inner dimension of the matrix product A * op(B). |
alpha | A scalar multiplier. |
a | A P-element array containing the diagonal elements of matrix A where P = MIN(m , k ) if lside is true; else, P = MIN(n , k ) if lside is false. |
b | The LDB-by-TDB matrix B where (LDB = leading dimension of B, and TDB = trailing dimension of B):
|
ldb | The leading dimension of matrix B. |
beta | A scalar multiplier. |
c | The m by n matrix C. |
ldc | The leading dimension of matrix C. |
ldb
, or ldc
are not correct.int la_diag_mtx_mult_cmplx | ( | bool | lside, |
int | opb, | ||
int | m, | ||
int | n, | ||
int | k, | ||
double complex | alpha, | ||
const double complex * | a, | ||
const double complex * | b, | ||
int | ldb, | ||
double complex | beta, | ||
double complex * | c, | ||
int | ldc ) |
Computes the matrix operation: \( C = \alpha A op(B) + \beta C \), or \( C = \alpha op(B) A + \beta * C \).
lside | Set to true to apply matrix A from the left; else, set to false to apply matrix A from the left. |
opb | Set to LA_TRANSPOSE to compute op(B) as a direct transpose of B, set to LA_HERMITIAN_TRANSPOSE to compute op(B) as the Hermitian transpose of B, otherwise, set to LA_NO_OPERATION to compute op(B) as B. |
m | The number of rows in the matrix C. |
n | The number of columns in the matrix C. |
k | The inner dimension of the matrix product A * op(B). |
alpha | A scalar multiplier. |
a | A P-element array containing the diagonal elements of matrix A where P = MIN(m , k ) if lside is true; else, P = MIN(n , k ) if lside is false. |
b | The LDB-by-TDB matrix B where (LDB = leading dimension of B, and TDB = trailing dimension of B):
|
ldb | The leading dimension of matrix B. |
beta | A scalar multiplier. |
c | The m by n matrix C. |
ldc | The leading dimension of matrix C. |
ldb
, or ldc
are not correct.int la_diag_mtx_mult_mixed | ( | bool | lside, |
int | opb, | ||
int | m, | ||
int | n, | ||
int | k, | ||
double complex | alpha, | ||
const double * | a, | ||
const double complex * | b, | ||
int | ldb, | ||
double complex | beta, | ||
double complex * | c, | ||
int | ldc ) |
Computes the matrix operation: \( C = \alpha A op(B) + \beta C \), or \( C = \alpha op(B) A + \beta C \).
lside | Set to true to apply matrix A from the left; else, set to false to apply matrix A from the left. |
opb | Set to LA_TRANSPOSE to compute op(B) as a direct transpose of B, set to LA_HERMITIAN_TRANSPOSE to compute op(B) as the Hermitian transpose of B, otherwise, set to LA_NO_OPERATION to compute op(B) as B. |
m | The number of rows in the matrix C. |
n | The number of columns in the matrix C. |
k | The inner dimension of the matrix product A * op(B). |
alpha | A scalar multiplier. |
a | A P-element array containing the diagonal elements of matrix A where P = MIN(m , k ) if lside is true; else, P = MIN(n , k ) if lside is false. |
b | The LDB-by-TDB matrix B where (LDB = leading dimension of B, and TDB = trailing dimension of B):
|
ldb | The leading dimension of matrix B. |
beta | A scalar multiplier. |
c | The m by n matrix C. |
ldc | The leading dimension of matrix C. |
ldb
, or ldc
are not correct.int la_eigen_asymm | ( | bool | vecs, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double complex * | vals, | ||
double complex * | v, | ||
int | ldv ) |
Computes the eigenvalues, and optionally the right eigenvectors of a square matrix.
vecs | Set to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues. |
n | The dimension of the matrix. |
a | On input, the N-by-N matrix on which to operate. On output, the contents of this matrix are overwritten. |
lda | The leading dimension of matrix A. |
vals | An N-element array containing the eigenvalues of the matrix. The eigenvalues are not sorted. |
v | An N-by-N matrix where the right eigenvectors will be written (one per column). |
lda
or ldv
is not correct.int la_eigen_cmplx | ( | bool | vecs, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
double complex * | vals, | ||
double complex * | v, | ||
int | ldv ) |
Computes the eigenvalues, and optionally the right eigenvectors of a square matrix.
vecs | Set to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues. |
n | The dimension of the matrix. |
a | On input, the N-by-N matrix on which to operate. On output, the contents of this matrix are overwritten. |
lda | The leading dimension of matrix A. |
vals | An N-element array containing the eigenvalues of the matrix. The eigenvalues are not sorted. |
v | An N-by-N matrix where the right eigenvectors will be written (one per column). |
lda
or ldv
is not correct.int la_eigen_gen | ( | bool | vecs, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double * | b, | ||
int | ldb, | ||
double complex * | alpha, | ||
double * | beta, | ||
double complex * | v, | ||
int | ldv ) |
Computes the eigenvalues, and optionally the right eigenvectors of a square matrix assuming the structure of the eigenvalue problem is \( A X = \lambda B X \).
vecs | Set to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues. |
n | The dimension of the matrix. |
a | On input, the N-by-N matrix A. On output, the contents of this matrix are overwritten. |
lda | The leading dimension of matrix A. |
b | On input, the N-by-N matrix B. On output, the contents of this matrix are overwritten. |
ldb | The leading dimension of matrix B. |
alpha | An N-element array a factor of the eigenvalues. The eigenvalues must be computed as ALPHA / BETA. This however, is not as trivial as it seems as it is entirely possible, and likely, that ALPHA / BETA can overflow or underflow. With that said, the values in ALPHA will always be less than and usually comparable with the NORM(A). |
beta | An N-element array that contains the denominator used to determine the eigenvalues as ALPHA / BETA. If used, the values in this array will always be less than and usually comparable with the NORM(B). |
v | An N-by-N matrix where the right eigenvectors will be written (one per column). |
lda
or ldv
is not correct.int la_eigen_symm | ( | bool | vecs, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double * | vals ) |
Computes the eigenvalues, and optionally the eigenvectors of a real, symmetric matrix.
vecs | Set to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues. |
n | The dimension of the matrix. |
a | On input, the N-by-N symmetric matrix on which to operate. On output, and if vecs is set to true, the matrix will contain the eigenvectors (one per column) corresponding to each eigenvalue in vals . If vecs is set to false, the lower triangular portion of the matrix is overwritten. |
lda | The leading dimension of matrix A. |
vals | An N-element array that will contain the eigenvalues sorted into ascending order. |
lda
is not correct.int la_form_lq | ( | int | m, |
int | n, | ||
double * | l, | ||
int | ldl, | ||
const double * | tau, | ||
double * | q, | ||
int | ldq ) |
Forms the matrix Q with orthonormal rows from the elementary reflectors returned by the base QR factorization algorithm.
m | The number of rows in R. |
n | The number of columns in R. |
l | On input, the M-by-N factored matrix as returned by the LQ factorization routine. On output, the lower triangular matrix L. |
ldl | The leading dimension of matrix L. |
tau | A MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r . |
q | An N-by-N matrix where the Q matrix will be written. |
ldq | The leading dimension of matrix Q. |
ldl
or ldq
are not correct.int la_form_lq_cmplx | ( | int | m, |
int | n, | ||
double complex * | l, | ||
int | ldl, | ||
const double complex * | tau, | ||
double complex * | q, | ||
int | ldq ) |
Forms the matrix Q with orthonormal rows from the elementary reflectors returned by the base QR factorization algorithm.
m | The number of rows in R. |
n | The number of columns in R. |
l | On input, the M-by-N factored matrix as returned by the LQ factorization routine. On output, the lower triangular matrix L. |
ldl | The leading dimension of matrix L. |
tau | A MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r . |
q | An N-by-N matrix where the Q matrix will be written. |
ldq | The leading dimension of matrix Q. |
ldl
or ldq
are not correct.int la_form_lu | ( | int | n, |
double * | a, | ||
int | lda, | ||
int * | ipvt, | ||
double * | u, | ||
int | ldu, | ||
double * | p, | ||
int | ldp ) |
Extracts the L, U, and P matrices from the LU factorization output from la_lu_factor.
n | The dimension of the input matrix. |
a | On input, the N-by-N matrix as output by la_lu_factor. On output, the N-by-N lower triangular matrix L. |
lda | The leading dimension of a . |
ipvt | The N-element pivot array as output by la_lu_factor. |
u | An N-by-N matrix where the U matrix will be written. |
ldu | The leading dimension of u . |
p | An N-by-N matrix where the row permutation matrix will be written. |
ldp | The leading dimension of p . |
lda
, ldu
, or ldp
is not correct. int la_form_lu_cmplx | ( | int | n, |
double complex * | a, | ||
int | lda, | ||
int * | ipvt, | ||
double complex * | u, | ||
int | ldu, | ||
double * | p, | ||
int | ldp ) |
Extracts the L, U, and P matrices from the LU factorization output from la_lu_factor_cmplx.
n | The dimension of the input matrix. |
a | On input, the N-by-N matrix as output by la_lu_factor_cmplx. On output, the N-by-N lower triangular matrix L. |
lda | The leading dimension of a . |
ipvt | The N-element pivot array as output by la_lu_factor_cmplx. |
u | An N-by-N matrix where the U matrix will be written. |
ldu | The leading dimension of u . |
p | An N-by-N matrix where the row permutation matrix will be written. |
ldp | The leading dimension of p . |
lda
, ldu
, or ldp
is not correct. int la_form_qr | ( | bool | fullq, |
int | m, | ||
int | n, | ||
double * | r, | ||
int | ldr, | ||
const double * | tau, | ||
double * | q, | ||
int | ldq ) |
Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm.
fullq | Set to true to always return the full Q matrix; else, set to false, and in the event that M > N, Q may be supplied as M-by-N, and therefore only return the useful submatrix Q1 (Q = [Q1, Q2]) as the factorization can be written as Q * R = [Q1, Q2] * [R1; 0]. |
m | The number of rows in R. |
n | The number of columns in R. |
r | On input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R. |
ldr | The leading dimension of matrix R. |
tau | A MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r . |
q | An M-by-M matrix where the full Q matrix will be written. In the event that fullq is set to false, and M > N, this matrix need only by M-by-N. |
ldq | The leading dimension of matrix Q. |
ldr
or ldq
are not correct.int la_form_qr_cmplx | ( | bool | fullq, |
int | m, | ||
int | n, | ||
double complex * | r, | ||
int | ldr, | ||
const double complex * | tau, | ||
double complex * | q, | ||
int | ldq ) |
Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm.
fullq | Set to true to always return the full Q matrix; else, set to false, and in the event that M > N, Q may be supplied as M-by-N, and therefore only return the useful submatrix Q1 (Q = [Q1, Q2]) as the factorization can be written as Q * R = [Q1, Q2] * [R1; 0]. |
m | The number of rows in R. |
n | The number of columns in R. |
r | On input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R. |
ldr | The leading dimension of matrix R. |
tau | A MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r . |
q | An M-by-M matrix where the full Q matrix will be written. In the event that fullq is set to false, and M > N, this matrix need only by M-by-N. |
ldq | The leading dimension of matrix Q. |
ldr
or ldq
are not correct.int la_form_qr_cmplx_pvt | ( | bool | fullq, |
int | m, | ||
int | n, | ||
double complex * | r, | ||
int | ldr, | ||
const double complex * | tau, | ||
const int * | pvt, | ||
double complex * | q, | ||
int | ldq, | ||
double complex * | p, | ||
int | ldp ) |
Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm. This routine also inflates the pivot array into an N-by-N matrix P such that \( A P = Q R \).
fullq | Set to true to always return the full Q matrix; else, set to false, and in the event that M > N, Q may be supplied as M-by-N, and therefore only return the useful submatrix Q1 (Q = [Q1, Q2]) as the factorization can be written as Q * R = [Q1, Q2] * [R1; 0]. |
m | The number of rows in R. |
n | The number of columns in R. |
r | On input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R. |
ldr | The leading dimension of matrix R. |
tau | A MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r . |
pvt | An N-element array containing the pivot information from the QR factorization. |
q | An M-by-M matrix where the full Q matrix will be written. In the event that fullq is set to false, and M > N, this matrix need only by M-by-N. |
ldq | The leading dimension of matrix Q. |
p | An N-by-N matrix where the pivot matrix P will be written. |
ldp | The leading dimension of matrix P. |
lda
is not correct.int la_form_qr_pvt | ( | bool | fullq, |
int | m, | ||
int | n, | ||
double * | r, | ||
int | ldr, | ||
const double * | tau, | ||
const int * | pvt, | ||
double * | q, | ||
int | ldq, | ||
double * | p, | ||
int | ldp ) |
Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm. This routine also inflates the pivot array into an N-by-N matrix P such that \( A P = Q R \).
fullq | Set to true to always return the full Q matrix; else, set to false, and in the event that M > N, Q may be supplied as M-by-N, and therefore only return the useful submatrix Q1 (Q = [Q1, Q2]) as the factorization can be written as Q * R = [Q1, Q2] * [R1; 0]. |
m | The number of rows in R. |
n | The number of columns in R. |
r | On input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R. |
ldr | The leading dimension of matrix R. |
tau | A MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r . |
pvt | An N-element array containing the pivot information from the QR factorization. |
q | An M-by-M matrix where the full Q matrix will be written. In the event that fullq is set to false, and M > N, this matrix need only by M-by-N. |
ldq | The leading dimension of matrix Q. |
p | An N-by-N matrix where the pivot matrix P will be written. |
ldp | The leading dimension of matrix P. |
lda
is not correct.int la_inverse | ( | int | n, |
double * | a, | ||
int | lda ) |
Computes the inverse of a square matrix.
n | The dimension of matrix A. |
a | On input, the N-by-N matrix to invert. On output, the inverted matrix. |
lda | The leading dimension of matrix A. |
lda
is not correct.int la_inverse_cmplx | ( | int | n, |
double complex * | a, | ||
int | lda ) |
Computes the inverse of a square matrix.
n | The dimension of matrix A. |
a | On input, the N-by-N matrix to invert. On output, the inverted matrix. |
lda | The leading dimension of matrix A. |
lda
is not correct.int la_lq_factor | ( | int | m, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double * | tau ) |
Computes the LQ factorization of an M-by-N matrix without pivoting.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau , represent the orthogonal matrix Q as a product of elementary reflectors. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors. |
lda
is not correct.int la_lq_factor_cmplx | ( | int | m, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
double complex * | tau ) |
Computes the LQ factorization of an M-by-N matrix without pivoting.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau , represent the orthogonal matrix Q as a product of elementary reflectors. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors. |
lda
is not correct.int la_lu_factor | ( | int | m, |
int | n, | ||
double * | a, | ||
int | lda, | ||
int * | ipvt ) |
Computes the LU factorization of an M-by-N matrix.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix on which to operate. On output, the LU factored matrix in the form [L\U] where the unit diagonal elements of L are not stored. |
lda | The leading dimension of matrix A. |
ipvt | An MIN(M, N)-element array used to track row-pivot operations. The array stored pivot information such that row I is interchanged with row IPVT(I). |
lda
is not correct.a
is found to be singular. int la_lu_factor_cmplx | ( | int | m, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
int * | ipvt ) |
Computes the LU factorization of an M-by-N matrix.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix on which to operate. On output, the LU factored matrix in the form [L\U] where the unit diagonal elements of L are not stored. |
lda | The leading dimension of matrix A. |
ipvt | An MIN(M, N)-element array used to track row-pivot operations. The array stored pivot information such that row I is interchanged with row IPVT(I). |
lda
is not correct.a
is found to be singular. int la_mtx_mult | ( | bool | transa, |
bool | transb, | ||
int | m, | ||
int | n, | ||
int | k, | ||
double | alpha, | ||
const double * | a, | ||
int | lda, | ||
const double * | b, | ||
int | ldb, | ||
double | beta, | ||
double * | c, | ||
int | ldc ) |
Computes the matrix operation \( C = \alpha op(A) op(B) + \beta C \).
transa | Set to true to compute op(A) as the transpose of A; else, set to false to compute op(A) as A. |
transb | Set to true to compute op(B) as the transpose of B; else, set to false to compute op(B) as B. |
m | The number of rows in c . |
n | The number of columns in c . |
k | The interior dimension of the product a and b . |
alpha | A scalar multiplier. |
a | If transa is true, this matrix must be k by m ; else, if transa is false, this matrix must be m by k . |
lda | The leading dimension of matrix a . |
b | If transb is true, this matrix must be n by k ; else, if transb is false, this matrix must be k by n . |
ldb | The leading dimension of matrix b . |
beta | A scalar multiplier. |
c | The m by n matrix C. |
ldc | The leading dimension of matrix c . |
lda
, ldb
, or ldc
are not correct. int la_mtx_mult_cmplx | ( | int | opa, |
int | opb, | ||
int | m, | ||
int | n, | ||
int | k, | ||
double complex | alpha, | ||
const double complex * | a, | ||
int | lda, | ||
const double complex * | b, | ||
int | ldb, | ||
double complex | beta, | ||
double complex * | c, | ||
int | ldc ) |
Computes the matrix operation \( C = \alpha op(A) op(B) + \beta C \).
opa | Set to LA_TRANSPOSE to compute op(A) as a direct transpose of A, set to LA_HERMITIAN_TRANSPOSE to compute op(A) as the Hermitian transpose of A, otherwise, set to LA_NO_OPERATION to compute op(A) as A. |
opb | Set to TLA_RANSPOSE to compute op(B) as a direct transpose of B, set to LA_HERMITIAN_TRANSPOSE to compute op(B) as the Hermitian transpose of B, otherwise, set to LA_NO_OPERATION to compute op(B) as B. |
mThe | number of rows in c . |
n | The number of columns in c . |
k | The interior dimension of the product a and b . |
alpha | A scalar multiplier. |
a | If opa is LA_TRANSPOSE or LA_HERMITIAN_TRANSPOSE, this matrix must be k by m ; else, this matrix must be m by k . |
lda | The leading dimension of matrix a . |
b | If opb is LA_TRANSPOSE or LA_HERMITIAN_TRANSPOSE, this matrix must be n by k ; else, this matrix must be k by n . |
ldb | The leading dimension of matrix b . |
beta | A scalar multiplier. |
c | The m by n matrix C. |
ldc | The leading dimension of matrix c . |
lda
, ldb
, or ldc
are not correct. int la_mult_lq | ( | bool | lside, |
bool | trans, | ||
int | m, | ||
int | n, | ||
int | k, | ||
const double * | a, | ||
int | lda, | ||
const double * | tau, | ||
double * | c, | ||
int | ldc ) |
Multiplies a general matrix by the orthogonal matrix Q from a LQ factorization such that: \( C = op(Q) C \), or \( C = C op(Q) \).
lside | Set to true to apply \( Q \) or \( Q^T \) from the left; else, set to false to apply \( Q \) or \( Q^T \) from the right. |
trans | Set to true to apply \( Q^T \); else, set to false. |
m | The number of rows in matrix C. |
n | The number of columns in matrix C. |
k | The number of elementary reflectors whose product defines the matrix Q. |
a | On input, an K-by-P matrix containing the elementary reflectors output from the LQ factorization. If lside is set to true, P = M; else, if lside is set to false, P = N. |
lda | The leading dimension of matrix A. |
tau | A K-element array containing the scalar factors of each elementary reflector defined in a . |
c | On input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C. |
ldc | The leading dimension of matrix C. |
lda
or ldc
are not correct.int la_mult_lq_cmplx | ( | bool | lside, |
bool | trans, | ||
int | m, | ||
int | n, | ||
int | k, | ||
const double complex * | a, | ||
int | lda, | ||
const double complex * | tau, | ||
double complex * | c, | ||
int | ldc ) |
Multiplies a general matrix by the orthogonal matrix Q from a LQ factorization such that: \( C = op(Q) C \), or \( C = C op(Q) \).
lside | Set to true to apply \( Q \) or \( Q^H \) from the left; else, set to false to apply \( Q \) or \( Q^H \) from the right. |
trans | Set to true to apply \( Q^H \); else, set to false. |
m | The number of rows in matrix C. |
n | The number of columns in matrix C. |
k | The number of elementary reflectors whose product defines the matrix Q. |
a | On input, an K-by-P matrix containing the elementary reflectors output from the LQ factorization. If lside is set to true, P = M; else, if lside is set to false, P = N. |
lda | The leading dimension of matrix A. |
tau | A K-element array containing the scalar factors of each elementary reflector defined in a . |
c | On input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C. |
ldc | The leading dimension of matrix C. |
lda
or ldc
are not correct.int la_mult_qr | ( | bool | lside, |
bool | trans, | ||
int | m, | ||
int | n, | ||
int | k, | ||
double * | a, | ||
int | lda, | ||
const double * | tau, | ||
double * | c, | ||
int | ldc ) |
Multiplies a general matrix by the orthogonal matrix Q from a QR factorization such that: \( C = op(Q) C \), or \( C = C op(Q) \).
lside | Set to true to apply \( Q \) or \( Q^T \) from the left; else, set to false to apply \( Q \) or \( Q^T \) from the right. |
trans | Set to true to apply \( Q^T \); else, set to false. |
m | The number of rows in matrix C. |
n | The number of columns in matrix C. |
k | The number of elementary reflectors whose product defines the matrix Q. |
a | On input, an LDA-by-K matrix containing the elementary reflectors output from the QR factorization. If lside is set to true, LDA = M, and M >= K >= 0; else, if lside is set to false, LDA = N, and N >= K >= 0. Notice, the contents of this matrix are restored on exit. |
lda | The leading dimension of matrix A. |
tau | A K-element array containing the scalar factors of each elementary reflector defined in a . |
c | On input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C. |
ldc | The leading dimension of matrix C. |
lda
or ldc
are not correct.int la_mult_qr_cmplx | ( | bool | lside, |
bool | trans, | ||
int | m, | ||
int | n, | ||
int | k, | ||
double complex * | a, | ||
int | lda, | ||
const double complex * | tau, | ||
double complex * | c, | ||
int | ldc ) |
Multiplies a general matrix by the orthogonal matrix Q from a QR factorization such that: \( C = op(Q) C \), or \( C = C op(Q) \).
lside | Set to true to apply \( Q \) or \( Q^H \) from the left; else, set to false to apply \( Q \) or \( Q^H \) from the right. |
trans | Set to true to apply \( Q^H \); else, set to false. |
m | The number of rows in matrix C. |
n | The number of columns in matrix C. |
k | The number of elementary reflectors whose product defines the matrix Q. |
a | On input, an LDA-by-K matrix containing the elementary reflectors output from the QR factorization. If lside is set to true, LDA = M, and M >= K >= 0; else, if lside is set to false, LDA = N, and N >= K >= 0. Notice, the contents of this matrix are restored on exit. |
lda | The leading dimension of matrix A. |
tau | A K-element array containing the scalar factors of each elementary reflector defined in a . |
c | On input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C. |
ldc | The leading dimension of matrix C. |
lda
or ldc
are not correct.int la_pinverse | ( | int | m, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double * | ainv, | ||
int | ldai ) |
Computes the Moore-Penrose pseudo-inverse of an M-by-N matrix by means of singular value decomposition.
m | The number of rows in the matrix. @parma n The number of columns in the matrix. |
a | On input, the M-by-N matrix to invert. The matrix is overwritten on output. |
lda | The leading dimension of matrix A. |
ainv | The N-by-M matrix where the pseudo-inverse of a will be written. |
ldai | The leading dimension of matrix AINV. |
lda
, or ldai
is not correct. int la_pinverse_cmplx | ( | int | m, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
double complex * | ainv, | ||
int | ldai ) |
Computes the Moore-Penrose pseudo-inverse of an M-by-N matrix by means of singular value decomposition.
m | The number of rows in the matrix. @parma n The number of columns in the matrix. |
a | On input, the M-by-N matrix to invert. The matrix is overwritten on output. |
lda | The leading dimension of matrix A. |
ainv | The N-by-M matrix where the pseudo-inverse of a will be written. |
ldai | The leading dimension of matrix AINV. |
lda
, or ldai
is not correct. int la_qr_factor | ( | int | m, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double * | tau ) |
Computes the QR factorization of an M-by-N matrix without pivoting.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau , represent the orthogonal matrix Q as a product of elementary reflectors. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors. |
lda
is not correct.int la_qr_factor_cmplx | ( | int | m, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
double complex * | tau ) |
Computes the QR factorization of an M-by-N matrix without pivoting.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau , represent the orthogonal matrix Q as a product of elementary reflectors. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors. |
lda
is not correct.int la_qr_factor_cmplx_pvt | ( | int | m, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
double complex * | tau, | ||
int * | jpvt ) |
Computes the QR factorization of an M-by-N matrix with column pivoting.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau , represent the orthogonal matrix Q as a product of elementary reflectors. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors. |
jpvt | On input, an N-element array that if JPVT(I) .ne. 0, the I-th column of A is permuted to the front of A * P; if JPVT(I) = 0, the I-th column of A is a free column. On output, if JPVT(I) = K, then the I-th column of A * P was the K-th column of A. |
lda
is not correct.int la_qr_factor_pvt | ( | int | m, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double * | tau, | ||
int * | jpvt ) |
Computes the QR factorization of an M-by-N matrix with column pivoting.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau , represent the orthogonal matrix Q as a product of elementary reflectors. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors. |
jpvt | On input, an N-element array that if JPVT(I) .ne. 0, the I-th column of A is permuted to the front of A * P; if JPVT(I) = 0, the I-th column of A is a free column. On output, if JPVT(I) = K, then the I-th column of A * P was the K-th column of A. |
lda
is not correct.int la_qr_rank1_update | ( | int | m, |
int | n, | ||
double * | q, | ||
int | ldq, | ||
double * | r, | ||
int | ldr, | ||
double * | u, | ||
double * | v ) |
Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where \( A = Q R \), and \( A1 = A + U V^T \) such that \( A1 = Q1 R1 \).
m | The number of rows in R. |
n | The number of columns in R. |
q | On input, the original M-by-K orthogonal matrix Q. On output, the updated matrix Q1. |
ldq | The leading dimension of matrix Q. |
r | On input, the M-by-N matrix R. On output, the updated matrix R1. |
ldr | The leading dimension of matrix R. |
u | On input, the M-element U update vector. On output, the original content of the array is overwritten. |
v | On input, the N-element V update vector. On output, the original content of the array is overwritten. |
ldq
or ldr
is not correct.int la_qr_rank1_update_cmplx | ( | int | m, |
int | n, | ||
double complex * | q, | ||
int | ldq, | ||
double complex * | r, | ||
int | ldr, | ||
double complex * | u, | ||
double complex * | v ) |
Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where \( A = Q R \), and \( A1 = A + U V^H \) such that \( A1 = Q1 R1 \).
m | The number of rows in R. |
n | The number of columns in R. |
q | On input, the original M-by-K orthogonal matrix Q. On output, the updated matrix Q1. |
ldq | The leading dimension of matrix Q. |
r | On input, the M-by-N matrix R. On output, the updated matrix R1. |
ldr | The leading dimension of matrix R. |
u | On input, the M-element U update vector. On output, the original content of the array is overwritten. |
v | On input, the N-element V update vector. On output, the original content of the array is overwritten. |
ldq
or ldr
is not correct.int la_rank | ( | int | m, |
int | n, | ||
double * | a, | ||
int | lda, | ||
int * | rnk ) |
Computes the rank of a matrix.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | The M-by-N matrix. The matrix is overwritten as part of this operation. |
lda | The leading dimension of matrix A. |
rnk | The rank of a . |
lda
is not correct.int la_rank1_update | ( | int | m, |
int | n, | ||
double | alpha, | ||
const double * | x, | ||
const double * | y, | ||
double * | a, | ||
int | lda ) |
Performs the rank-1 update to matrix A such that: \( A = \alpha X Y^T + A \), where \( A \) is an M-by-N matrix, \( \alpha \) is a scalar, \( X \) is an M-element array, and \( Y \) is an N-element array.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
alpha | The scalar multiplier. |
x | An M-element array. |
y | An N-element array. |
a | On input, the M-by-N matrix to update. On output, the updated M-by-N matrix. |
lda | The leading dimension of matrix A. |
lda
is not correct. int la_rank1_update_cmplx | ( | int | m, |
int | n, | ||
double complex | alpha, | ||
const double complex * | x, | ||
const double complex * | y, | ||
double complex * | a, | ||
int | lda ) |
Performs the rank-1 update to matrix A such that: \( A = \alpha X Y^H + A \), where \( A \) is an M-by-N matrix, \( \alpha \) is a scalar, \( X \) is an M-element array, and \( Y \) is an N-element array.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
alpha | The scalar multiplier. |
x | An M-element array. |
y | An N-element array. |
a | On input, the M-by-N matrix to update. On output, the updated M-by-N matrix. |
lda | The leading dimension of matrix A. |
lda
is not correct. int la_rank_cmplx | ( | int | m, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
int * | rnk ) |
Computes the rank of a matrix.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | The M-by-N matrix. The matrix is overwritten as part of this operation. |
lda | The leading dimension of matrix A. |
rnk | The rank of a . |
lda
is not correct.int la_solve_cholesky | ( | bool | upper, |
int | m, | ||
int | n, | ||
const double * | a, | ||
int | lda, | ||
double * | b, | ||
int | ldb ) |
Solves a system of Cholesky factored equations.
upper | Set to true if the original matrix A was factored such that \( A = U^T U \); else, set to false if the factorization of A was \( A = L^T L \). |
m | The number of rows in matrix B. |
n | The number of columns in matrix B. |
a | The M-by-M Cholesky factored matrix. |
lda | The leading dimension of matrix A. |
b | On input, the M-by-N right-hand-side matrix B. On output, the M-by-N solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct. int la_solve_cholesky_cmplx | ( | bool | upper, |
int | m, | ||
int | n, | ||
const double complex * | a, | ||
int | lda, | ||
double complex * | b, | ||
int | ldb ) |
Solves a system of Cholesky factored equations.
upper | Set to true if the original matrix A was factored such that \( A = U^T U \); else, set to false if the factorization of A was \( A = L^T L \). |
m | The number of rows in matrix B. |
n | The number of columns in matrix B. |
a | The M-by-M Cholesky factored matrix. |
lda | The leading dimension of matrix A. |
b | On input, the M-by-N right-hand-side matrix B. On output, the M-by-N solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct. int la_solve_least_squares | ( | int | m, |
int | n, | ||
int | k, | ||
double * | a, | ||
int | lda, | ||
double * | b, | ||
int | ldb ) |
Solves the overdetermined or underdetermined system ( \( A X = B \)) of M equations of N unknowns using a QR or LQ factorization of the matrix A. Notice, it is assumed that matrix A has full rank.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N matrix A. On output, if M >= N, the QR factorization of A in the form as output by qr_factor; else, if M < N, the LQ factorization of A. |
lda | The leading dimension of matrix A. |
b | If M >= N, the M-by-NRHS matrix B. On output, the first N rows contain the N-by-NRHS solution matrix X. If M < N, an N-by-NRHS matrix with the first M rows containing the matrix B. On output, the N-by-NRHS solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct.a
is not of full rank. int la_solve_least_squares_cmplx | ( | int | m, |
int | n, | ||
int | k, | ||
double complex * | a, | ||
int | lda, | ||
double complex * | b, | ||
int | ldb ) |
Solves the overdetermined or underdetermined system ( \( A X = B \)) of M equations of N unknowns using a QR or LQ factorization of the matrix A. Notice, it is assumed that matrix A has full rank.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N matrix A. On output, if M >= N, the QR factorization of A in the form as output by qr_factor; else, if M < N, the LQ factorization of A. |
lda | The leading dimension of matrix A. |
b | If M >= N, the M-by-NRHS matrix B. On output, the first N rows contain the N-by-NRHS solution matrix X. If M < N, an N-by-NRHS matrix with the first M rows containing the matrix B. On output, the N-by-NRHS solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct.a
is not of full rank. int la_solve_lq | ( | int | m, |
int | n, | ||
int | k, | ||
const double * | a, | ||
int | lda, | ||
const double * | tau, | ||
double * | b, | ||
int | ldb ) |
Solves a system of M QR-factored equations of N unknowns where N >= M.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N QR factored matrix as returned by lq_factor. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by lq_factor. |
b | On input, an N-by-K matrix containing the first M rows of the right-hand-side matrix. On output, the solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct, or if m
is less than n
.int la_solve_lq_cmplx | ( | int | m, |
int | n, | ||
int | k, | ||
const double complex * | a, | ||
int | lda, | ||
const double complex * | tau, | ||
double complex * | b, | ||
int | ldb ) |
Solves a system of M QR-factored equations of N unknowns where N >= M.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N QR factored matrix as returned by lq_factor. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by lq_factor. |
b | On input, an N-by-K matrix containing the first M rows of the right-hand-side matrix. On output, the solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct, or if m
is less than n
.int la_solve_lu | ( | int | m, |
int | n, | ||
const double * | a, | ||
int | lda, | ||
const int * | ipvt, | ||
double * | b, | ||
int | ldb ) |
Solves a system of LU-factored equations.
m | The number of rows in matrix B. |
n | The number of columns in matrix B. |
a | The M-by-M LU factored matrix. |
lda | The leading dimension of matrix A. |
ipvt | The M-element pivot array from the LU factorization. |
b | On input, the M-by-N right-hand-side. On output, the M-by-N solution. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct. int la_solve_lu_cmplx | ( | int | m, |
int | n, | ||
const double complex * | a, | ||
int | lda, | ||
const int * | ipvt, | ||
double complex * | b, | ||
int | ldb ) |
Solves a system of LU-factored equations.
m | The number of rows in matrix B. |
n | The number of columns in matrix B. |
a | The M-by-M LU factored matrix. |
lda | The leading dimension of matrix A. |
ipvt | The M-element pivot array from the LU factorization. |
b | On input, the M-by-N right-hand-side. On output, the M-by-N solution. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct. int la_solve_qr | ( | int | m, |
int | n, | ||
int | k, | ||
double * | a, | ||
int | lda, | ||
const double * | tau, | ||
double * | b, | ||
int | ldb ) |
Solves a system of M QR-factored equations of N unknowns where M >= N.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor. |
b | On input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct, or if m
is less than n
.int la_solve_qr_cmplx | ( | int | m, |
int | n, | ||
int | k, | ||
double complex * | a, | ||
int | lda, | ||
const double complex * | tau, | ||
double complex * | b, | ||
int | ldb ) |
Solves a system of M QR-factored equations of N unknowns where M >= N.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor. |
b | On input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct, or if m
is less than n
.int la_solve_qr_cmplx_pvt | ( | int | m, |
int | n, | ||
int | k, | ||
double complex * | a, | ||
int | lda, | ||
const double complex * | tau, | ||
const int * | jpvt, | ||
double complex * | b, | ||
int | ldb ) |
Solves a system of M QR-factored equations of N unknowns.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor. |
jpvt | The N-element array that was used to track the column pivoting operations in the QR factorization. |
b | On input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct.int la_solve_qr_pvt | ( | int | m, |
int | n, | ||
int | k, | ||
double * | a, | ||
int | lda, | ||
const double * | tau, | ||
const int * | jpvt, | ||
double * | b, | ||
int | ldb ) |
Solves a system of M QR-factored equations of N unknowns.
m | The number of equations (rows in matrix A). |
n | The number of unknowns (columns in matrix A). |
k | The number of columns in the right-hand-side matrix. |
a | On input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored. |
lda | The leading dimension of matrix A. |
tau | A MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor. |
jpvt | The N-element array that was used to track the column pivoting operations in the QR factorization. |
b | On input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct.int la_solve_tri_mtx | ( | bool | lside, |
bool | upper, | ||
bool | trans, | ||
bool | nounit, | ||
int | m, | ||
int | n, | ||
double | alpha, | ||
const double * | a, | ||
int | lda, | ||
double * | b, | ||
int | ldb ) |
Solves one of the matrix equations: \( op(A) X = \alpha B \), or \( X op(A) = \alpha B \), where \( A \) is a triangular matrix.
lside | Set to true to solve \( op(A) X = \alpha B \); else, set to false to solve \( X op(A) = \alpha B \). |
upper | Set to true if A is an upper triangular matrix; else, set to false if A is a lower triangular matrix. |
trans | Set to true if \( op(A) = A^T \); else, set to false if \( op(A) = A \). |
nounit | Set to true if A is not a unit-diagonal matrix (ones on every diagonal element); else, set to false if A is a unit-diagonal matrix. |
m | The number of rows in matrix B. |
n | The number of columns in matrix B. |
alpha | The scalar multiplier to B. |
a | If lside is true, the M-by-M triangular matrix on which to operate; else, if lside is false, the N-by-N triangular matrix on which to operate. |
lda | The leading dimension of matrix A. |
b | On input, the M-by-N right-hand-side. On output, the M-by-N solution. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct. int la_solve_tri_mtx_cmplx | ( | bool | lside, |
bool | upper, | ||
bool | trans, | ||
bool | nounit, | ||
int | m, | ||
int | n, | ||
double complex | alpha, | ||
const double complex * | a, | ||
int | lda, | ||
double complex * | b, | ||
int | ldb ) |
Solves one of the matrix equations: \( op(A) X = \alpha B \), or \( X op(A) = \alpha B \), where \( A \) is a triangular matrix.
lside | Set to true to solve \( op(A) X = \alpha B \); else, set to false to solve \( X op(A) = \alpha B \). |
upper | Set to true if A is an upper triangular matrix; else, set to false if A is a lower triangular matrix. |
trans | Set to true if \( op(A) = A^H \); else, set to false if \( op(A) = A \). |
nounit | Set to true if A is not a unit-diagonal matrix (ones on every diagonal element); else, set to false if A is a unit-diagonal matrix. |
m | The number of rows in matrix B. |
n | The number of columns in matrix B. |
alpha | The scalar multiplier to B. |
a | If lside is true, the M-by-M triangular matrix on which to operate; else, if lside is false, the N-by-N triangular matrix on which to operate. |
lda | The leading dimension of matrix A. |
b | On input, the M-by-N right-hand-side. On output, the M-by-N solution. |
ldb | The leading dimension of matrix B. |
lda
, or ldb
is not correct. int la_sort_eigen | ( | bool | ascend, |
int | n, | ||
double * | vals, | ||
double * | vecs, | ||
int | ldv ) |
A sorting routine specifically tailored for sorting of eigenvalues and their associated eigenvectors using a quick-sort approach.
ascend | |
n | The number of eigenvalues. |
vals | On input, an N-element array containing the eigenvalues. On output, the sorted eigenvalues. |
vecs | On input, an N-by-N matrix containing the eigenvectors associated with vals (one vector per column). On output, the sorted eigenvector matrix. |
ldv | The leading dimension of vecs . |
ldv
is not correct.int la_sort_eigen_cmplx | ( | bool | ascend, |
int | n, | ||
double complex * | vals, | ||
double complex * | vecs, | ||
int | ldv ) |
A sorting routine specifically tailored for sorting of eigenvalues and their associated eigenvectors using a quick-sort approach.
ascend | |
n | The number of eigenvalues. |
vals | On input, an N-element array containing the eigenvalues. On output, the sorted eigenvalues. |
vecs | On input, an N-by-N matrix containing the eigenvectors associated with vals (one vector per column). On output, the sorted eigenvector matrix. |
ldv | The leading dimension of vecs . |
ldv
is not correct.int la_svd | ( | int | m, |
int | n, | ||
double * | a, | ||
int | lda, | ||
double * | s, | ||
double * | u, | ||
int | ldu, | ||
double * | vt, | ||
int | ldv ) |
Computes the singular value decomposition of a matrix \( A \). The SVD is defined as: \( A = U S V^T \), where \( U \) is an M-by-M orthogonal matrix, \( S \) is an M-by-N diagonal matrix, and \( V \) is an N-by-N orthogonal matrix.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. The matrix is overwritten on output. |
lda | The leading dimension of matrix A. |
s | A MIN(M, N)-element array containing the singular values of a sorted in descending order. |
u | An M-by-M matrix where the orthogonal U matrix will be written. |
ldu | The leading dimension of matrix U. |
vt | An N-by-N matrix where the transpose of the right singular vector matrix V. |
ldv | The leading dimension of matrix V. |
lda
, ldu
, or ldv
is not correct.int la_svd_cmplx | ( | int | m, |
int | n, | ||
double complex * | a, | ||
int | lda, | ||
double * | s, | ||
double complex * | u, | ||
int | ldu, | ||
double complex * | vt, | ||
int | ldv ) |
Computes the singular value decomposition of a matrix \( A \). The SVD is defined as: \( A = U S V^H \), where \( U \) is an M-by-M orthogonal matrix, \( S \) is an M-by-N diagonal matrix, and \( V \) is an N-by-N orthogonal matrix.
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | On input, the M-by-N matrix to factor. The matrix is overwritten on output. |
lda | The leading dimension of matrix A. |
s | A MIN(M, N)-element array containing the singular values of a sorted in descending order. |
u | An M-by-M matrix where the orthogonal U matrix will be written. |
ldu | The leading dimension of matrix U. |
vt | An N-by-N matrix where the conjugate transpose of the right singular vector matrix V. |
ldv | The leading dimension of vt . |
lda
, ldu
, or ldv
is not correct.int la_trace | ( | int | m, |
int | n, | ||
const double * | a, | ||
int | lda, | ||
double * | rst ) |
Computes the trace of a matrix (the sum of the main diagonal elements).
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | The M-by-N matrix on which to operate. |
lda | The leading dimension of the matrix. |
rst | The results of the operation. |
lda
is not correct. int la_trace_cmplx | ( | int | m, |
int | n, | ||
const double complex * | a, | ||
int | lda, | ||
double complex * | rst ) |
Computes the trace of a matrix (the sum of the main diagonal elements).
m | The number of rows in the matrix. |
n | The number of columns in the matrix. |
a | The M-by-N matrix on which to operate. |
lda | The leading dimension of the matrix. |
rst | The results of the operation. |
lda
is not correct. int la_tri_mtx_mult | ( | bool | upper, |
double | alpha, | ||
int | n, | ||
const double * | a, | ||
int | lda, | ||
double | beta, | ||
double * | b, | ||
int | ldb ) |
Computes the triangular matrix operation: \( B = \alpha A^T A + \beta B \), or \( B = \alpha A A^T + \beta B \), where \( A \) is a triangular matrix.
upper | Set to true if matrix \( A \) is upper triangular, and \( B = \alpha A^T A + \beta B \) is to be calculated; else, set to false if \( A \) is lower triangular, and \( B = \alpha A A^T + \beta B \) is to be computed. |
alpha | A scalar multiplier. |
n | The dimension of the matrix. |
a | The n by n triangular matrix A. Notice, if upper is true, only the upper triangular portion of this matrix is referenced; else, if upper is false, only the lower triangular portion of this matrix is referenced. |
lda | The leading dimension of matrix A. |
beta | A scalar multiplier. |
b | On input, the n by n matrix B. On output, the n by n resulting matrix. |
ldb | The leading dimension of matrix B. |
lda
or ldb
are not correct. int la_tri_mtx_mult_cmplx | ( | bool | upper, |
double complex | alpha, | ||
int | n, | ||
const double complex * | a, | ||
int | lda, | ||
double complex | beta, | ||
double complex * | b, | ||
int | ldb ) |
Computes the triangular matrix operation: \( B = \alpha A^T A + \beta B \), or \( B = \alpha A A^T + \beta B \), where \( A \) is a triangular matrix.
upper | Set to true if matrix \( A \) is upper triangular, and \( B = \alpha A^T A + \beta B \) is to be calculated; else, set to false if \( A \) is lower triangular, and \( B = \alpha A A^T + \beta B \) is to be computed. |
alpha | A scalar multiplier. |
n | The dimension of the matrix. |
a | The n by n triangular matrix A. Notice, if upper is true, only the upper triangular portion of this matrix is referenced; else, if upper is false, only the lower triangular portion of this matrix is referenced. |
lda | The leading dimension of matrix A. |
beta | A scalar multiplier. |
b | On input, the n by n matrix B. On output, the n by n resulting matrix. |
ldb | The leading dimension of matrix B. |
lda
or ldb
are not correct.