LAMA
lama::OpenMPLAPACK Class Reference

#include <OpenMPLAPACK.hpp>

Public Member Functions

template<>
IndexType getrf (const enum CBLAS_ORDER order, const int m, const int n, float *const A, const int lda, int *const ipiv)
template<>
IndexType getrf (const enum CBLAS_ORDER order, const int m, const int n, double *const A, const int lda, int *const ipiv)
template<>
void getinv (const IndexType n, float *a, const IndexType lda)
template<>
void getinv (const IndexType n, double *a, const IndexType lda)
template<>
int getri (const enum CBLAS_ORDER UNUSED(order), const int UNUSED(n), float *const UNUSED(a), const int UNUSED(lda), int *const UNUSED(ipiv))
template<>
int getri (const enum CBLAS_ORDER UNUSED(order), const int UNUSED(n), double *const UNUSED(a), const int UNUSED(lda), int *const UNUSED(ipiv))
template<>
int trtrs (const enum CBLAS_ORDER order, const enum CBLAS_UPLO uplo, const enum CBLAS_TRANSPOSE trans, const enum CBLAS_DIAG diag, const int n, const int nrhs, const float *A, const int lda, float *B, const int ldb)
template<>
int trtrs (const enum CBLAS_ORDER order, const enum CBLAS_UPLO uplo, const enum CBLAS_TRANSPOSE trans, const enum CBLAS_DIAG diag, const int n, const int nrhs, const double *A, const int lda, double *B, const int ldb)
template<>
int tptrs (const enum CBLAS_ORDER order, const enum CBLAS_UPLO uplo, const enum CBLAS_TRANSPOSE trans, const enum CBLAS_DIAG diag, const int n, const int nrhs, const float *AP, float *B, const int ldb)
template<>
int tptrs (const enum CBLAS_ORDER order, const enum CBLAS_UPLO uplo, const enum CBLAS_TRANSPOSE trans, const enum CBLAS_DIAG diag, const int n, const int nrhs, const double *AP, double *B, const int ldb)
template<>
void laswp (const enum CBLAS_ORDER order, const int N, float *A, const int LDA, const int K1, const int K2, const int *ipiv, const int INCX, SyncToken *syncToken)
template<>
void laswp (const enum CBLAS_ORDER order, const int N, double *A, const int LDA, const int K1, const int K2, const int *ipiv, const int INCX, SyncToken *syncToken)

Static Public Member Functions

template<typename T >
static IndexType getrf (const enum CBLAS_ORDER order, const IndexType m, const IndexType n, T *const a, const IndexType lda, IndexType *const ipiv)
 getrf_cpu computes the LU factorization of a general m-by-n matrix A in floating point single precision as A = P*L*U, where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n) and U is upper triangular (upper trapezoidal if m < n).
template<typename T >
static IndexType getri (const enum CBLAS_ORDER order, const IndexType n, T *const A, const IndexType lda, IndexType *const ipiv)
template<typename T >
static void getinv (const IndexType n, T *a, const IndexType lda)
 Implementation of LAPACKInterface::getinv vi LAPACK.
template<typename T >
static IndexType trtrs (const enum CBLAS_ORDER order, const enum CBLAS_UPLO uplo, const enum CBLAS_TRANSPOSE trans, const enum CBLAS_DIAG diag, const IndexType n, const IndexType nrhs, const T *A, const IndexType lda, T *B, const IndexType ldb)
 trtrs solves the following equation system: op(A)*X = B where op(A) is either A, AT or AH; and B is a matrix of right hand sides and will contain the solution of all equations on output.
template<typename T >
static IndexType tptrs (const enum CBLAS_ORDER order, const enum CBLAS_UPLO uplo, const enum CBLAS_TRANSPOSE trans, const enum CBLAS_DIAG diag, const IndexType n, const IndexType nrhs, const T *AP, T *B, const IndexType ldb)
 tptrs solves the following equation system: op(A)*X = B where op(A) is either A, AT or AH; and B is a matrix of right hand sides and will contain the solution of all equations on output.
template<typename T >
static void laswp (const enum CBLAS_ORDER order, const IndexType n, T *A, const IndexType lda, const IndexType k1, const IndexType k2, const IndexType *ipiv, const IndexType incx, SyncToken *syncToken)
 laswp_cpu performs a series of row interchanges on the matrix A.

Private Member Functions

 LAMA_LOG_DECL_STATIC_LOGGER (logger)

Member Function Documentation

template<typename T >
static void lama::OpenMPLAPACK::getinv ( const IndexType  n,
T *  a,
const IndexType  lda 
) [static]

Implementation of LAPACKInterface::getinv vi LAPACK.

template<>
void lama::OpenMPLAPACK::getinv ( const IndexType  n,
float *  a,
const IndexType  lda 
)
template<>
void lama::OpenMPLAPACK::getinv ( const IndexType  n,
double *  a,
const IndexType  lda 
)
template<typename T >
static IndexType lama::OpenMPLAPACK::getrf ( const enum CBLAS_ORDER  order,
const IndexType  m,
const IndexType  n,
T *const  a,
const IndexType  lda,
IndexType *const  ipiv 
) [static]

