linalg 1.8.2
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
linalg.h File Reference
#include <stdbool.h>
#include <complex.h>
Include dependency graph for linalg.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)
 

Macro Definition Documentation

◆ LA_ARRAY_SIZE_ERROR

#define LA_ARRAY_SIZE_ERROR   102

Definition at line 196 of file linalg.h.

◆ LA_CONVERGENCE_ERROR

#define LA_CONVERGENCE_ERROR   106

Definition at line 200 of file linalg.h.

◆ LA_HERMITIAN_TRANSPOSE

#define LA_HERMITIAN_TRANSPOSE   2

Definition at line 193 of file linalg.h.

◆ LA_INVALID_INPUT_ERROR

#define LA_INVALID_INPUT_ERROR   101

Definition at line 195 of file linalg.h.

◆ LA_INVALID_OPERATION_ERROR

#define LA_INVALID_OPERATION_ERROR   107

Definition at line 201 of file linalg.h.

◆ LA_MATRIX_FORMAT_ERROR

#define LA_MATRIX_FORMAT_ERROR   104

Definition at line 198 of file linalg.h.

◆ LA_NO_ERROR

#define LA_NO_ERROR   0

Definition at line 194 of file linalg.h.

◆ LA_NO_OPERATION

#define LA_NO_OPERATION   0

Definition at line 191 of file linalg.h.

◆ LA_OUT_OF_MEMORY_ERROR

#define LA_OUT_OF_MEMORY_ERROR   105

Definition at line 199 of file linalg.h.

◆ LA_SINGULAR_MATRIX_ERROR

#define LA_SINGULAR_MATRIX_ERROR   103

Definition at line 197 of file linalg.h.

◆ LA_TRANSPOSE

#define LA_TRANSPOSE   1

Definition at line 192 of file linalg.h.

Function Documentation

◆ la_band_diag_mtx_mult()

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.

Parameters
leftSet to true to compute A = alpha * A * B; else, set to false to compute A = alpha * B * A.
mThe number of rows in matrix A.
klThe number of subdiagonals. Must be at least 0.
kuThe number of superdiagonals. Must be at least 0.
alphaThe scalar multiplier.
aThe 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
bAn array containing the diagonal elements of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if b and f are not compatible in size.
  • LA_ARRAY_SIZE_ERROR: Occurs if a and b are not compatible in terms of internal dimensions.
  • LA_INVALID_INPUT_ERROR: Occurs if either ku or kl are not zero or greater.

◆ la_band_diag_mtx_mult_cmplx()

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.

Parameters
leftSet to true to compute A = alpha * A * B; else, set to false to compute A = alpha * B * A.
mThe number of rows in matrix A.
klThe number of subdiagonals. Must be at least 0.
kuThe number of superdiagonals. Must be at least 0.
alphaThe scalar multiplier.
aThe 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
bAn array containing the diagonal elements of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if b and f are not compatible in size.
  • LA_ARRAY_SIZE_ERROR: Occurs if a and b are not compatible in terms of internal dimensions.
  • LA_INVALID_INPUT_ERROR: Occurs if either ku or kl are not zero or greater.

◆ la_band_mtx_vec_mult()

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.

Parameters
transSet to true for op(A) == A^T; else, false for op(A) == A.
klThe number of subdiagonals. Must be at least 0.
kuThe number of superdiagonals. Must be at least 0.
alphaA scalar multiplier.
aThe M-by-N matrix A storing the banded matrix in a compressed form supplied column by column.
xIf trans is true, this is an M-element vector; else, if trans is false, this is an N-element vector.
betaA scalar multiplier.
yOn 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.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input arrays are not sized appropriately.
  • LA_INVALID_INPUT_ERROR: Occurs if either ku or kl are not zero or greater.

◆ la_band_mtx_vec_mult_cmplx()

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.

Parameters
transSet 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.
klThe number of subdiagonals. Must be at least 0.
kuThe number of superdiagonals. Must be at least 0.
alphaA scalar multiplier.
aThe M-by-N matrix A storing the banded matrix in a compressed form supplied column by column.
xIf trans is true, this is an M-element vector; else, if trans is false, this is an N-element vector.
betaA scalar multiplier.
yOn 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.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input arrays are not sized appropriately.
  • LA_INVALID_INPUT_ERROR: Occurs if either ku or kl are not zero or greater.

