![]() |
LAPACK LQ decomposition functions (gelqf) | |
void | blaze::gelqf (int m, int n, float *A, int lda, float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the LQ decomposition of the given dense single precision column-major matrix. More... | |
void | blaze::gelqf (int m, int n, double *A, int lda, double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the LQ decomposition of the given dense double precision column-major matrix. More... | |
void | blaze::gelqf (int m, int n, complex< float > *A, int lda, complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the LQ decomposition of the given dense single precision complex column-major matrix. More... | |
void | blaze::gelqf (int m, int n, complex< double > *A, int lda, complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the LQ decomposition of the given dense double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::gelqf (DenseMatrix< MT, SO > &A, typename MT::ElementType *tau) |
LAPACK kernel for the LQ decomposition of the given dense matrix. More... | |
LAPACK QL decomposition functions (geqlf) | |
void | blaze::geqlf (int m, int n, float *A, int lda, float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the QL decomposition of the given dense single precision column-major matrix. More... | |
void | blaze::geqlf (int m, int n, double *A, int lda, double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the QL decomposition of the given dense double precision column-major matrix. More... | |
void | blaze::geqlf (int m, int n, complex< float > *A, int lda, complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the QL decomposition of the given dense single precision complex column-major matrix. More... | |
void | blaze::geqlf (int m, int n, complex< double > *A, int lda, complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the QL decomposition of the given dense double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::geqlf (DenseMatrix< MT, SO > &A, typename MT::ElementType *tau) |
LAPACK kernel for the QL decomposition of the given dense matrix. More... | |
LAPACK QR decomposition functions (geqrf) | |
void | blaze::geqrf (int m, int n, float *A, int lda, float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the QR decomposition of the given dense single precision column-major matrix. More... | |
void | blaze::geqrf (int m, int n, double *A, int lda, double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the QR decomposition of the given dense double precision column-major matrix. More... | |
void | blaze::geqrf (int m, int n, complex< float > *A, int lda, complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the QR decomposition of the given dense single precision complex column-major matrix. More... | |
void | blaze::geqrf (int m, int n, complex< double > *A, int lda, complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the QR decomposition of the given dense double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::geqrf (DenseMatrix< MT, SO > &A, typename MT::ElementType *tau) |
LAPACK kernel for the QR decomposition of the given dense matrix. More... | |
LAPACK RQ decomposition functions (gerqf) | |
void | blaze::gerqf (int m, int n, float *A, int lda, float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the RQ decomposition of the given dense single precision column-major matrix. More... | |
void | blaze::gerqf (int m, int n, double *A, int lda, double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the RQ decomposition of the given dense single precision column-major matrix. More... | |
void | blaze::gerqf (int m, int n, complex< float > *A, int lda, complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the RQ decomposition of the given dense single precision complex column-major matrix. More... | |
void | blaze::gerqf (int m, int n, complex< double > *A, int lda, complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the RQ decomposition of the given dense double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::gerqf (DenseMatrix< MT, SO > &A, typename MT::ElementType *tau) |
LAPACK kernel for the RQ decomposition of the given dense matrix. More... | |
LAPACK LU decomposition functions (getrf) | |
void | blaze::getrf (int m, int n, float *A, int lda, int *ipiv, int *info) |
LAPACK kernel for the LU decomposition of the given dense general single precision column-major matrix. More... | |
void | blaze::getrf (int m, int n, double *A, int lda, int *ipiv, int *info) |
LAPACK kernel for the LU decomposition of the given dense general double precision column-major matrix. More... | |
void | blaze::getrf (int m, int n, complex< float > *A, int lda, int *ipiv, int *info) |
LAPACK kernel for the LU decomposition of the given dense general single precision complex column-major matrix. More... | |
void | blaze::getrf (int m, int n, complex< double > *A, int lda, int *ipiv, int *info) |
LAPACK kernel for the LU decomposition of the given dense general double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::getrf (DenseMatrix< MT, SO > &A, int *ipiv) |
LAPACK LDLH decomposition functions (hetrf) | |
void | blaze::hetrf (char uplo, int n, complex< float > *A, int lda, int *ipiv, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the decomposition of the given dense Hermitian indefinite single precision complex column-major matrix. More... | |
void | blaze::hetrf (char uplo, int n, complex< double > *A, int lda, int *ipiv, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the decomposition of the given dense Hermitian indefinite double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::hetrf (DenseMatrix< MT, SO > &A, char uplo, int *ipiv) |
LAPACK kernel for the decomposition of the given dense Hermitian indefinite matrix. More... | |
LAPACK functions to reconstruct Q from a LQ decomposition (orglq) | |
void | blaze::orglq (int m, int n, int k, float *A, int lda, const float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition. More... | |
void | blaze::orglq (int m, int n, int k, double *A, int lda, const double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::orglq (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition. More... | |
LAPACK functions to reconstruct Q from a QL decomposition (orgql) | |
void | blaze::orgql (int m, int n, int k, float *A, int lda, const float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition. More... | |
void | blaze::orgql (int m, int n, int k, double *A, int lda, const double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::orgql (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition. More... | |
LAPACK functions to reconstruct Q from a QR decomposition (orgqr) | |
void | blaze::orgqr (int m, int n, int k, float *A, int lda, const float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition. More... | |
void | blaze::orgqr (int m, int n, int k, double *A, int lda, const double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::orgqr (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition. More... | |
LAPACK functions to reconstruct Q from a RQ decomposition (orgrq) | |
void | blaze::orgrq (int m, int n, int k, float *A, int lda, const float *tau, float *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition. More... | |
void | blaze::orgrq (int m, int n, int k, double *A, int lda, const double *tau, double *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::orgrq (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition. More... | |
LAPACK LLH (Cholesky) decomposition functions (potrf) | |
void | blaze::potrf (char uplo, int n, float *A, int lda, int *info) |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite single precision column-major matrix. More... | |
void | blaze::potrf (char uplo, int n, double *A, int lda, int *info) |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite double precision column-major matrix. More... | |
void | blaze::potrf (char uplo, int n, complex< float > *A, int lda, int *info) |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite single precision complex column-major matrix. More... | |
void | blaze::potrf (char uplo, int n, complex< double > *A, int lda, int *info) |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::potrf (DenseMatrix< MT, SO > &A, char uplo) |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite matrix. More... | |
LAPACK LDLT decomposition functions (sytrf) | |
void | blaze::sytrf (char uplo, int n, float *A, int lda, int *ipiv, float *work, int lwork, int *info) |
LAPACK kernel for the decomposition of the given dense symmetric indefinite single precision column-major matrix. More... | |
void | blaze::sytrf (char uplo, int n, double *A, int lda, int *ipiv, double *work, int lwork, int *info) |
LAPACK kernel for the decomposition of the given dense symmetric indefinite double precision column-major matrix. More... | |
void | blaze::sytrf (char uplo, int n, complex< float > *A, int lda, int *ipiv, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the decomposition of the given dense symmetric indefinite single precision complex column-major matrix. More... | |
void | blaze::sytrf (char uplo, int n, complex< double > *A, int lda, int *ipiv, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the decomposition of the given dense symmetric indefinite double precision complex column-major matrix. More... | |
template<typename MT , bool SO> | |
void | blaze::sytrf (DenseMatrix< MT, SO > &A, char uplo, int *ipiv) |
LAPACK kernel for the decomposition of the given dense symmetric indefinite matrix. More... | |
LAPACK functions to reconstruct Q from a LQ decomposition (unglq) | |
void | blaze::unglq (int m, int n, int k, complex< float > *A, int lda, const complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition. More... | |
void | blaze::unglq (int m, int n, int k, complex< double > *A, int lda, const complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::unglq (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition. More... | |
LAPACK functions to reconstruct Q from a QL decomposition (ungql) | |
void | blaze::ungql (int m, int n, int k, complex< float > *A, int lda, const complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition. More... | |
void | blaze::ungql (int m, int n, int k, complex< double > *A, int lda, const complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::ungql (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition. More... | |
LAPACK functions to reconstruct Q from a QR decomposition (ungqr) | |
void | blaze::ungqr (int m, int n, int k, complex< float > *A, int lda, const complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition. More... | |
void | blaze::ungqr (int m, int n, int k, complex< double > *A, int lda, const complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::ungqr (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition. More... | |
LAPACK functions to reconstruct Q from a RQ decomposition (ungrq) | |
void | blaze::ungrq (int m, int n, int k, complex< float > *A, int lda, const complex< float > *tau, complex< float > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition. More... | |
void | blaze::ungrq (int m, int n, int k, complex< double > *A, int lda, const complex< double > *tau, complex< double > *work, int lwork, int *info) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition. More... | |
template<typename MT , bool SO> | |
void | blaze::ungrq (DenseMatrix< MT, SO > &A, const typename MT::ElementType *tau) |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition. More... | |
|
inline |
LAPACK kernel for the LQ decomposition of the given dense single precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix LQ decomposition of a general m-by-n single precision column-major matrix based on the LAPACK sgelqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:n)
is stored on exit in A(i,i+1:n)
, and tau
in tau(i)
. Thus on exit the elements on and below the diagonal of the matrix contain the m-by-min(m,n) lower trapezoidal matrix L
(L
is lower triangular if m <= n); the elements above the diagonal, with the array tau
, represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgelqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LQ decomposition of the given dense double precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix LQ decomposition of a general m-by-n double precision column-major matrix based on the LAPACK sgelqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:n)
is stored on exit in A(i,i+1:n)
, and tau
in tau(i)
. Thus on exit the elements on and below the diagonal of the matrix contain the m-by-min(m,n) lower trapezoidal matrix L
(L
is lower triangular if m <= n); the elements above the diagonal, with the array tau
, represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgelqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LQ decomposition of the given dense single precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix LQ decomposition of a general m-by-n single precision complex column-major matrix based on the LAPACK sgelqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:n)
is stored on exit in A(i,i+1:n)
, and tau
in tau(i)
. Thus on exit the elements on and below the diagonal of the matrix contain the m-by-min(m,n) lower trapezoidal matrix L
(L
is lower triangular if m <= n); the elements above the diagonal, with the array tau
, represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgelqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LQ decomposition of the given dense double precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix LQ decomposition of a general m-by-n double precision complex column-major matrix based on the LAPACK sgelqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:n)
is stored on exit in A(i,i+1:n)
, and tau
in tau(i)
. Thus on exit the elements on and below the diagonal of the matrix contain the m-by-min(m,n) lower trapezoidal matrix L
(L
is lower triangular if m <= n); the elements above the diagonal, with the array tau
, represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgelqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LQ decomposition of the given dense matrix.
A | The matrix to be decomposed. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function performs the dense matrix LQ decomposition of a general m-by-n matrix based on the LAPACK gelqf() functions. Note that this function can only be used for general, non-adapted matrices with float
, double
, complex<float>
, or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
In case of a column-major matrix, the resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:n)
is stored on exit in A(i,i+1:n)
, and tau
in tau(i)
. Thus on exit the elements on and below the diagonal of the matrix contain the m-by-min(m,n) lower trapezoidal matrix L
(L
is lower triangular if m <= n); the elements above the diagonal, with the array tau
, represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
In case of a row-major matrix, the resulting decomposition is transposed, i.e. the elementary reflectors are stored below the diagonal and the elements on and above the diagonal contain the min(m,n)-by-m upper trapezoidal matrix L
.
For more information on the gelqf() functions (i.e. sgelqf(), dgelqf(), cgelqf(), and zgelqf()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QL decomposition of the given dense single precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QL decomposition of a general m-by-n single precision column-major matrix based on the LAPACK sgeqlf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(m-k+i+1:m) = 0
and v(m-k+i) = 1
. v(1:m-k+i-1)
is stored on exit in A(1:m-k+i-1,n-k+i)
, and tau
in tau(i)
. Thus in case m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L
and in case m <= n, the elements on and below the (n-m)-th subdiagonal contain the m-by-n lower trapezoidal matrix L
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgeqlf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QL decomposition of the given dense double precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QL decomposition of a general m-by-n double precision column-major matrix based on the LAPACK sgeqlf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(m-k+i+1:m) = 0
and v(m-k+i) = 1
. v(1:m-k+i-1)
is stored on exit in A(1:m-k+i-1,n-k+i)
, and tau
in tau(i)
. Thus in case m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L
and in case m <= n, the elements on and below the (n-m)-th subdiagonal contain the m-by-n lower trapezoidal matrix L
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgeqlf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QL decomposition of the given dense single precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QL decomposition of a general m-by-n single precision complex column-major matrix based on the LAPACK sgeqlf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(m-k+i+1:m) = 0
and v(m-k+i) = 1
. v(1:m-k+i-1)
is stored on exit in A(1:m-k+i-1,n-k+i)
, and tau
in tau(i)
. Thus in case m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L
and in case m <= n, the elements on and below the (n-m)-th subdiagonal contain the m-by-n lower trapezoidal matrix L
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgeqlf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QL decomposition of the given dense double precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QL decomposition of a general m-by-n double precision complex column-major matrix based on the LAPACK sgeqlf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(m-k+i+1:m) = 0
and v(m-k+i) = 1
. v(1:m-k+i-1)
is stored on exit in A(1:m-k+i-1,n-k+i)
, and tau
in tau(i)
. Thus in case m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L
and in case m <= n, the elements on and below the (n-m)-th subdiagonal contain the m-by-n lower trapezoidal matrix L
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgeqlf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QL decomposition of the given dense matrix.
A | The matrix to be decomposed. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function performs the dense matrix QL decomposition of a general m-by-n matrix based on the LAPACK geqlf() functions. Note that this function can only be used for general, non-adapted matrices with float
, double
, complex<float>
, or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
In case of a column-major matrix, the resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(m-k+i+1:m) = 0
and v(m-k+i) = 1
. v(1:m-k+i-1)
is stored on exit in A(1:m-k+i-1,n-k+i)
, and tau
in tau(i)
. Thus in case m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L
and in case m <= n, the elements on and below the (n-m)-th subdiagonal contain the m-by-n lower trapezoidal matrix L
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
In case of a row-major matrix, the resulting decomposition is transposed, i.e. the elementary reflectors are stored in the lower part of the matrix and the elements in the upper part contain the n-by-min(m,n) upper trapezoidal matrix L
.
For more information on the geqlf() functions (i.e. sgeqlf(), dgeqlf(), cgeqlf(), and zgeqlf()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QR decomposition of the given dense single precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QR decomposition of a general m-by-n single precision column-major matrix based on the LAPACK sgeqrf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:m)
is stored on exit in A(i+1:m,i)
, and tau
in tau(i)
. Thus on exit the elements on and above the diagonal of the matrix 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 orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgeqrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QR decomposition of the given dense double precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QR decomposition of a general m-by-n double precision column-major matrix based on the LAPACK dgeqrf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:m)
is stored on exit in A(i+1:m,i)
, and tau
in tau(i)
. Thus on exit the elements on and above the diagonal of the matrix 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 orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the dgeqrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QR decomposition of the given dense single precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QR decomposition of a general m-by-n single precision complex column-major matrix based on the LAPACK cgeqrf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:m)
is stored on exit in A(i+1:m,i)
, and tau
in tau(i)
. Thus on exit the elements on and above the diagonal of the matrix 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 orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the cgeqrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QR decomposition of the given dense double precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix QR decomposition of a general m-by-n double precision complex column-major matrix based on the LAPACK zgeqrf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:m)
is stored on exit in A(i+1:m,i)
, and tau
in tau(i)
. Thus on exit the elements on and above the diagonal of the matrix 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 orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the zgeqrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the QR decomposition of the given dense matrix.
A | The matrix to be decomposed. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function performs the dense matrix QR decomposition of a general m-by-n matrix based on the LAPACK geqrf() functions. Note that this function can only be used for general, non-adapted matrices with float
, double
, complex<float>
, or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
In case of a column-major matrix, the resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(0:i-1) = 0
and v(i) = 1
. v(i+1:m)
is stored on exit in A(i+1:m,i)
, and tau
in tau(i)
. Thus on exit the elements on and above the diagonal of the matrix 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 orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
In case of a row-major matrix, the resulting decomposition is transposed, i.e. the elementary reflectors are stored above the diagonal and the elements on and below the diagonal contain the n-by-min(m,n) lower trapezoidal matrix R
.
For more information on the geqrf() functions (i.e. sgeqrf(), dgeqrf(), cgeqrf(), and zgeqrf()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the RQ decomposition of the given dense single precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix RQ decomposition of a general m-by-n single precision column-major matrix based on the LAPACK sgerqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(n-k+i+1:n) = 0
and v(n-k+i) = 1
. v(1:n-k+i-1)
is stored on exit in A(m-k+i,1:n-k+i-1)
, and tau
in tau(i)
. Thus in case m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the m-by-m upper triangular matrix R
and in case m >= n, the elements on and above the (m-n)-th subdiagonal contain the m-by-n upper trapezoidal matrix R
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the sgerqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the RQ decomposition of the given dense single precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix RQ decomposition of a general m-by-n double precision column-major matrix based on the LAPACK dgerqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(n-k+i+1:n) = 0
and v(n-k+i) = 1
. v(1:n-k+i-1)
is stored on exit in A(m-k+i,1:n-k+i-1)
, and tau
in tau(i)
. Thus in case m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the m-by-m upper triangular matrix R
and in case m >= n, the elements on and above the (m-n)-th subdiagonal contain the m-by-n upper trapezoidal matrix R
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the dgerqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the RQ decomposition of the given dense single precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix RQ decomposition of a general m-by-n single precision complex column-major matrix based on the LAPACK cgerqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(n-k+i+1:n) = 0
and v(n-k+i) = 1
. v(1:n-k+i-1)
is stored on exit in A(m-k+i,1:n-k+i-1)
, and tau
in tau(i)
. Thus in case m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the m-by-m upper triangular matrix R
and in case m >= n, the elements on and above the (m-n)-th subdiagonal contain the m-by-n upper trapezoidal matrix R
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the cgerqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the RQ decomposition of the given dense double precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix RQ decomposition of a general m-by-n double precision complex column-major matrix based on the LAPACK zgerqf() function. The resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(n-k+i+1:n) = 0
and v(n-k+i) = 1
. v(1:n-k+i-1)
is stored on exit in A(m-k+i,1:n-k+i-1)
, and tau
in tau(i)
. Thus in case m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the m-by-m upper triangular matrix R
and in case m >= n, the elements on and above the (m-n)-th subdiagonal contain the m-by-n upper trapezoidal matrix R
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
The info argument provides feedback on the success of the function call:
For more information on the zgerqf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the RQ decomposition of the given dense matrix.
A | The matrix to be decomposed. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function performs the dense matrix RQ decomposition of a general m-by-n matrix based on the LAPACK gerqf() functions. Note that this function can only be used for general, non-adapted matrices with float
, double
, complex<float>
, or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
In case of a column-major matrix, the resulting decomposition has the form
where the Q
is represented as a product of elementary reflectors
Each H(i) has the form
where tau
is a real scalar, and v
is a real vector with v(n-k+i+1:n) = 0
and v(n-k+i) = 1
. v(1:n-k+i-1)
is stored on exit in A(m-k+i,1:n-k+i-1)
, and tau
in tau(i)
. Thus in case m <= n, the upper triangle of the subarray A(1:m,n-m+1:n) contains the m-by-m upper triangular matrix R
and in case m >= n, the elements on and above the (m-n)-th subdiagonal contain the m-by-n upper trapezoidal matrix R
; the remaining elements in combination with the array tau
represent the orthogonal matrix Q
as a product of min(m,n) elementary reflectors.
In case of a row-major matrix, the resulting decomposition is transposed, i.e. the elementary reflectors are stored in the upper part of the matrix and the elements in the lower part contain the min(m,n)-by-m lower trapezoidal matrix R
.
For more information on the gerqf() functions (i.e. sgerqf(), dgerqf(), cgerqf(), and zgerqf()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LU decomposition of the given dense general single precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array for the pivot indices; size >= min( m, n ). |
info | Return code of the function call. |
This function performs the dense matrix LU decomposition of a general m-by-n single precision column-major matrix based on the LAPACK sgetrf() function, which uses partial pivoting with row interchanges. The resulting decomposition has the form
where P
is a permutation matrix, L
is a lower unitriangular matrix (lower trapezoidal if m > n), and U
is an upper triangular matrix (upper trapezoidal if m < n). The resulting decomposition is stored within the matrix A: L
is stored in the lower part of A and U
is stored in the upper part. The unit diagonal elements of L
are not stored.
The info argument provides feedback on the success of the function call:
For more information on the sgetrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LU decomposition of the given dense general double precision column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array for the pivot indices; size >= min( m, n ). |
info | Return code of the function call. |
This function performs the dense matrix LU decomposition of a general m-by-n double precision column-major matrix based on the LAPACK dgetrf() function, which uses partial pivoting with row interchanges. The resulting decomposition has the form
where P
is a permutation matrix, L
is a lower unitriangular matrix (lower trapezoidal if m > n), and U
is an upper triangular matrix (upper trapezoidal if m < n). The resulting decomposition is stored within the matrix A: L
is stored in the lower part of A and U
is stored in the upper part. The unit diagonal elements of L
are not stored.
The info argument provides feedback on the success of the function call:
For more information on the dgetrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LU decomposition of the given dense general single precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array for the pivot indices; size >= min( m, n ). |
info | Return code of the function call. |
This function performs the dense matrix LU decomposition of a general m-by-n single precision complex column-major matrix based on the LAPACK cgetrf() function, which uses partial pivoting with row interchanges. The resulting decomposition has the form
where P
is a permutation matrix, L
is a lower unitriangular matrix (lower trapezoidal if m > n), and U
is an upper triangular matrix (upper trapezoidal if m < n). The resulting decomposition is stored within the matrix A: L
is stored in the lower part of A and U
is stored in the upper part. The unit diagonal elements of L
are not stored.
The info argument provides feedback on the success of the function call:
For more information on the cgetrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the LU decomposition of the given dense general double precision complex column-major matrix.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array for the pivot indices; size >= min( m, n ). |
info | Return code of the function call. |
This function performs the dense matrix LU decomposition of a general m-by-n double precision complex column-major matrix based on the LAPACK zgetrf() function, which uses partial pivoting with row interchanges. The resulting decomposition has the form
where P
is a permutation matrix, L
is a lower unitriangular matrix (lower trapezoidal if m > n), and U
is an upper triangular matrix (upper trapezoidal if m < n). The resulting decomposition is stored within the matrix A: L
is stored in the lower part of A and U
is stored in the upper part. The unit diagonal elements of L
are not stored.
The info argument provides feedback on the success of the function call:
For more information on the zgetrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense Hermitian indefinite single precision complex column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the Hermitian matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix decomposition of a Hermitian indefinite single precision column-major matrix based on the LAPACK chetrf() function, which uses the Bunch-Kaufman diagonal pivoting method. The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
If the function exits successfully (i.e. info = 0) then the first element of the work array returns the optimal lwork. For optimal performance lwork >= n*NB, where NB is the optimal blocksize returned by the LAPACK function ilaenv(). If lwork = -1 then a workspace query is assumed. The function only calculates the optimal size of the work array and returns this value as the first entry of the work array.
For more information on the chetrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense Hermitian indefinite double precision complex column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the Hermitian matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix decomposition of a Hermitian indefinite double precision column-major matrix based on the LAPACK zhetrf() function, which uses the Bunch-Kaufman diagonal pivoting method. The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
If the function exits successfully (i.e. info = 0) then the first element of the work array returns the optimal lwork. For optimal performance lwork >= n*NB, where NB is the optimal blocksize returned by the LAPACK function ilaenv(). If lwork = -1 then a workspace query is assumed. The function only calculates the optimal size of the work array and returns this value as the first entry of the work array.
For more information on the zhetrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense Hermitian indefinite matrix.
A | The matrix to be decomposed. |
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
ipiv | Auxiliary array of size n for the pivot indices. |
std::invalid_argument | Invalid non-square matrix provided. |
std::invalid_argument | Invalid uplo argument provided. |
This function performs the dense matrix decomposition of a Hermitian indefinite matrix based on the LAPACK hetrf() functions, which use the Bunch-Kaufman diagonal pivoting method. Note that the function only works for general, non-adapted matrices with complex<float>
or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The function fails if ...
'L'
nor 'U'
.In all failure cases a std::invalid_argument exception is thrown.
For more information on the hetrf() functions (i.e. chetrf() and zhetrf()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a LQ decomposition based on the LAPACK sorglq() function for single precision column-major matrices that have already been factorized by the sgelqf() function. The info argument provides feedback on the success of the function call:
For more information on the sorglq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a LQ decomposition based on the LAPACK dorglq() function for double precision column-major matrices that have already been factorized by the dgelqf() function. The info argument provides feedback on the success of the function call:
For more information on the dorglq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a LQ decomposition based on the LAPACK orglq() functions from matrices that have already been LQ factorized by the gelqf() functions. Note that this function can only be used for general, non-adapted matrices with float
or double
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major m-by-min(m,n) or column-major min(m,n)-by-n Q matrix is stored in the within the given matrix A:
For more information on the orglq() functions (i.e. sorglq() and dorglq()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QL decomposition based on the LAPACK sorgql() function for single precision column-major matrices that have already been factorized by the sgeqlf() function. The info argument provides feedback on the success of the function call:
For more information on the sorgql() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QL decomposition based on the LAPACK dorgql() function for double precision column-major matrices that have already been factorized by the dgeqlf() function. The info argument provides feedback on the success of the function call:
For more information on the dorgql() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a QL decomposition based on the LAPACK orgql() functions from matrices that have already been QL factorized by the geqlf() functions. Note that this function can only be used for general, non-adapted matrices with float
or double
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major min(m,n)-by-n or column-major m-by-min(m,n) Q matrix is stored in the within the given matrix A:
For more information on the orgql() functions (i.e. sorgql() and dorgql()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QR decomposition based on the LAPACK sorgqr() function for single precision column-major matrices that have already been factorized by the sgeqrf() function. The info argument provides feedback on the success of the function call:
For more information on the sorgqr() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QR decomposition based on the LAPACK dorgqr() function for double precision column-major matrices that have already been factorized by the dgeqrf() function. The info argument provides feedback on the success of the function call:
For more information on the dorgqr() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a QR decomposition based on the LAPACK orgqr() functions from matrices that have already been QR factorized by the geqrf() functions. Note that this function can only be used for general, non-adapted matrices with float
or double
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major min(m,n)-by-n or column-major m-by-min(m,n) Q matrix is stored in the within the given matrix A:
For more information on the orgqr() functions (i.e. sorgqr() and dorgqr()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a RQ decomposition based on the LAPACK sorgrq() function for single precision column-major matrices that have already been factorized by the sgerqf() function. The info argument provides feedback on the success of the function call:
For more information on the sorgrq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a RQ decomposition based on the LAPACK dorgrq() function for double precision column-major matrices that have already been factorized by the dgerqf() function. The info argument provides feedback on the success of the function call:
For more information on the dorgrq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a RQ decomposition based on the LAPACK orgrq() functions from matrices that have already been RQ factorized by the gerqf() functions. Note that this function can only be used for general, non-adapted matrices with float
or double
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major m-by-min(m,n) or column-major min(m,n)-by-n Q matrix is stored in the within the given matrix A:
For more information on the orgrq() functions (i.e. sorgrq() and dorgrq()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite single precision column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
info | Return code of the function call. |
This function performs the dense matrix Cholesky decomposition of a symmetric positive definite single precision column-major matrix based on the LAPACK spotrf() function. The decomposition has the form
where U
is an upper triangular matrix and L
is a lower triangular matrix. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
For more information on the spotrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite double precision column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
info | Return code of the function call. |
This function performs the dense matrix Cholesky decomposition of a symmetric positive definite double precision column-major matrix based on the LAPACK dpotrf() function. The decomposition has the form
where U
is an upper triangular matrix and L
is a lower triangular matrix. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
For more information on the dpotrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite single precision complex column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
info | Return code of the function call. |
This function performs the dense matrix Cholesky decomposition of a symmetric positive definite single precision complex column-major matrix based on the LAPACK cpotrf() function. The decomposition has the form
where U
is an upper triangular matrix and L
is a lower triangular matrix. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
For more information on the cpotrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite double precision complex column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
info | Return code of the function call. |
This function performs the dense matrix Cholesky decomposition of a symmetric positive definite double precision complex column-major matrix based on the LAPACK zpotrf() function. The decomposition has the form
where U
is an upper triangular matrix and L
is a lower triangular matrix. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
For more information on the zpotrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the Cholesky decomposition of the given dense positive definite matrix.
A | The matrix to be decomposed. |
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
std::invalid_argument | Invalid non-square matrix provided. |
std::invalid_argument | Invalid uplo argument provided. |
std::invalid_argument | Decomposition of singular matrix failed. |
This function performs the dense matrix Cholesky decomposition of a symmetric positive definite matrix based on the LAPACK potrf() functions. Note that the function only works for general, non-adapted matrices with float
, double
, complex<float>
, or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The decomposition has the form
where U
is an upper triangular matrix and L
is a lower triangular matrix. The Cholesky decomposition fails if ...
'L'
nor 'U'
.In all failure cases a std::invalid_argument exception is thrown.
For more information on the potrf() functions (i.e. spotrf(), dpotrf(), cpotrf(), and zpotrf()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense symmetric indefinite single precision column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the symmetric matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix decomposition of a symmetric indefinite single precision column-major matrix based on the LAPACK ssytrf() function, which uses the Bunch-Kaufman diagonal pivoting method. The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
If the function exits successfully (i.e. info = 0) then the first element of the work array returns the optimal lwork. For optimal performance lwork >= n*NB, where NB is the optimal blocksize returned by the LAPACK function ilaenv(). If lwork = -1 then a workspace query is assumed. The function only calculates the optimal size of the work array and returns this value as the first entry of the work array.
For more information on the ssytrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense symmetric indefinite double precision column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the symmetric matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix decomposition of a symmetric indefinite double precision column-major matrix based on the LAPACK dsytrf() function, which uses the Bunch-Kaufman diagonal pivoting method. The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
If the function exits successfully (i.e. info = 0) then the first element of the work array returns the optimal lwork. For optimal performance lwork >= n*NB, where NB is the optimal blocksize returned by the LAPACK function ilaenv(). If lwork = -1 then a workspace query is assumed. The function only calculates the optimal size of the work array and returns this value as the first entry of the work array.
For more information on the dsytrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense symmetric indefinite single precision complex column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the symmetric matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix decomposition of a symmetric indefinite single precision column-major matrix based on the LAPACK csytrf() function, which uses the Bunch-Kaufman diagonal pivoting method. The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
If the function exits successfully (i.e. info = 0) then the first element of the work array returns the optimal lwork. For optimal performance lwork >= n*NB, where NB is the optimal blocksize returned by the LAPACK function ilaenv(). If lwork = -1 then a workspace query is assumed. The function only calculates the optimal size of the work array and returns this value as the first entry of the work array.
For more information on the csytrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense symmetric indefinite double precision complex column-major matrix.
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of the symmetric matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function performs the dense matrix decomposition of a symmetric indefinite double precision column-major matrix based on the LAPACK zsytrf() function, which uses the Bunch-Kaufman diagonal pivoting method. The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The info argument provides feedback on the success of the function call:
If the function exits successfully (i.e. info = 0) then the first element of the work array returns the optimal lwork. For optimal performance lwork >= n*NB, where NB is the optimal blocksize returned by the LAPACK function ilaenv(). If lwork = -1 then a workspace query is assumed. The function only calculates the optimal size of the work array and returns this value as the first entry of the work array.
For more information on the zsytrf() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the decomposition of the given dense symmetric indefinite matrix.
A | The matrix to be decomposed. |
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
ipiv | Auxiliary array of size n for the pivot indices. |
std::invalid_argument | Invalid non-square matrix provided. |
std::invalid_argument | Invalid uplo argument provided. |
This function performs the dense matrix decomposition of a symmetric indefinite matrix based on the LAPACK sytrf() functions, which use the Bunch-Kaufman diagonal pivoting method. Note that the function only works for general, non-adapted matrices with float
, double
, complex<float>
, or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The decomposition has the form
where U
(or L
) is a product of permutation and unit upper (lower) triangular matrices, and D
is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The resulting decomposition is stored within A: In case uplo is set to 'L'
the result is stored in the lower part of the matrix and the upper part remains untouched, in case uplo is set to 'U'
the result is stored in the upper part and the lower part remains untouched.
The function fails if ...
'L'
nor 'U'
.In all failure cases a std::invalid_argument exception is thrown.
For more information on the sytrf() functions (i.e. ssytrf(), dsytrf(), csytrf(), and zsytrf()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a LQ decomposition based on the LAPACK cunglq() function for single precision complex column-major matrices that have already been factorized by the cgelqf() function. The info argument provides feedback on the success of the function call:
For more information on the cunglq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a LQ decomposition based on the LAPACK zunglq() function for double precision complex column-major matrices that have already been factorized by the zgelqf() function. The info argument provides feedback on the success of the function call:
For more information on the zunglq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a LQ decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a LQ decomposition based on the LAPACK unglq() functions from matrices that have already been LQ factorized by the gelqf() functions. Note that this function can only be used for general, non-adapted matrices with complex<float>
or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major m-by-min(m,n) or column-major min(m,n)-by-n Q matrix is stored in the within the given matrix A:
For more information on the unglq() functions (i.e. cunglq() and zunglq()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QL decomposition based on the LAPACK cungql() function for single precision column-major matrices that have already been factorized by the sgeqlf() function. The info argument provides feedback on the success of the function call:
For more information on the cungql() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QL decomposition based on the LAPACK zungql() function for double precision column-major matrices that have already been factorized by the dgeqlf() function. The info argument provides feedback on the success of the function call:
For more information on the zungql() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QL decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a QL decomposition based on the LAPACK ungql() functions from matrices that have already been QL factorized by the geqlf() functions. Note that this function can only be used for general, non-adapted matrices with float
or double
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major m-by-min(m,n) or column-major min(m,n)-by-n Q matrix is stored in the within the given matrix A:
For more information on the ungql() functions (i.e. cungql() and zungql()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QR decomposition based on the LAPACK cungqr() function for single precision complex column-major matrices that have already been factorized by the cgeqrf() function. The info argument provides feedback on the success of the function call:
For more information on the cungqr() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision complex column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a QR decomposition based on the LAPACK zungqr() function for double precision complex column-major matrices that have already been factorized by the zgeqrf() function. The info argument provides feedback on the success of the function call:
For more information on the zungqr() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a QR decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a QR decomposition based on the LAPACK ungqr() functions from matrices that have already been QR factorized by the geqrf() functions. Note that this function can only be used for general, non-adapted matrices with complex<float>
or complex<double>
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major min(m,n)-by-n or column-major m-by-min(m,n) Q matrix is stored in the within the given matrix A:
For more information on the ungqr() functions (i.e. cungqr() and zungqr()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the single precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a RQ decomposition based on the LAPACK cungrq() function for single precision column-major matrices that have already been factorized by the sgerqf() function. The info argument provides feedback on the success of the function call:
For more information on the cungrq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
m | The number of rows of the given matrix ![]() |
n | The number of columns of the given matrix ![]() |
k | The number of elementary reflectors, whose product defines the matrix ![]() |
A | Pointer to the first element of the double precision column-major matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
work | Auxiliary array; size >= max( 1, lwork ). |
lwork | The dimension of the array work; size >= max( 1, n ). |
info | Return code of the function call. |
This function generates all or part of the orthogonal matrix Q from a RQ decomposition based on the LAPACK zungrq() function for double precision column-major matrices that have already been factorized by the dgerqf() function. The info argument provides feedback on the success of the function call:
For more information on the zungrq() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for the reconstruction of the orthogonal matrix Q from a RQ decomposition.
A | The decomposed matrix. |
tau | Array for the scalar factors of the elementary reflectors; size >= min( m, n ). |
This function reconstructs the orthogonal matrix Q of a RQ decomposition based on the LAPACK ungrq() functions from matrices that have already been RQ factorized by the gerqf() functions. Note that this function can only be used for general, non-adapted matrices with float
or double
element type. The attempt to call the function with any adapted matrix or matrices of any other element type results in a compile time error!
The row-major m-by-min(m,n) or column-major min(m,n)-by-n Q matrix is stored in the within the given matrix A:
For more information on the ungrq() functions (i.e. cungrq() and zungrq()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/