getrf_cpu computes the LU factorization of a general m-by-n matrix A in floating point single precision as A = P*L*U, where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n) and U is upper triangular (upper trapezoidal if m < n).

The routine uses partial pivoting, with row interchanges. L and U will be stored within A whereby the diagonal elements of L will not be stored.

Parameters:
[in]orderSpecifies, whether the matrix is stored in column major order (i.e. CblasColMajor) or in row major order (i.e. CblasRowMajor).
[in]mNumber of rows of matrix A; m must be at least zero. m specifies, how many rows of A will be touched by this function.
[in]nNumber of columns of matrix A; n must be at least zero. n specifies, how many columns of A will be touched by this function.
[in,out]aOn input, a is holding the array, containing the matrix A. On output it will be overwritten by L and U, whereby the diagonal elements of L will not be stored.
[in]ldaThe first dimension of array a. lda specifies the actual number of rows of A. lda is used to compute the position of the next column, i.e. if position (r,c) is going to be accessed its position will be computed by (c*lda + r).
[in,out]ipivThe array holding the permutation of the matrix. Its size must be at least min(m,n). It contains the info that row i has been changed with ipiv[i]. ipiv is assumed to be stored in C-Style, i.e the values, representing the indexes of the matrix are assumed to be starting with zero and ending with m-1. It will also leave with this assumption, but in between, it is first incremented and afterwards decremented, to fit the Fortran interface.
Returns:
info If info=0, the execution is successful. If info = -i, the i-th parameter had an illegal value. If info = i, uii is 0. The factorization has been completed, but U is exactly singular. Division by 0 will occur if you use the factor U for solving a system of linear equations.

Referenced by lama::LUSolver::computeLUFactorization().

template<>
IndexType lama::OpenMPLAPACK::getrf ( const enum CBLAS_ORDER  order,
const int  m,
const int  n,
float *const  A,
const int  lda,
int *const  ipiv 
)
template<>
IndexType lama::OpenMPLAPACK::getrf ( const enum CBLAS_ORDER  order,
const int  m,
const int  n,
double *const  A,
const int  lda,
int *const  ipiv 
)
template<typename T >
static IndexType lama::OpenMPLAPACK::getri ( const enum CBLAS_ORDER  order,
const IndexType  n,
T *const  A,
const IndexType  lda,
IndexType *const  ipiv 
) [static]
template<>
int lama::OpenMPLAPACK::getri ( const enum CBLAS_ORDER   UNUSEDorder,
const int   UNUSEDn,
float *const   UNUSEDa,
const int   UNUSEDlda,
int *const   UNUSEDipiv 
)

References LAMA_THROWEXCEPTION.

template<>
int lama::OpenMPLAPACK::getri ( const enum CBLAS_ORDER   UNUSEDorder,
const int   UNUSEDn,
double *const   UNUSEDa,
const int   UNUSEDlda,
int *const   UNUSEDipiv 
)

References LAMA_THROWEXCEPTION.

template<typename T >
static void lama::OpenMPLAPACK::laswp ( const enum CBLAS_ORDER  order,
const IndexType  n,
T *  A,
const IndexType  lda,
const IndexType  k1,
const IndexType  k2,
const IndexType ipiv,
const IndexType  incx,
SyncToken syncToken 
) [static]

laswp_cpu performs a series of row interchanges on the matrix A.

One row interchange is initiated for each of rows k1 through k2 of A.

Parameters:
[in]orderSpecifies, whether the matrix is stored in column major order (i.e. CblasColMajor) or in row major order (i.e. CblasRowMajor). Since a translation of the data would be too expensiv, if it was stored in row major order, the BLAS level1 function SSWAP will be called instead. The beginning column of the vector in A will then be LDA-N.
[in]nThe number of columns of the matrix A.
[in,out]AArray of dimension (LDA,N). On entry, the matrix of column dimension N to which the row interchanges will be applied. On exit, the permuted matrix.
[in]ldaIf the matrix is stored in column major order, lda specifies the actual number of rows of A. If else the matrix is stored in row major order, lda specifies the actual number of columns of A.
[in]k1The first element of ipiv for which a row interchange will be done.
[in]k2The last element of ipiv for which a row interchange will be done.
[in]ipivArray of dimension (k2*abs(incx)). The vector of pivot indices. Only the elements in positions k1 through k2 of ipiv are accessed. ipiv(k) = l implies rows k and l are to be interchanged.
[in]incxThe increment between successive values of ipiv. If ipiv is negative, the pivots are applied in reverse order.

Referenced by lama::LUSolver::computeLUFactorization().