◆ la_band_to_full_mtx()

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.

Parameters
klThe number of subdiagonals. Must be at least 0.
kuThe number of superdiagonals. Must be at least 0.
bThe banded matrix to convert, stored in dense form. See band_mtx_vec_mult for details on this storage method.
fThe M-by-N element full matrix.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if b and f are not compatible in size.
  • LA_INVALID_INPUT_ERROR: Occurs if either ku or kl are not zero or greater.

◆ la_band_to_full_mtx_cmplx()

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.

Parameters
klThe number of subdiagonals. Must be at least 0.
kuThe number of superdiagonals. Must be at least 0.
bThe banded matrix to convert, stored in dense form. See band_mtx_vec_mult for details on this storage method.
fThe M-by-N element full matrix.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if b and f are not compatible in size.
  • LA_INVALID_INPUT_ERROR: Occurs if either ku or kl are not zero or greater.

◆ la_cholesky_factor()

int la_cholesky_factor ( bool upper,
int n,
double * a,
int lda )

Computes the Cholesky factorization of a symmetric, positive definite matrix.

Parameters
upperSet 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 \).
nThe dimension of matrix A.
aOn 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.
ldaThe leading dimension of matrix A.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_MATRIX_FORMAT_ERROR: Occurs if a is not positive definite.

◆ la_cholesky_factor_cmplx()

int la_cholesky_factor_cmplx ( bool upper,
int n,
double complex * a,
int lda )

Computes the Cholesky factorization of a symmetric, positive definite matrix.

Parameters
upperSet 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 \).
nThe dimension of matrix A.
aOn 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.
ldaThe leading dimension of matrix A.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_MATRIX_FORMAT_ERROR: Occurs if a is not positive definite.

◆ la_cholesky_rank1_downdate()

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

Parameters
nThe dimension of the matrix.
rOn input, the N-by-N upper triangular matrix R. On output, the updated matrix R1.
ldrThe leading dimension of matrix R.
uOn input, the N-element update vector U. On output, the rotation sines used to transform R to R1.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldr is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_MATRIX_FORMAT_ERROR: Occurs if the downdated matrix is not positive definite.

◆ la_cholesky_rank1_downdate_cmplx()

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

Parameters
nThe dimension of the matrix.
rOn input, the N-by-N upper triangular matrix R. On output, the updated matrix R1.
ldrThe leading dimension of matrix R.
uOn input, the N-element update vector U. On output, the rotation sines used to transform R to R1.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldr is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_MATRIX_FORMAT_ERROR: Occurs if the downdated matrix is not positive definite.

◆ la_cholesky_rank1_update()

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

Parameters
nThe dimension of the matrix.
rOn input, the N-by-N upper triangular matrix R. On output, the updated matrix R1.
ldrThe leading dimension of matrix R.
uOn input, the N-element update vector U. On output, the rotation sines used to transform R to R1.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldr is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_cholesky_rank1_update_cmplx()

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

WARNING
Untested code that needs more work before being ready for normal use.
Parameters
nThe dimension of the matrix.
rOn input, the N-by-N upper triangular matrix R. On output, the updated matrix R1.
ldrThe leading dimension of matrix R.
uOn input, the N-element update vector U. On output, the rotation sines used to transform R to R1.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldr is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_det()

int la_det ( int n,
double * a,
int lda,
double * d )

Computes the determinant of a square matrix.

Parameters
nThe dimension of the matrix.
aThe N-by-N matrix. The matrix is overwritten on output.
ldaThe leading dimension of the matrix.
dThe determinant of a.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input arrays are not sized appropriately.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_det_cmplx()

int la_det_cmplx ( int n,
double complex * a,
int lda,
double complex * d )

Computes the determinant of a square matrix.

Parameters
nThe dimension of the matrix.
aThe N-by-N matrix. The matrix is overwritten on output.
ldaThe leading dimension of the matrix.
dThe determinant of a.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input arrays are not sized appropriately.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_diag_mtx_mult()

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

Parameters
lsideSet to true to apply matrix A from the left; else, set to false to apply matrix A from the left.
transSet to true if op(B) == B**T; else, set to false if op(B) == B.
mThe number of rows in the matrix C.
nThe number of columns in the matrix C.
kThe inner dimension of the matrix product A * op(B).
alphaA scalar multiplier.
aA 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.
bThe LDB-by-TDB matrix B where (LDB = leading dimension of B, and TDB = trailing dimension of B):
  • lside == true & trans == true: LDB = n, TDB = k
  • lside == true & trans == false: LDB = k, TDB = n
  • lside == false & trans == true: LDB = k, TDB = m
  • lside == false & trans == false: LDB = m, TDB = k
ldbThe leading dimension of matrix B.
betaA scalar multiplier.
cThe m by n matrix C.
ldcThe leading dimension of matrix C.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldb, or ldc are not correct.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input array sizes are incorrect.

◆ la_diag_mtx_mult_cmplx()

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

Parameters
lsideSet to true to apply matrix A from the left; else, set to false to apply matrix A from the left.
opbSet 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.
mThe number of rows in the matrix C.
nThe number of columns in the matrix C.
kThe inner dimension of the matrix product A * op(B).
alphaA scalar multiplier.
aA 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.
bThe LDB-by-TDB matrix B where (LDB = leading dimension of B, and TDB = trailing dimension of B):
  • lside == true & trans == true: LDB = n, TDB = k
  • lside == true & trans == false: LDB = k, TDB = n
  • lside == false & trans == true: LDB = k, TDB = m
  • lside == false & trans == false: LDB = m, TDB = k
ldbThe leading dimension of matrix B.
betaA scalar multiplier.
cThe m by n matrix C.
ldcThe leading dimension of matrix C.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldb, or ldc are not correct.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input array sizes are incorrect.

◆ la_diag_mtx_mult_mixed()

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

Parameters
lsideSet to true to apply matrix A from the left; else, set to false to apply matrix A from the left.
opbSet 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.
mThe number of rows in the matrix C.
nThe number of columns in the matrix C.
kThe inner dimension of the matrix product A * op(B).
alphaA scalar multiplier.
aA 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.
bThe LDB-by-TDB matrix B where (LDB = leading dimension of B, and TDB = trailing dimension of B):
  • lside == true & trans == true: LDB = n, TDB = k
  • lside == true & trans == false: LDB = k, TDB = n
  • lside == false & trans == true: LDB = k, TDB = m
  • lside == false & trans == false: LDB = m, TDB = k
ldbThe leading dimension of matrix B.
betaA scalar multiplier.
cThe m by n matrix C.
ldcThe leading dimension of matrix C.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldb, or ldc are not correct.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input array sizes are incorrect.

◆ la_eigen_asymm()

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.

Parameters
vecsSet to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues.
nThe dimension of the matrix.
aOn input, the N-by-N matrix on which to operate. On output, the contents of this matrix are overwritten.
ldaThe leading dimension of matrix A.
valsAn N-element array containing the eigenvalues of the matrix. The eigenvalues are not sorted.
vAn N-by-N matrix where the right eigenvectors will be written (one per column).
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldv is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs if the algorithm failed to converge.

◆ la_eigen_cmplx()

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.

Parameters
vecsSet to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues.
nThe dimension of the matrix.
aOn input, the N-by-N matrix on which to operate. On output, the contents of this matrix are overwritten.
ldaThe leading dimension of matrix A.
valsAn N-element array containing the eigenvalues of the matrix. The eigenvalues are not sorted.
vAn N-by-N matrix where the right eigenvectors will be written (one per column).
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldv is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs if the algorithm failed to converge.

◆ la_eigen_gen()

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

Parameters
vecsSet to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues.
nThe dimension of the matrix.
aOn input, the N-by-N matrix A. On output, the contents of this matrix are overwritten.
ldaThe leading dimension of matrix A.
bOn input, the N-by-N matrix B. On output, the contents of this matrix are overwritten.
ldbThe leading dimension of matrix B.
alphaAn 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).
betaAn 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).
vAn N-by-N matrix where the right eigenvectors will be written (one per column).
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldv is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs if the algorithm failed to converge.

◆ la_eigen_symm()

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.

Parameters
vecsSet to true to compute the eigenvectors as well as the eigenvalues; else, set to false to just compute the eigenvalues.
nThe dimension of the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
valsAn N-element array that will contain the eigenvalues sorted into ascending order.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs if the algorithm failed to converge.

◆ la_form_lq()

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.

Parameters
mThe number of rows in R.
nThe number of columns in R.
lOn input, the M-by-N factored matrix as returned by the LQ factorization routine. On output, the lower triangular matrix L.
ldlThe leading dimension of matrix L.
tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
qAn N-by-N matrix where the Q matrix will be written.
ldqThe leading dimension of matrix Q.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldl or ldq are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_form_lq_cmplx()

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.

WARNING
Untested code that needs more work before being ready for normal use.
Parameters
mThe number of rows in R.
nThe number of columns in R.
lOn input, the M-by-N factored matrix as returned by the LQ factorization routine. On output, the lower triangular matrix L.
ldlThe leading dimension of matrix L.
tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
qAn N-by-N matrix where the Q matrix will be written.
ldqThe leading dimension of matrix Q.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldl or ldq are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_form_lu()

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.

Parameters
nThe dimension of the input matrix.
aOn input, the N-by-N matrix as output by la_lu_factor. On output, the N-by-N lower triangular matrix L.
ldaThe leading dimension of a.
ipvtThe N-element pivot array as output by la_lu_factor.
uAn N-by-N matrix where the U matrix will be written.
lduThe leading dimension of u.
pAn N-by-N matrix where the row permutation matrix will be written.
ldpThe leading dimension of p.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, ldu, or ldp is not correct.

◆ la_form_lu_cmplx()

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.

Parameters
nThe dimension of the input matrix.
aOn input, the N-by-N matrix as output by la_lu_factor_cmplx. On output, the N-by-N lower triangular matrix L.
ldaThe leading dimension of a.
ipvtThe N-element pivot array as output by la_lu_factor_cmplx.
uAn N-by-N matrix where the U matrix will be written.
lduThe leading dimension of u.
pAn N-by-N matrix where the row permutation matrix will be written.
ldpThe leading dimension of p.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, ldu, or ldp is not correct.

◆ la_form_qr()

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.

Parameters
fullqSet 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].
mThe number of rows in R.
nThe number of columns in R.
rOn input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R.
ldrThe leading dimension of matrix R.
tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
qAn 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.
ldqThe leading dimension of matrix Q.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldr or ldq are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_form_qr_cmplx()

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.

Parameters
fullqSet 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].
mThe number of rows in R.
nThe number of columns in R.
rOn input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R.
ldrThe leading dimension of matrix R.
tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
qAn 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.
ldqThe leading dimension of matrix Q.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldr or ldq are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_form_qr_cmplx_pvt()

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

Parameters
fullqSet 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].
mThe number of rows in R.
nThe number of columns in R.
rOn input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R.
ldrThe leading dimension of matrix R.
tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
pvtAn N-element array containing the pivot information from the QR factorization.
qAn 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.
ldqThe leading dimension of matrix Q.
pAn N-by-N matrix where the pivot matrix P will be written.
ldpThe leading dimension of matrix P.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_form_qr_pvt()

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

Parameters
fullqSet 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].
mThe number of rows in R.
nThe number of columns in R.
rOn input, the M-by-N factored matrix as returned by the QR factorization routine. On output, the upper triangular matrix R.
ldrThe leading dimension of matrix R.
tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
pvtAn N-element array containing the pivot information from the QR factorization.
qAn 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.
ldqThe leading dimension of matrix Q.
pAn N-by-N matrix where the pivot matrix P will be written.
ldpThe leading dimension of matrix P.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_inverse()

int la_inverse ( int n,
double * a,
int lda )

Computes the inverse of a square matrix.

Parameters
nThe dimension of matrix A.
aOn input, the N-by-N matrix to invert. On output, the inverted matrix.
ldaThe leading dimension of matrix A.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_SINGULAR_MATRIX_ERROR: Occurs if the input matrix is singular.

◆ la_inverse_cmplx()

int la_inverse_cmplx ( int n,
double complex * a,
int lda )

Computes the inverse of a square matrix.

Parameters
nThe dimension of matrix A.
aOn input, the N-by-N matrix to invert. On output, the inverted matrix.
ldaThe leading dimension of matrix A.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_SINGULAR_MATRIX_ERROR: Occurs if the input matrix is singular.

◆ la_lq_factor()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_lq_factor_cmplx()

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.

WARNING
Untested code that needs more work before being ready for normal use.
Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_lu_factor()

int la_lu_factor ( int m,
int n,
double * a,
int lda,
int * ipvt )

Computes the LU factorization of an M-by-N matrix.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
ipvtAn 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).
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_SINGULAR_MATRIX_ERROR: Occurs as a warning if a is found to be singular.

◆ la_lu_factor_cmplx()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
ipvtAn 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).
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_SINGULAR_MATRIX_ERROR: Occurs as a warning if a is found to be singular.

◆ la_mtx_mult()

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

Parameters
transaSet to true to compute op(A) as the transpose of A; else, set to false to compute op(A) as A.
transbSet to true to compute op(B) as the transpose of B; else, set to false to compute op(B) as B.
mThe number of rows in c.
nThe number of columns in c.
kThe interior dimension of the product a and b.
alphaA scalar multiplier.
aIf transa is true, this matrix must be k by m; else, if transa is false, this matrix must be m by k.
ldaThe leading dimension of matrix a.
bIf transb is true, this matrix must be n by k; else, if transb is false, this matrix must be k by n.
ldbThe leading dimension of matrix b.
betaA scalar multiplier.
cThe m by n matrix C.
ldcThe leading dimension of matrix c.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, ldb, or ldc are not correct.

◆ la_mtx_mult_cmplx()

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

Parameters
opaSet 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.
opbSet 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.
mThenumber of rows in c.
nThe number of columns in c.
kThe interior dimension of the product a and b.
alphaA scalar multiplier.
aIf opa is LA_TRANSPOSE or LA_HERMITIAN_TRANSPOSE, this matrix must be k by m; else, this matrix must be m by k.
ldaThe leading dimension of matrix a.
bIf opb is LA_TRANSPOSE or LA_HERMITIAN_TRANSPOSE, this matrix must be n by k; else, this matrix must be k by n.
ldbThe leading dimension of matrix b.
betaA scalar multiplier.
cThe m by n matrix C.
ldcThe leading dimension of matrix c.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, ldb, or ldc are not correct.

◆ la_mult_lq()

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

Parameters
lsideSet to true to apply \( Q \) or \( Q^T \) from the left; else, set to false to apply \( Q \) or \( Q^T \) from the right.
transSet to true to apply \( Q^T \); else, set to false.
mThe number of rows in matrix C.
nThe number of columns in matrix C.
kThe number of elementary reflectors whose product defines the matrix Q.
aOn 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.
ldaThe leading dimension of matrix A.
tauA K-element array containing the scalar factors of each elementary reflector defined in a.
cOn input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C.
ldcThe leading dimension of matrix C.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldc are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_mult_lq_cmplx()

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

WARNING
Untested code that needs more work before being ready for normal use.
Parameters
lsideSet to true to apply \( Q \) or \( Q^H \) from the left; else, set to false to apply \( Q \) or \( Q^H \) from the right.
transSet to true to apply \( Q^H \); else, set to false.
mThe number of rows in matrix C.
nThe number of columns in matrix C.
kThe number of elementary reflectors whose product defines the matrix Q.
aOn 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.
ldaThe leading dimension of matrix A.
tauA K-element array containing the scalar factors of each elementary reflector defined in a.
cOn input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C.
ldcThe leading dimension of matrix C.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldc are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_mult_qr()

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

