<T>LAPACK 0.1.2
C++ Template Linear Algebra PACKage
|
Functions | |
template<typename TA , typename TT > | |
int | tlapack::legacy::geqr2 (idx_t m, idx_t n, TA *A, idx_t lda, TT *tau) |
Computes a QR factorization of a matrix A. | |
template<typename TA , typename TB > | |
void | tlapack::legacy::lacpy (MatrixType matrixtype, idx_t m, idx_t n, const TA *A, idx_t lda, TB *B, idx_t ldb) |
Copies a real matrix from A to B where A is either a full, upper triangular or lower triangular matrix. | |
template<class norm_t , typename TA > | |
real_type< TA > | tlapack::legacy::lange (norm_t normType, idx_t m, idx_t n, const TA *A, idx_t lda) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value. | |
template<class norm_t , typename TA > | |
real_type< TA > | tlapack::legacy::lanhe (norm_t normType, Uplo uplo, idx_t n, const TA *A, idx_t lda) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value of a symmetric matrix. | |
template<class norm_t , typename TA > | |
real_type< TA > | tlapack::legacy::lansy (norm_t normType, Uplo uplo, idx_t n, const TA *A, idx_t lda) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value of a symmetric matrix. | |
template<typename TA > | |
real_type< TA > | tlapack::legacy::lantr (Norm normType, Uplo uplo, Diag diag, idx_t m, idx_t n, const TA *A, idx_t lda) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value of a triangular matrix. | |
template<class side_t , typename TV , typename TC > | |
void | tlapack::legacy::larf (side_t side, idx_t m, idx_t n, TV const *v, int_t incv, scalar_type< TV, TC > tau, TC *C, idx_t ldC) |
Applies an elementary reflector H to a m-by-n matrix C. | |
template<class side_t , class trans_t , class direction_t , class storev_t , typename TV , typename TC > | |
int | tlapack::legacy::larfb (side_t side, trans_t trans, direction_t direction, storev_t storev, idx_t m, idx_t n, idx_t k, TV const *V, idx_t ldv, TV const *T, idx_t ldt, TC *C, idx_t ldc) |
Applies a block reflector \(H\) or its conjugate transpose \(H^H\) to a m-by-n matrix C, from either the left or the right. | |
template<typename T > | |
void | tlapack::legacy::larfg (idx_t n, T *alpha, T *x, int_t incx, T *tau) |
Generates a elementary Householder reflection. | |
template<class direction_t , class storev_t , typename scalar_t > | |
int | tlapack::legacy::larft (direction_t direction, storev_t storev, idx_t n, idx_t k, const scalar_t *V, idx_t ldV, const scalar_t *tau, scalar_t *T, idx_t ldT) |
Forms the triangular factor T of a block reflector H of order n, which is defined as a product of k elementary reflectors. | |
template<typename T > | |
void | tlapack::legacy::larnv (idx_t idist, idx_t *iseed, idx_t n, T *x) |
Returns a vector of n random numbers from a uniform or normal distribution. | |
template<class matrixtype_t , typename T > | |
int | tlapack::legacy::lascl (matrixtype_t matrixtype, idx_t kl, idx_t ku, const real_type< T > &b, const real_type< T > &a, idx_t m, idx_t n, T *A, idx_t lda) |
Multiplies a matrix A by the real scalar a/b. | |
template<typename TA > | |
void | tlapack::legacy::laset (MatrixType matrixtype, idx_t m, idx_t n, TA alpha, TA beta, TA *A, idx_t lda) |
Initializes a matrix to diagonal and off-diagonal values. | |
template<typename TX > | |
void | tlapack::legacy::lassq (idx_t n, TX const *x, int_t incx, real_type< TX > &scl, real_type< TX > &sumsq) |
Updates a sum of squares represented in scaled form. | |
template<class uplo_t , typename T > | |
int | tlapack::legacy::potrf (uplo_t uplo, idx_t n, T *A, idx_t lda) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using a blocked algorithm. | |
template<class uplo_t , typename T > | |
int | tlapack::legacy::potrs (uplo_t uplo, idx_t n, idx_t nrhs, const T *A, idx_t lda, T *B, idx_t ldb) |
Apply the Cholesky factorization to solve a linear system. | |
template<typename TA , typename TT > | |
int | tlapack::legacy::ung2r (idx_t m, idx_t n, idx_t k, TA *A, idx_t lda, const TT *tau) |
Generates a m-by-n matrix Q with orthogonal columns. | |
template<class side_t , class trans_t , typename TA , typename TC > | |
int | tlapack::legacy::unm2r (side_t side, trans_t trans, idx_t m, idx_t n, idx_t k, TA *A, idx_t lda, const real_type< TA, TC > *tau, TC *C, idx_t ldc) |
Applies unitary matrix Q to a matrix C. | |
template<class side_t , class trans_t , typename TA , typename TC > | |
int | tlapack::legacy::unmqr (side_t side, trans_t trans, idx_t m, idx_t n, idx_t k, TA const *A, idx_t lda, TA const *tau, TC *C, idx_t ldc) |
Multiplies the general m-by-n matrix C by Q from geqrf() using a blocked code as follows: | |
Computes a QR factorization of a matrix A.
[in] | m | The number of rows of the matrix A. |
[in] | n | The number of columns of the matrix A. |
[in,out] | A | m-by-n matrix. On exit, the elements on and above the diagonal of the array contain the min(m,n)-by-n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal, with the array tau, represent the unitary matrix Q as a product of elementary reflectors. |
[in] | lda | The leading dimension of A. lda >= max(1,m). |
[out] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. The subarray tau[1:n-1] is used as the workspace. |
void tlapack::legacy::lacpy | ( | MatrixType | matrixtype, |
idx_t | m, | ||
idx_t | n, | ||
const TA * | A, | ||
idx_t | lda, | ||
TB * | B, | ||
idx_t | ldb | ||
) |
Copies a real matrix from A to B where A is either a full, upper triangular or lower triangular matrix.
[in] | matrixtype | : 'U': A is assumed to be upper triangular; elements below thediagonal are not referenced. 'L': A is assumed to be lower triangular; elements above the diagonal are not referenced. otherwise, A is assumed to be a full matrix. |
[in] | m | Number of rows of A. |
[in] | n | Number of columns of A. |
[in] | A | m-by-n matrix. |
[in] | lda | Leading dimension of A. |
[out] | B | Matrix with at least m rows and at least n columns. |
[in] | ldb | Leading dimension of B. |
real_type< TA > tlapack::legacy::lange | ( | norm_t | normType, |
idx_t | m, | ||
idx_t | n, | ||
const TA * | A, | ||
idx_t | lda | ||
) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value.
normType | Type should be specified as follows: Norm::Max = maximum absolute value over all elements in A. Note: this is not a consistent matrix norm. Norm::One = one norm of the matrix A, the maximum value of the sumsof each column. Norm::Inf = the infinity norm of the matrix A, the maximum value of the sum of each row. Norm::Fro = the Frobenius norm of the matrix A. This the square root of the sum of the squares of each element in A. |
m | Number of rows to be included in the norm. m >= 0 |
n | Number of columns to be included in the norm. n >= 0 |
A | matrix size m-by-n. |
lda | Column length of the matrix A. ldA >= m |
real_type< TA > tlapack::legacy::lanhe | ( | norm_t | normType, |
Uplo | uplo, | ||
idx_t | n, | ||
const TA * | A, | ||
idx_t | lda | ||
) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value of a symmetric matrix.
normType | Type should be specified as follows: Norm::Max = maximum absolute value over all elements in A. Note: this is not a consistent matrix norm. Norm::One = one norm of the matrix A, the maximum value of the sumsof each column. Norm::Inf = the infinity norm of the matrix A, the maximum value of the sum of each row. Norm::Fro = the Frobenius norm of the matrix A. This the square root of the sum of the squares of each element in A. |
uplo | Indicates whether the symmetric matrix A is stored as upper triangular or lower triangular. The other strict triangular part of A is not referenced. |
n | Number of columns to be included in the norm. n >= 0 |
A | symmetric matrix size lda-by-n. |
lda | Leading dimension of matrix A. ldA >= m |
real_type< TA > tlapack::legacy::lansy | ( | norm_t | normType, |
Uplo | uplo, | ||
idx_t | n, | ||
const TA * | A, | ||
idx_t | lda | ||
) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value of a symmetric matrix.
normType | Type should be specified as follows: Norm::Max = maximum absolute value over all elements in A. Note: this is not a consistent matrix norm. Norm::One = one norm of the matrix A, the maximum value of the sumsof each column. Norm::Inf = the infinity norm of the matrix A, the maximum value of the sum of each row. Norm::Fro = the Frobenius norm of the matrix A. This the square root of the sum of the squares of each element in A. |
uplo | Indicates whether the symmetric matrix A is stored as upper triangular or lower triangular. The other triangular part of A is not referenced. |
n | Number of columns to be included in the norm. n >= 0 |
A | symmetric matrix size lda-by-n. |
lda | Leading dimension of matrix A. ldA >= m |
real_type< TA > tlapack::legacy::lantr | ( | Norm | normType, |
Uplo | uplo, | ||
Diag | diag, | ||
idx_t | m, | ||
idx_t | n, | ||
const TA * | A, | ||
idx_t | lda | ||
) |
Calculates the value of the one norm, Frobenius norm, infinity norm, or element of largest absolute value of a triangular matrix.
normType | Type should be specified as follows: Norm::Max = maximum absolute value over all elements in A. Note: this is not a consistent matrix norm. Norm::One = one norm of the matrix A, the maximum value of the sumsof each column. Norm::Inf = the infinity norm of the matrix A, the maximum value of the sum of each row. Norm::Fro = the Frobenius norm of the matrix A. This the square root of the sum of the squares of each element in A. | |
uplo | Indicates whether A is upper or lower triangular. The other strict triangular part of A is not referenced. | |
[in] | diag | Whether A has a unit or non-unit diagonal:
|
m | Number of rows to be included in the norm. m >= 0 | |
n | Number of columns to be included in the norm. n >= 0 | |
A | symmetric matrix size lda-by-n. | |
lda | Leading dimension of matrix A. ldA >= m |
void tlapack::legacy::larf | ( | side_t | side, |
idx_t | m, | ||
idx_t | n, | ||
TV const * | v, | ||
int_t | incv, | ||
scalar_type< TV, TC > | tau, | ||
TC * | C, | ||
idx_t | ldC | ||
) |
Applies an elementary reflector H to a m-by-n matrix C.
int tlapack::legacy::larfb | ( | side_t | side, |
trans_t | trans, | ||
direction_t | direction, | ||
storev_t | storev, | ||
idx_t | m, | ||
idx_t | n, | ||
idx_t | k, | ||
TV const * | V, | ||
idx_t | ldv, | ||
TV const * | T, | ||
idx_t | ldt, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Applies a block reflector \(H\) or its conjugate transpose \(H^H\) to a m-by-n matrix C, from either the left or the right.
[in] | side |
|
[in] | trans |
|
[in] | direction | Indicates how H is formed from a product of elementary reflectors
|
[in] | storev | Indicates how the vectors which define the elementary reflectors are stored: |
[in] | m | The number of rows of the matrix C. |
[in] | n | The number of columns of the matrix C. |
[in] | k | The order of the matrix T (= the number of elementary reflectors whose product defines the block reflector).
|
[in] | V |
|
[in] | ldv | The leading dimension of the array V.
|
[in] | T | The k-by-k matrix T, stored in an ldt-by-k array. The triangular k-by-k matrix T in the representation of the block reflector. |
[in] | ldt | The leading dimension of the array T. ldt >= k. |
[in,out] | C | The m-by-n matrix C, stored in an ldc-by-n array. On entry, the m-by-n matrix C. On exit, C is overwritten by \(H C\) or \(H^H C\) or \(C H\) or \(C H^H\). |
[in] | ldc | The leading dimension of the array C. ldc >= max(1,m). |
The shape of the matrix V and the storage of the vectors which define the H(i) is best illustrated by the following example with n = 5 and k = 3. The elements equal to 1 are not stored. The rest of the array is not used.
direction = Forward and direction = Forward and storev = Columnwise: storev = Rowwise: V = ( 1 ) V = ( 1 v1 v1 v1 v1 ) ( v1 1 ) ( 1 v2 v2 v2 ) ( v1 v2 1 ) ( 1 v3 v3 ) ( v1 v2 v3 ) ( v1 v2 v3 ) direction = Backward and direction = Backward and storev = Columnwise: storev = Rowwise: V = ( v1 v2 v3 ) V = ( v1 v1 1 ) ( v1 v2 v3 ) ( v2 v2 v2 1 ) ( 1 v2 v3 ) ( v3 v3 v3 v3 1 ) ( 1 v3 ) ( 1 )
void tlapack::legacy::larfg | ( | idx_t | n, |
T * | alpha, | ||
T * | x, | ||
int_t | incx, | ||
T * | tau | ||
) |
Generates a elementary Householder reflection.
int tlapack::legacy::larft | ( | direction_t | direction, |
storev_t | storev, | ||
idx_t | n, | ||
idx_t | k, | ||
const scalar_t * | V, | ||
idx_t | ldV, | ||
const scalar_t * | tau, | ||
scalar_t * | T, | ||
idx_t | ldT | ||
) |
Forms the triangular factor T of a block reflector H of order n, which is defined as a product of k elementary reflectors.
If direction = Direction::Forward, H = H_1 H_2 . . . H_k
and T is upper triangular. If direction = Direction::Backward, H = H_k . . . H_2 H_1 and T is lower triangular.
If storev = StoreV::Columnwise, the vector which defines the elementary reflector H(i) is stored in the i-th column of the array V, and
H = I - V * T * V'
If storev = StoreV::Rowwise, the vector which defines the elementary reflector H(i) is stored in the i-th row of the array V, and
H = I - V' * T * V
The shape of the matrix V and the storage of the vectors which define the H(i) is best illustrated by the following example with n = 5 and k = 3. The elements equal to 1 are not stored.
direction=Direction::Forward & storev=StoreV::Columnwise
direction=Direction::Forward & storev=StoreV::Rowwise
V = ( 1 ) V = ( 1 v1 v1 v1 v1 ) ( v1 1 ) ( 1 v2 v2 v2 ) ( v1 v2 1 ) ( 1 v3 v3 ) ( v1 v2 v3 ) ( v1 v2 v3 )
direction=Direction::Backward & storev=StoreV::Columnwise direction=Direction::Backward & storev=StoreV::Rowwise
V = ( v1 v2 v3 ) V = ( v1 v1 1 ) ( v1 v2 v3 ) ( v2 v2 v2 1 ) ( 1 v2 v3 ) ( v3 v3 v3 v3 1 ) ( 1 v3 ) ( 1 )
direction | Specifies the direction in which the elementary reflectors are multiplied to form the block reflector. Direction::Forward Direction::Backward | |
storev | Specifies how the vectors which define the elementary reflectors are stored. StoreV::Columnwise StoreV::Rowwise | |
n | The order of the block reflector H. n >= 0. | |
k | The order of the triangular factor T, or the number of elementary reflectors. k >= 1. | |
[in] | V | Real matrix containing the vectors defining the elementary reflector H. If stored columnwise, V is n-by-k. If stored rowwise, V is k-by-n. |
ldV | Column length of the matrix V. If stored columnwise, ldV >= n. If stored rowwise, ldV >= k. | |
[in] | tau | Real vector of length k containing the scalar factors of the elementary reflectors H. |
[out] | T | Real matrix of size k-by-k containing the triangular factor of the block reflector. If the direction of the elementary reflectors is forward, T is upper triangular; if the direction of the elementary reflectors is backward, T is lower triangular. |
ldT | Column length of the matrix T. ldT >= k. |
Returns a vector of n random numbers from a uniform or normal distribution.
This implementation uses the Mersenne Twister 19937 generator (class std::mt19937), which is a Mersenne Twister pseudo-random generator of 32-bit numbers with a state size of 19937 bits.
Requires ISO C++ 2011 random number generators.
[in] | idist | Specifies the distribution: 1: real and imaginary parts each uniform (0,1) 2: real and imaginary parts each uniform (-1,1) 3: real and imaginary parts each normal (0,1) 4: uniformly distributed on the disc abs(z) < 1 5: uniformly distributed on the circle abs(z) = 1 |
[in,out] | iseed | Seed for the random number generator. The seed is updated inside the routine ( Currently: seed_out := seed_in
|
[in] | n | Length of vector x. |
[out] | x | Pointer to real vector of length n. |
int tlapack::legacy::lascl | ( | matrixtype_t | matrixtype, |
idx_t | kl, | ||
idx_t | ku, | ||
const real_type< T > & | b, | ||
const real_type< T > & | a, | ||
idx_t | m, | ||
idx_t | n, | ||
T * | A, | ||
idx_t | lda | ||
) |
Multiplies a matrix A by the real scalar a/b.
Multiplication of a matrix A by scalar a/b is done without over/underflow as long as the final result \(a A/b\) does not over/underflow. The parameter type specifies that A may be full, upper triangular, lower triangular, upper Hessenberg, or banded.
[in] | matrixtype | Specifies the type of matrix A. MatrixType::General: A is a full matrix. MatrixType::Lower: A is a lower triangular matrix. MatrixType::Upper: A is an upper triangular matrix. MatrixType::Hessenberg: A is an upper Hessenberg matrix. MatrixType::LowerBand: A is a symmetric band matrix with lower bandwidth kl and upperbandwidth ku and with the only the lower half stored. MatrixType::UpperBand: A is a symmetric band matrix with lower bandwidth kl and upper bandwidth ku and with the only the upper half stored. MatrixType::Band: A is a band matrix with lower bandwidth kl and upper bandwidth ku. |
[in] | kl | The lower bandwidth of A, used only for banded matrix types B, Q and Z. |
[in] | ku | The upper bandwidth of A, used only for banded matrix types B, Q and Z. |
[in] | b | The denominator of the scalar a/b. |
[in] | a | The numerator of the scalar a/b. |
[in] | m | The number of rows of the matrix A. m>=0 |
[in] | n | The number of columns of the matrix A. n>=0 |
[in,out] | A | Pointer to the matrix A [in/out]. |
[in] | lda | The column length of the matrix A. |
void tlapack::legacy::laset | ( | MatrixType | matrixtype, |
idx_t | m, | ||
idx_t | n, | ||
TA | alpha, | ||
TA | beta, | ||
TA * | A, | ||
idx_t | lda | ||
) |
Initializes a matrix to diagonal and off-diagonal values.
[in] | matrixtype | : 'U': A is assumed to be upper triangular; elements below thediagonal are not referenced. 'L': A is assumed to be lower triangular; elements above the diagonal are not referenced. otherwise, A is assumed to be a full matrix. |
[in] | m | Number of rows of A. |
[in] | n | Number of columns of A. |
[in] | alpha | Value to be assigned to the off-diagonal elements of A. |
[in] | beta | Value to assign to the diagonal elements of A. |
[in] | A | m-by-n matrix. |
[in] | lda | Leading dimension of A. |
void tlapack::legacy::lassq | ( | idx_t | n, |
TX const * | x, | ||
int_t | incx, | ||
real_type< TX > & | scl, | ||
real_type< TX > & | sumsq | ||
) |
Updates a sum of squares represented in scaled form.
\[ scl_{[OUT]}^2 sumsq_{[OUT]} = \sum_{i = 0}^n x_i^2 + scl_{[IN]}^2 sumsq_{[IN]}, \]
The value of sumsq is assumed to be non-negative.
If (scale * sqrt( sumsq )) > tbig on entry then we require: scale >= sqrt( TINY*EPS ) / sbig on entry, and if 0 < (scale * sqrt( sumsq )) < tsml on entry then we require: scale <= sqrt( HUGE ) / ssml on entry, where tbig – upper threshold for values whose square is representable; sbig – scaling constant for big numbers;
[in] | n | The number of elements to be used from the vector x. |
[in] | x | Array of dimension \((1+(n-1)*|incx|)\). |
[in] | incx | The increment between successive values of the vector x. If incx > 0, X(i*incx) = x_i for 0 <= i < n If incx < 0, X((n-i-1)*(-incx)) = x_i for 0 <= i < n If incx = 0, x isn't a vector so there is no need to call this subroutine. If you call it anyway, it will count x_0 in the vector norm n times. |
[in] | scl | |
[in] | sumsq |
int tlapack::legacy::potrs | ( | uplo_t | uplo, |
idx_t | n, | ||
idx_t | nrhs, | ||
const T * | A, | ||
idx_t | lda, | ||
T * | B, | ||
idx_t | ldb | ||
) |
Apply the Cholesky factorization to solve a linear system.
Generates a m-by-n matrix Q with orthogonal columns.
\[ Q = H_1 H_2 ... H_k \]
[in] | m | The number of rows of the matrix A. m>=0 |
[in] | n | The number of columns of the matrix A. n>=0 |
[in] | k | The number of elementary reflectors whose product defines the matrix Q. n>=k>=0 |
[in,out] | A | m-by-n matrix. On entry, the i-th column must contains the vector which defines the elementary reflector \(H_i\), for \(i=0,1,...,k-1\), as returned by GEQRF in the first k columns of its array argument A. On exit, the m-by-n matrix \(Q = H_1 H_2 ... H_k\). |
[in] | lda | The leading dimension of A. lda >= max(1,m). |
[in] | tau | Real vector of length min(m,n). The scalar factors of the elementary reflectors. |
int tlapack::legacy::unm2r | ( | side_t | side, |
trans_t | trans, | ||
idx_t | m, | ||
idx_t | n, | ||
idx_t | k, | ||
TA * | A, | ||
idx_t | lda, | ||
const real_type< TA, TC > * | tau, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Applies unitary matrix Q to a matrix C.
[in] | side | Specifies which side Q is to be applied. 'L': apply Q or Q' from the Left; 'R': apply Q or Q' from the Right. |
[in] | trans | Specifies whether Q or Q' is applied. 'N': No transpose, apply Q; 'T': Transpose, apply Q'. |
[in] | m | The number of rows of the matrix C. |
[in] | n | The number of columns of the matrix C. |
[in] | k | The number of elementary reflectors whose product defines the matrix Q. If side='L', m>=k>=0; if side='R', n>=k>=0. |
[in] | A | Matrix containing the elementary reflectors H. If side='L', A is k-by-m; if side='R', A is k-by-n. |
[in] | lda | The column length of the matrix A. ldA>=k. |
[in] | tau | Real vector of length k containing the scalar factors of the elementary reflectors. |
[in,out] | C | m-by-n matrix. On exit, C is replaced by one of the following: If side='L' & trans='N': C <- Q * C If side='L' & trans='T': C <- Q'* C If side='R' & trans='T': C <- C * Q' If side='R' & trans='N': C <- C * Q |
ldc | The column length the matrix C. ldC>=m. |
int tlapack::legacy::unmqr | ( | side_t | side, |
trans_t | trans, | ||
idx_t | m, | ||
idx_t | n, | ||
idx_t | k, | ||
TA const * | A, | ||
idx_t | lda, | ||
TA const * | tau, | ||
TC * | C, | ||
idx_t | ldc | ||
) |
Multiplies the general m-by-n matrix C by Q from geqrf() using a blocked code as follows:
[in] | side |
|
[in] | trans |
|
[in] | m | The number of rows of the matrix C. m >= 0. |
[in] | n | The number of columns of the matrix C. n >= 0. |
[in] | k | The number of elementary reflectors whose product defines the matrix Q.
|
[in] | A |
|
[in] | lda | The leading dimension of the array A.
|
[in] | tau | The vector tau of length k. tau[i] must contain the scalar factor of the elementary reflector H(i), as returned by geqrf(). |
[in,out] | C | The m-by-n matrix C, stored in an ldc-by-n array. On entry, the m-by-n matrix C. On exit, C is overwritten by \(Q C\) or \(Q^H C\) or \(C Q^H\) or \(C Q\). |
[in] | ldc | The leading dimension of the array C. ldc >= max(1,m). |