![]() |
LAPACK general linear system functions (gesv) | |
void | blaze::gesv (int n, int nrhs, float *A, int lda, int *ipiv, float *B, int ldb, int *info) |
LAPACK kernel for solving a general single precision linear system of equations ( ![]() | |
void | blaze::gesv (int n, int nrhs, double *A, int lda, int *ipiv, double *B, int ldb, int *info) |
LAPACK kernel for solving a general double precision linear system of equations ( ![]() | |
void | blaze::gesv (int n, int nrhs, complex< float > *A, int lda, int *ipiv, complex< float > *B, int ldb, int *info) |
LAPACK kernel for solving a general single precision complex linear system of equations ( ![]() | |
void | blaze::gesv (int n, int nrhs, complex< double > *A, int lda, int *ipiv, complex< double > *B, int ldb, int *info) |
LAPACK kernel for solving a general double precision complex linear system of equations ( ![]() | |
template<typename MT , bool SO, typename VT , bool TF> | |
void | blaze::gesv (DenseMatrix< MT, SO > &A, DenseVector< VT, TF > &b, int *ipiv) |
LAPACK kernel for solving a general linear system of equations ( ![]() | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
void | blaze::gesv (DenseMatrix< MT1, SO1 > &A, DenseMatrix< MT2, SO2 > &B, int *ipiv) |
LAPACK kernel for solving a general linear system of equations ( ![]() | |
LAPACK Hermitian indefinite linear system functions (hesv) | |
void | blaze::hesv (char uplo, int n, int nrhs, complex< float > *A, int lda, int *ipiv, complex< float > *B, int ldb, complex< float > *work, int lwork, int *info) |
LAPACK kernel for solving a Hermitian indefinite single precision complex linear system of equations ( ![]() | |
void | blaze::hesv (char uplo, int n, int nrhs, complex< double > *A, int lda, int *ipiv, complex< double > *B, int ldb, complex< double > *work, int lwork, int *info) |
LAPACK kernel for solving a Hermitian indefinite double precision complex linear system of equations ( ![]() | |
template<typename MT , bool SO, typename VT , bool TF> | |
void | blaze::hesv (DenseMatrix< MT, SO > &A, DenseVector< VT, TF > &b, char uplo, int *ipiv) |
LAPACK kernel for solving a Hermitian indefinite linear system of equations ( ![]() | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
void | blaze::hesv (DenseMatrix< MT1, SO1 > &A, DenseMatrix< MT2, SO2 > &B, char uplo, int *ipiv) |
LAPACK kernel for solving a Hermitian indefinite linear system of equations ( ![]() | |
LAPACK positive definite linear system functions (posv) | |
void | blaze::posv (char uplo, int n, int nrhs, float *A, int lda, float *B, int ldb, int *info) |
LAPACK kernel for solving a positive definite single precision linear system of equations ( ![]() | |
void | blaze::posv (char uplo, int n, int nrhs, double *A, int lda, double *B, int ldb, int *info) |
LAPACK kernel for solving a positive definite double precision linear system of equations ( ![]() | |
void | blaze::posv (char uplo, int n, int nrhs, complex< float > *A, int lda, complex< float > *B, int ldb, int *info) |
LAPACK kernel for solving a positive definite single precision complex linear system of equations ( ![]() | |
void | blaze::posv (char uplo, int n, int nrhs, complex< double > *A, int lda, complex< double > *B, int ldb, int *info) |
LAPACK kernel for solving a positive definite double precision complex linear system of equations ( ![]() | |
template<typename MT , bool SO, typename VT , bool TF> | |
void | blaze::posv (DenseMatrix< MT, SO > &A, DenseVector< VT, TF > &b, char uplo) |
LAPACK kernel for solving a positive definite linear system of equations ( ![]() | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
void | blaze::posv (DenseMatrix< MT1, SO1 > &A, DenseMatrix< MT2, SO2 > &B, char uplo) |
LAPACK kernel for solving a positive definite linear system of equations ( ![]() | |
LAPACK symmetric indefinite linear system functions (sysv) | |
void | blaze::sysv (char uplo, int n, int nrhs, float *A, int lda, int *ipiv, float *B, int ldb, float *work, int lwork, int *info) |
LAPACK kernel for solving a symmetric indefinite single precision linear system of equations ( ![]() | |
void | blaze::sysv (char uplo, int n, int nrhs, double *A, int lda, int *ipiv, double *B, int ldb, double *work, int lwork, int *info) |
LAPACK kernel for solving a symmetric indefinite double precision linear system of equations ( ![]() | |
void | blaze::sysv (char uplo, int n, int nrhs, complex< float > *A, int lda, int *ipiv, complex< float > *B, int ldb, complex< float > *work, int lwork, int *info) |
LAPACK kernel for solving a symmetric indefinite single precision complex linear system of equations ( ![]() | |
void | blaze::sysv (char uplo, int n, int nrhs, complex< double > *A, int lda, int *ipiv, complex< double > *B, int ldb, complex< double > *work, int lwork, int *info) |
LAPACK kernel for solving a symmetric indefinite double precision complex linear system of equations ( ![]() | |
template<typename MT , bool SO, typename VT , bool TF> | |
void | blaze::sysv (DenseMatrix< MT, SO > &A, DenseVector< VT, TF > &b, char uplo, int *ipiv) |
LAPACK kernel for solving a symmetric indefinite linear system of equations ( ![]() | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
void | blaze::sysv (DenseMatrix< MT1, SO1 > &A, DenseMatrix< MT2, SO2 > &B, char uplo, int *ipiv) |
LAPACK kernel for solving a symmetric indefinite linear system of equations ( ![]() | |
LAPACK triangular linear system functions (trsv) | |
void | blaze::trsv (char uplo, char trans, char diag, int n, const float *A, int lda, float *x, int incX) |
LAPACK kernel for solving a triangular single precision linear system of equations ( ![]() | |
void | blaze::trsv (char uplo, char trans, char diag, int n, const double *A, int lda, double *x, int incX) |
LAPACK kernel for solving a triangular double precision linear system of equations ( ![]() | |
void | blaze::trsv (char uplo, char trans, char diag, int n, const complex< float > *A, int lda, complex< float > *x, int incX) |
LAPACK kernel for solving a triangular single precision complex linear system of equations ( ![]() | |
void | blaze::trsv (char uplo, char trans, char diag, int n, const complex< double > *A, int lda, complex< double > *x, int incX) |
LAPACK kernel for solving a triangular double precision complex linear system of equations ( ![]() | |
template<typename MT , bool SO, typename VT , bool TF> | |
void | blaze::trsv (const DenseMatrix< MT, SO > &A, DenseVector< VT, TF > &b, char uplo, char trans, char diag) |
LAPACK kernel for solving a triangular linear system of equations ( ![]() ![]() | |
|
inline |
LAPACK kernel for solving a general linear system of equations ( ).
A | The column-major system matrix. |
b | The right-hand side vector. |
ipiv | Auxiliary array of size n for the pivot indices. |
std::invalid_argument | Invalid non-square matrix provided. |
std::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK gesv() functions to compute the solution to the system of general linear equations:
In this context the general system matrix A n-by-n matrix and x and b are n-dimensional vectors. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the vector b contains the solution of the linear system of equations and A has been decomposed by means of an LU decomposition with partial pivoting and row interchanges. The decomposition has the form
where P
is a permutation matrix, L
is a lower unitriangular matrix, and U
is an upper triangular matrix. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
In all failure cases an exception is thrown.
For more information on the gesv() functions (i.e. sgesv(), dgesv(), cgesv(), and zgesv()), see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a general linear system of equations ( ).
A | The column-major system matrix. |
B | The matrix of right-hand sides. |
ipiv | Auxiliary array of size n for the pivot indices. |
std::invalid_argument | Invalid non-square matrix provided. |
std::invalid_argument | Matrix sizes do not match. |
std::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK gesv() functions to compute the solution to the general system of linera equations:
In this context the general system matrix A is a n-by-n matrix and X and B are either row-major m-by-n matrices or column-major n-by-m matrices. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the matrix B contains the solutions of the linear system of equations and A has been decomposed by means of an LU decomposition with partial pivoting and row interchanges. The decomposition has the form
where P
is a permutation matrix, L
is a lower unitriangular matrix, and U
is an upper triangular matrix. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
In all failure cases an exception is thrown.
For more information on the gesv() functions (i.e. sgesv(), dgesv(), cgesv(), and zgesv()), see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a general single precision linear system of equations ( ).
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the single precision column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK sgesv() function to compute the solution to the general system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The LU decomposition with partial pivoting and row interchanges is used to factor A as
where P
is a permutation matrix, L
is a lower unitriangular matrix, and U
is an upper triangular matrix. The resulting decomposition is stored within 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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the sgesv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a general double precision linear system of equations ( ).
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the double precision column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK dgesv() function to compute the solution to the general system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The LU decomposition with partial pivoting and row interchanges is used to factor A as
where P
is a permutation matrix, L
is a lower unitriangular matrix, and U
is an upper triangular matrix. The resulting decomposition is stored within 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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the dgesv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a general single precision complex linear system of equations ( ).
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the single precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK cgesv() function to compute the solution to the general system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The LU decomposition with partial pivoting and row interchanges is used to factor A as
where P
is a permutation matrix, L
is a lower unitriangular matrix, and U
is an upper triangular matrix. The resulting decomposition is stored within 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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the cgesv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a general double precision complex linear system of equations ( ).
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the double precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK zgesv() function to compute the solution to the general system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The LU decomposition with partial pivoting and row interchanges is used to factor A as
where P
is a permutation matrix, L
is a lower unitriangular matrix, and U
is an upper triangular matrix. The resulting decomposition is stored within 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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the zgesv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a Hermitian indefinite linear system of equations ( ).
A | The column-major system matrix. |
b | The right-hand side vector. |
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. |
std::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK hesv() functions to compute the solution to the Hermitian indefinite system of linear equations;
In this context the Hermitian indefinite system matrix A is a n-by-n matrix and x and b are n-dimensional vectors. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the vector b contains the solution of the linear system of equations and A has been decomposed by means of the Bunch-Kaufman decomposition. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
'L'
nor 'U'
;In all failure cases an exception is thrown.
For more information on the hesv() functions (i.e. chesv() and zhesv()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a Hermitian indefinite linear system of equations ( ).
A | The system matrix. |
B | The matrix of right-hand sides. |
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. |
std::invalid_argument | Matrix sizes do not match. |
std::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK hesv() functions to compute the solution to the Hermitian indefinite system of linear equations:
In this context the Hermitian indefinite system matrix A is a n-by-n matrix and X and / B are either row-major m-by-n matrices or column-major n-by-m matrices. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the matrix B contains the solution of the linear system of equations and A has been decomposed by means of a Bunch-Kaufman decomposition. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
'L'
nor 'U'
;In all failure cases an exception is thrown.
For more information on the hesv() functions (i.e. chesv() and zhesv()) see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a Hermitian indefinite single precision complex linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the single precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
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 uses the LAPACK chesv() function to compute the solution to the Hermitian indefinite system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The Bunch-Kaufman decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the chesv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a Hermitian indefinite double precision complex linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the double precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
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 uses the LAPACK zhesv() function to compute the solution to the Hermitian indefinite system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The Bunch-Kaufman decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the zhesv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a positive definite linear system of equations ( ).
A | The column-major system matrix. |
b | The right-hand side vector. |
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::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK posv() functions to compute the solution to the positive definite system of linear equations:
In this context the positive definite system matrix A is a n-by-n matrix and x and b are n-dimensional vectors. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the vector b contains the solution of the linear system of equations and A has been decomposed by means of a Cholesky decomposition. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
'L'
nor 'U'
;In all failure cases an exception is thrown.
For more information on the posv() functions (i.e. sposv(), dposv(), cposv(), and zposv()), see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a positive definite linear system of equations ( ).
A | The system matrix. |
B | The matrix of right-hand sides. |
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 | Matrix sizes do not match. |
std::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK posv() functions to compute the solution to the positive definite system of linear equations:
In this context the positive definite system matrix A is a n-by-n matrix and X and B are either row-major m-by-n matrices or column-major n-by-m matrices. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the matrix B contains the solution of the linear system of equations and A has been decomposed by means of a Cholesky decomposition. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
'L'
nor 'U'
;In all failure cases an exception is thrown.
For more information on the posv() functions (i.e. sposv(), dposv(), cposv(), and zposv()), see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a positive definite single precision linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the single precision column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK sposv() function to compute the solution to the positive definite system of linear equations , where A is a n-by-n positive definite matrix and X and B are n-by-nrhs matrices.
The Cholesky decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the sposv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a positive definite double precision linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the double precision column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK dposv() function to compute the solution to the positive definite system of linear equations , where A is a n-by-n positive definite matrix and X and B are n-by-nrhs matrices.
The Cholesky decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the dposv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a positive definite single precision complex linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the single precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK cposv() function to compute the solution to the positive definite system of linear equations , where A is a n-by-n positive definite matrix and X and B are n-by-nrhs matrices.
The Cholesky decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the cposv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a positive definite double precision complex linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the double precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
info | Return code of the function call. |
This function uses the LAPACK zposv() function to compute the solution to the positive definite system of linear equations , where A is a n-by-n positive definite matrix and X and B are n-by-nrhs matrices.
The Cholesky decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the zposv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a symmetric indefinite linear system of equations ( ).
A | The column-major system matrix. |
b | The right-hand side vector. |
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. |
std::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK sysv() functions to compute the solution to the symmetric indefinite system of linear equations;
In this context the symmetric indefinite system matrix A is a n-by-n matrix and x and b are n-dimensional vectors. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the vector b contains the solution of the linear system of equations and A has been decomposed by means of the Bunch-Kaufman decomposition. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
'L'
nor 'U'
;In all failure cases an exception is thrown.
For more information on the sysv() functions (i.e. ssysv(), dsysv(), csysv(), and zsysv()), see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a symmetric indefinite linear system of equations ( ).
A | The system matrix. |
B | The matrix of right-hand sides. |
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. |
std::invalid_argument | Matrix sizes do not match. |
std::runtime_error | Inversion of singular matrix failed. |
This function uses the LAPACK sysv() functions to compute the solution to the symmetric indefinite system of linear equations:
In this context the symmetric indefinite system matrix A is a n-by-n matrix and X and / B are either row-major m-by-n matrices or column-major n-by-m matrices. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the matrix B contains the solution of the linear system of equations and A has been decomposed by means of a Bunch-Kaufman decomposition. 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 factored form of A is then used to solve the system of equations.
The function fails if ...
'L'
nor 'U'
;In all failure cases an exception is thrown.
For more information on the sysv() functions (i.e. ssysv(), dsysv(), csysv(), and zsysv()), see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a symmetric indefinite single precision linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the single precision column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
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 uses the LAPACK ssysv() function to compute the solution to the symmetric indefinite system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The Bunch-Kaufman decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the ssysv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a symmetric indefinite double precision linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the double precision column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
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 uses the LAPACK dsysv() function to compute the solution to the symmetric indefinite system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The Bunch-Kaufman decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the dsysv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a symmetric indefinite single precision complex linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the single precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
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 uses the LAPACK csysv() function to compute the solution to the symmetric indefinite system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The Bunch-Kaufman decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the csysv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a symmetric indefinite double precision complex linear system of equations ( ).
uplo | 'L' to use the lower part of the matrix, 'U' to use the upper part. |
n | The number of rows/columns of matrix A ![]() |
nrhs | The number of right-hand side vectors ![]() |
A | Pointer to the first element of the double precision complex column-major square matrix. |
lda | The total number of elements between two columns of matrix A ![]() |
ipiv | Auxiliary array of size n for the pivot indices. |
B | Pointer to the first element of the column-major matrix. |
ldb | The total number of elements between two columns of matrix B ![]() |
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 uses the LAPACK zsysv() function to compute the solution to the symmetric indefinite system of linear equations , where A is a n-by-n matrix and X and B are n-by-nrhs matrices.
The Bunch-Kaufman decomposition is used to factor A as
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 factored form of A is then used to solve the system of equations.
The info argument provides feedback on the success of the function call:
For more information on the zsysv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a triangular linear system of equations ( ). (
).
A | The system matrix. |
b | The right-hand side vector. |
uplo | 'L' in case of a lower matrix, 'U' in case of an upper matrix. |
trans | 'N' for ![]() 'T' for ![]() C for ![]() |
diag | 'U' in case of a unitriangular matrix, 'N' otherwise. |
std::invalid_argument | Invalid non-square matrix provided. |
std::invalid_argument | Invalid uplo argument provided. |
std::invalid_argument | Invalid trans argument provided. |
std::invalid_argument | Invalid diag argument provided. |
This function uses the LAPACK trsv() functions to compute the solution to the triangular system of linear equations:
In this context the positive definite system matrix A is a n-by-n matrix and x and b are n-dimensional vectors. 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 adaptors or matrices of any other element type results in a compile time error!
If the function exits successfully, the vector x contains the solution of the linear system of equations. The function fails if ...
'L'
nor 'U'
;'N'
nor 'T'
nor 'C'
;'U'
nor 'N'
.In all failure cases a std::invalid_argument exception is thrown.
Examples:
For more information on the trsv() functions (i.e. strsv(), dtrsv(), ctrsv(), and ztrsv()), see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a triangular single precision linear system of equations ( ).
uplo | 'L' in case of a lower matrix, 'U' in case of an upper matrix. |
trans | 'N' for ![]() 'T' for ![]() C for ![]() |
diag | 'U' in case of a unitriangular matrix, 'N' otherwise. |
n | The number of rows/columns of the column-major triangular matrix ![]() |
A | Pointer to the first element of the single precision column-major square matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
x | Pointer to the first element of vector x. |
incX | The stride within vector x. |
This function uses the LAPACK strsv() function to compute the solution to the triangular system of linear equations , where A is a n-by-n triangular matrix and x and b are n-dimensional vectors. The trans argument specifies the form of the linear system of equations:
For more information on the strsv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a triangular double precision linear system of equations ( ).
uplo | 'L' in case of a lower matrix, 'U' in case of an upper matrix. |
trans | 'N' for ![]() 'T' for ![]() C for ![]() |
diag | 'U' in case of a unitriangular matrix, 'N' otherwise. |
n | The number of rows/columns of the column-major triangular matrix ![]() |
A | Pointer to the first element of the double precision column-major square matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
x | Pointer to the first element of vector x. |
incX | The stride within vector x. |
This function uses the LAPACK dtrsv() function to compute the solution to the triangular system of linear equations , where A is a n-by-n triangular matrix and x and b are n-dimensional vectors. The trans argument specifies the form of the linear system of equations:
For more information on the dtrsv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a triangular single precision complex linear system of equations ( ).
uplo | 'L' in case of a lower matrix, 'U' in case of an upper matrix. |
trans | 'N' for ![]() 'T' for ![]() C for ![]() |
diag | 'U' in case of a unitriangular matrix, 'N' otherwise. |
n | The number of rows/columns of the column-major triangular matrix ![]() |
A | Pointer to the first element of the single precision complex column-major square matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
x | Pointer to the first element of vector x. |
incX | The stride within vector x. |
This function uses the LAPACK ctrsv() function to compute the solution to the triangular system of linear equations , where A is a n-by-n triangular matrix and x and b are n-dimensional vectors. The trans argument specifies the form of the linear system of equations:
For more information on the ctrsv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/
|
inline |
LAPACK kernel for solving a triangular double precision complex linear system of equations ( ).
uplo | 'L' in case of a lower matrix, 'U' in case of an upper matrix. |
trans | 'N' for ![]() 'T' for ![]() C for ![]() |
diag | 'U' in case of a unitriangular matrix, 'N' otherwise. |
n | The number of rows/columns of the column-major triangular matrix ![]() |
A | Pointer to the first element of the double precision complex column-major square matrix. |
lda | The total number of elements between two columns of the matrix ![]() |
x | Pointer to the first element of vector x. |
incX | The stride within vector x. |
This function uses the LAPACK ztrsv() function to compute the solution to the triangular system of linear equations , where A is a n-by-n triangular matrix and x and b are n-dimensional vectors. The trans argument specifies the form of the linear system of equations:
For more information on the ztrsv() function, see the LAPACK online documentation browser:
http://www.netlib.org/lapack/explore-html/