Parameters
lsideSet to true to apply \( Q \) or \( Q^T \) from the left; else, set to false to apply \( Q \) or \( Q^T \) from the right.
transSet to true to apply \( Q^T \); else, set to false.
mThe number of rows in matrix C.
nThe number of columns in matrix C.
kThe number of elementary reflectors whose product defines the matrix Q.
aOn 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.
ldaThe leading dimension of matrix A.
tauA K-element array containing the scalar factors of each elementary reflector defined in a.
cOn input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C.
ldcThe leading dimension of matrix C.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldc are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_mult_qr_cmplx()

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

Parameters
lsideSet to true to apply \( Q \) or \( Q^H \) from the left; else, set to false to apply \( Q \) or \( Q^H \) from the right.
transSet to true to apply \( Q^H \); else, set to false.
mThe number of rows in matrix C.
nThe number of columns in matrix C.
kThe number of elementary reflectors whose product defines the matrix Q.
aOn 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.
ldaThe leading dimension of matrix A.
tauA K-element array containing the scalar factors of each elementary reflector defined in a.
cOn input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C.
ldcThe leading dimension of matrix C.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldc are not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_pinverse()

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.

Parameters
mThe number of rows in the matrix. @parma n The number of columns in the matrix.
aOn input, the M-by-N matrix to invert. The matrix is overwritten on output.
ldaThe leading dimension of matrix A.
ainvThe N-by-M matrix where the pseudo-inverse of a will be written.
ldaiThe leading dimension of matrix AINV.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldai is not correct.

◆ la_pinverse_cmplx()

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.

Parameters
mThe number of rows in the matrix. @parma n The number of columns in the matrix.
aOn input, the M-by-N matrix to invert. The matrix is overwritten on output.
ldaThe leading dimension of matrix A.
ainvThe N-by-M matrix where the pseudo-inverse of a will be written.
ldaiThe leading dimension of matrix AINV.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldai is not correct.

◆ la_qr_factor()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_qr_factor_cmplx()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_qr_factor_cmplx_pvt()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
jpvtOn 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.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_qr_factor_pvt()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn 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.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
jpvtOn 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.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_qr_rank1_update()

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

Parameters
mThe number of rows in R.
nThe number of columns in R.
qOn input, the original M-by-K orthogonal matrix Q. On output, the updated matrix Q1.
ldqThe leading dimension of matrix Q.
rOn input, the M-by-N matrix R. On output, the updated matrix R1.
ldrThe leading dimension of matrix R.
uOn input, the M-element U update vector. On output, the original content of the array is overwritten.
vOn input, the N-element V update vector. On output, the original content of the array is overwritten.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldq or ldr is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_qr_rank1_update_cmplx()

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

WARNING
Untested code that needs more work before being ready for normal use.
Parameters
mThe number of rows in R.
nThe number of columns in R.
qOn input, the original M-by-K orthogonal matrix Q. On output, the updated matrix Q1.
ldqThe leading dimension of matrix Q.
rOn input, the M-by-N matrix R. On output, the updated matrix R1.
ldrThe leading dimension of matrix R.
uOn input, the M-element U update vector. On output, the original content of the array is overwritten.
vOn input, the N-element V update vector. On output, the original content of the array is overwritten.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldq or ldr is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if there is insufficient memory available.

◆ la_rank()

int la_rank ( int m,
int n,
double * a,
int lda,
int * rnk )

Computes the rank of a matrix.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aThe M-by-N matrix. The matrix is overwritten as part of this operation.
ldaThe leading dimension of matrix A.
rnkThe rank of a.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs as a warning if the QR iteration process could not converge to a zero value.

◆ la_rank1_update()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
alphaThe scalar multiplier.
xAn M-element array.
yAn N-element array.
aOn input, the M-by-N matrix to update. On output, the updated M-by-N matrix.
ldaThe leading dimension of matrix A.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.

◆ la_rank1_update_cmplx()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
alphaThe scalar multiplier.
xAn M-element array.
yAn N-element array.
aOn input, the M-by-N matrix to update. On output, the updated M-by-N matrix.
ldaThe leading dimension of matrix A.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.

◆ la_rank_cmplx()

int la_rank_cmplx ( int m,
int n,
double complex * a,
int lda,
int * rnk )

Computes the rank of a matrix.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aThe M-by-N matrix. The matrix is overwritten as part of this operation.
ldaThe leading dimension of matrix A.
rnkThe rank of a.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs as a warning if the QR iteration process could not converge to a zero value.

◆ la_solve_cholesky()

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.

Parameters
upperSet 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 \).
mThe number of rows in matrix B.
nThe number of columns in matrix B.
aThe M-by-M Cholesky factored matrix.
ldaThe leading dimension of matrix A.
bOn input, the M-by-N right-hand-side matrix B. On output, the M-by-N solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.

◆ la_solve_cholesky_cmplx()

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.

Parameters
upperSet 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 \).
mThe number of rows in matrix B.
nThe number of columns in matrix B.
aThe M-by-M Cholesky factored matrix.
ldaThe leading dimension of matrix A.
bOn input, the M-by-N right-hand-side matrix B. On output, the M-by-N solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.

◆ la_solve_least_squares()

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.

Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn 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.
ldaThe leading dimension of matrix A.
bIf 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.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_INVALID_OPERATION_ERROR: Occurs if a is not of full rank.

◆ la_solve_least_squares_cmplx()

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.

Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn 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.
ldaThe leading dimension of matrix A.
bIf 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.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_INVALID_OPERATION_ERROR: Occurs if a is not of full rank.

◆ la_solve_lq()

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.

Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn input, the M-by-N QR factored matrix as returned by lq_factor.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by lq_factor.
bOn input, an N-by-K matrix containing the first M rows of the right-hand-side matrix. On output, the solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct, or if m is less than n.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_solve_lq_cmplx()

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.

WARNING
Untested code that needs more work before being ready for normal use.
Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn input, the M-by-N QR factored matrix as returned by lq_factor.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by lq_factor.
bOn input, an N-by-K matrix containing the first M rows of the right-hand-side matrix. On output, the solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct, or if m is less than n.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_solve_lu()

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.

Parameters
mThe number of rows in matrix B.
nThe number of columns in matrix B.
aThe M-by-M LU factored matrix.
ldaThe leading dimension of matrix A.
ipvtThe M-element pivot array from the LU factorization.
bOn input, the M-by-N right-hand-side. On output, the M-by-N solution.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.

◆ la_solve_lu_cmplx()

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.

Parameters
mThe number of rows in matrix B.
nThe number of columns in matrix B.
aThe M-by-M LU factored matrix.
ldaThe leading dimension of matrix A.
ipvtThe M-element pivot array from the LU factorization.
bOn input, the M-by-N right-hand-side. On output, the M-by-N solution.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.

◆ la_solve_qr()

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.

Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor.
bOn input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct, or if m is less than n.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_solve_qr_cmplx()

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.

Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor.
bOn input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct, or if m is less than n.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_solve_qr_cmplx_pvt()

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.

Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor.
jpvtThe N-element array that was used to track the column pivoting operations in the QR factorization.
bOn input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_solve_qr_pvt()

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.

Parameters
mThe number of equations (rows in matrix A).
nThe number of unknowns (columns in matrix A).
kThe number of columns in the right-hand-side matrix.
aOn input, the M-by-N QR factored matrix as returned by qr_factor. On output, the contents of this matrix are restored.
ldaThe leading dimension of matrix A.
tauA MIN(M, N)-element array containing the scalar factors of the elementary reflectors as returned by qr_factor.
jpvtThe N-element array that was used to track the column pivoting operations in the QR factorization.
bOn input, the M-by-K right-hand-side matrix. On output, the first N rows are overwritten by the solution matrix X.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_solve_tri_mtx()

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.

Parameters
lsideSet to true to solve \( op(A) X = \alpha B \); else, set to false to solve \( X op(A) = \alpha B \).
upperSet to true if A is an upper triangular matrix; else, set to false if A is a lower triangular matrix.
transSet to true if \( op(A) = A^T \); else, set to false if \( op(A) = A \).
nounitSet 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.
mThe number of rows in matrix B.
nThe number of columns in matrix B.
alphaThe scalar multiplier to B.
aIf 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.
ldaThe leading dimension of matrix A.
bOn input, the M-by-N right-hand-side. On output, the M-by-N solution.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.

◆ la_solve_tri_mtx_cmplx()

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.

Parameters
lsideSet to true to solve \( op(A) X = \alpha B \); else, set to false to solve \( X op(A) = \alpha B \).
upperSet to true if A is an upper triangular matrix; else, set to false if A is a lower triangular matrix.
transSet to true if \( op(A) = A^H \); else, set to false if \( op(A) = A \).
nounitSet 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.
mThe number of rows in matrix B.
nThe number of columns in matrix B.
alphaThe scalar multiplier to B.
aIf 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.
ldaThe leading dimension of matrix A.
bOn input, the M-by-N right-hand-side. On output, the M-by-N solution.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, or ldb is not correct.

◆ la_sort_eigen()

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.

Parameters
ascend
nThe number of eigenvalues.
valsOn input, an N-element array containing the eigenvalues. On output, the sorted eigenvalues.
vecsOn input, an N-by-N matrix containing the eigenvectors associated with vals (one vector per column). On output, the sorted eigenvector matrix.
ldvThe leading dimension of vecs.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldv is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_sort_eigen_cmplx()

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.

Parameters
ascend
nThe number of eigenvalues.
valsOn input, an N-element array containing the eigenvalues. On output, the sorted eigenvalues.
vecsOn input, an N-by-N matrix containing the eigenvectors associated with vals (one vector per column). On output, the sorted eigenvector matrix.
ldvThe leading dimension of vecs.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if ldv is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.

◆ la_svd()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn input, the M-by-N matrix to factor. The matrix is overwritten on output.
ldaThe leading dimension of matrix A.
sA MIN(M, N)-element array containing the singular values of a sorted in descending order.
uAn M-by-M matrix where the orthogonal U matrix will be written.
lduThe leading dimension of matrix U.
vtAn N-by-N matrix where the transpose of the right singular vector matrix V.
ldvThe leading dimension of matrix V.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, ldu, or ldv is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs as a warning if the QR iteration process could not converge to a zero value.

◆ la_svd_cmplx()

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.

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aOn input, the M-by-N matrix to factor. The matrix is overwritten on output.
ldaThe leading dimension of matrix A.
sA MIN(M, N)-element array containing the singular values of a sorted in descending order.
uAn M-by-M matrix where the orthogonal U matrix will be written.
lduThe leading dimension of matrix U.
vtAn N-by-N matrix where the conjugate transpose of the right singular vector matrix V.
ldvThe leading dimension of vt.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda, ldu, or ldv is not correct.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_CONVERGENCE_ERROR: Occurs as a warning if the QR iteration process could not converge to a zero value.

◆ la_trace()

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

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aThe M-by-N matrix on which to operate.
ldaThe leading dimension of the matrix.
rstThe results of the operation.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.

◆ la_trace_cmplx()

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

Parameters
mThe number of rows in the matrix.
nThe number of columns in the matrix.
aThe M-by-N matrix on which to operate.
ldaThe leading dimension of the matrix.
rstThe results of the operation.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda is not correct.

◆ la_tri_mtx_mult()

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.

Parameters
upperSet 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.
alphaA scalar multiplier.
nThe dimension of the matrix.
aThe 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.
ldaThe leading dimension of matrix A.
betaA scalar multiplier.
bOn input, the n by n matrix B. On output, the n by n resulting matrix.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldb are not correct.

◆ la_tri_mtx_mult_cmplx()

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.

Parameters
upperSet 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.
alphaA scalar multiplier.
nThe dimension of the matrix.
aThe 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.
ldaThe leading dimension of matrix A.
betaA scalar multiplier.
bOn input, the n by n matrix B. On output, the n by n resulting matrix.
ldbThe leading dimension of matrix B.
Returns
An error code. The following codes are possible.
  • LA_NO_ERROR: No error occurred. Successful operation.
  • LA_INVALID_INPUT_ERROR: Occurs if lda or ldb are not correct.