template<>
void lama::OpenMPLAPACK::laswp ( const enum CBLAS_ORDER  order,
const int  N,
float *  A,
const int  LDA,
const int  K1,
const int  K2,
const int *  ipiv,
const int  INCX,
SyncToken syncToken 
)
template<>
void lama::OpenMPLAPACK::laswp ( const enum CBLAS_ORDER  order,
const int  N,
double *  A,
const int  LDA,
const int  K1,
const int  K2,
const int *  ipiv,
const int  INCX,
SyncToken syncToken 
)
template<typename T >
static IndexType lama::OpenMPLAPACK::tptrs ( const enum CBLAS_ORDER  order,
const enum CBLAS_UPLO  uplo,
const enum CBLAS_TRANSPOSE  trans,
const enum CBLAS_DIAG  diag,
const IndexType  n,
const IndexType  nrhs,
const T *  AP,
T *  B,
const IndexType  ldb 
) [static]

tptrs solves the following equation system: op(A)*X = B where op(A) is either A, AT or AH; and B is a matrix of right hand sides and will contain the solution of all equations on output.

Parameters:
[in]orderSpecifies, whether the matrix is stored in column major order (i.e. CblasColMajor) or in row major order (i.e. CblasRowMajor).
[in]uploSpecifies, whether matrix A is upper triangular (i.e. CblasUpper) or lower triangular (i.e. CblasLower).
[in]transSpecifies op(A). if trans == CblasNoTrans, op(A) = A; if trans == CblasTrans, op(A) = AT; if trans == CblasConjTrans, op(A) = AH;
[in]diagSpecifies, whether the triangualr matrix A is a unit triangular matrix, i.e. the diagonal elements of A are one. if diag == CblasNonUnit, the diagonal elements of A are not assumed to be one; if diag == CblasUnit, the diagonal elements of A are assumed to be one and therefor not referenced;
[in]nThe number of columns of A and rows of B; n must be at least 0.
[in]nrhsThe number of columns of B; nrhs must be at least 0.
[in]APThe array containing matrix A.
[in,out]BOn input B is the array holding the right hand sides of the equations, on output, it will hold the solution for each equation system.
[in]ldbThe first dimension of B, i.e. the actual number of columns of B.
template<>
int lama::OpenMPLAPACK::tptrs ( const enum CBLAS_ORDER  order,
const enum CBLAS_UPLO  uplo,
const enum CBLAS_TRANSPOSE  trans,
const enum CBLAS_DIAG  diag,
const int  n,
const int  nrhs,
const float *  AP,
float *  B,
const int  ldb 
)
template<>
int lama::OpenMPLAPACK::tptrs ( const enum CBLAS_ORDER  order,
const enum CBLAS_UPLO  uplo,
const enum CBLAS_TRANSPOSE  trans,
const enum CBLAS_DIAG  diag,
const int  n,
const int  nrhs,
const double *  AP,
double *  B,
const int  ldb 
)
template<typename T >
static IndexType lama::OpenMPLAPACK::trtrs ( const enum CBLAS_ORDER  order,
const enum CBLAS_UPLO  uplo,
const enum CBLAS_TRANSPOSE  trans,
const enum CBLAS_DIAG  diag,
const IndexType  n,
const IndexType  nrhs,
const T *  A,
const IndexType  lda,
T *  B,
const IndexType  ldb 
) [static]

trtrs solves the following equation system: op(A)*X = B where op(A) is either A, AT or AH; and B is a matrix of right hand sides and will contain the solution of all equations on output.

Parameters:
[in]orderSpecifies, whether the matrix is stored in column major order (i.e. CblasColMajor) or in row major order (i.e. CblasRowMajor).
[in]uploSpecifies, whether matrix A is upper triangular (i.e. CblasUpper) or lower triangular (i.e. CblasLower).
[in]transSpecifies op(A). if trans == CblasNoTrans, op(A) = A; if trans == CblasTrans, op(A) = AT; if trans == CblasConjTrans, op(A) = AH;
[in]diagSpecifies, whether the triangualr matrix A is a unit triangular matrix, i.e. the diagonal elements of A are one. if diag == CblasNonUnit, the diagonal elements of A are not assumed to be one; if diag == CblasUnit, the diagonal elements of A are assumed to be one and therefor not referenced;
[in]nThe number of columns of A and rows of B; n must be at least 0.
[in]nrhsThe number of columns of B; nrhs must be at least 0.
[in]AThe array containing matrix A.
[in]ldaThe first dimension of A, i.e. the actual number of rows of A.
[in,out]BOn input B is the array holding the right hand sides of the equations, on output, it will hold the solution for each equation system.
[in]ldbThe first dimension of B, i.e. the actual number of columns of B.

Referenced by lama::LUSolver::solve().

template<>
int lama::OpenMPLAPACK::trtrs ( const enum CBLAS_ORDER  order,
const enum CBLAS_UPLO  uplo,
const enum CBLAS_TRANSPOSE  trans,
const enum CBLAS_DIAG  diag,
const int  n,
const int  nrhs,
const float *  A,
const int  lda,
float *  B,
const int  ldb 
)
template<>
int lama::OpenMPLAPACK::trtrs ( const enum CBLAS_ORDER  order,
const enum CBLAS_UPLO  uplo,
const enum CBLAS_TRANSPOSE  trans,
const enum CBLAS_DIAG  diag,
const int  n,
const int  nrhs,
const double *  A,
const int  lda,
double *  B,
const int  ldb 
)

The documentation for this class was generated from the following files: