Mentions légales du service

Skip to content
Snippets Groups Projects
Commit 390cf7c8 authored by PRUVOST Florent's avatar PRUVOST Florent Committed by Mathieu Faverge
Browse files

Blas/Lapack API: add lapack routines for Cholesky family functions

parent 83e00175
No related branches found
No related tags found
No related merge requests found
Showing
with 583 additions and 19 deletions
......@@ -863,10 +863,11 @@
* *Map functions*
* map: apply a user operator on each tile of the matrix
In addition, all *BLAS 3 routines* (gemm, hemm, her2k, herk, lauum, symm,
syr2k, syrk, trmm, trsm) can be called using an equivalent of the
(C)BLAS/LAPACK(E) API. The parameters are the same and the user just has
to add *CHAMELEON_* to the standard name of the routine. For example, in C
In addition, all *BLAS 3 routines* ~gemm~, ~hemm~, ~her2k~, ~herk~, ~lauum~, ~symm~,
~syr2k~, ~syrk~, ~trmm~, ~trsm~ and *LAPACK* ~lacpy~, ~laset~, ~posv~, ~potrf~, ~potri~,
~potrs~, ~trtri~ can be called using an equivalent of the (C)BLAS/LAPACK(E)
API. The parameters are the same and the user just has to add *CHAMELEON_*
to the standard name of the routine. For example, in C
#+begin_src
CHAMELEON_Init(4,0);
CHAMELEON_cblas_dgemm(CblasColMajor, CblasNoTrans, CblasNoTrans,
......@@ -874,7 +875,7 @@
CHAMELEON_Finalize();
#+end_src
In Fortran, the function names are for example: ~CHAMELEON_blas_dgemm~
instead of ~DGEMM~ and ~CHAMELEON_lapack_dlauum~ instead of ~DLAUUM~.
instead of ~DGEMM~ and ~CHAMELEON_lapack_dposv~ instead of ~DPOSV~.
**** Options routines
<<sec:options_routines>>
......
......@@ -63,11 +63,16 @@ set(ZSRC
src/lapack_zlacpy.c
src/lapack_zlaset.c
src/lapack_zlauum.c
src/lapack_zposv.c
src/lapack_zpotrf.c
src/lapack_zpotri.c
src/lapack_zpotrs.c
src/lapack_zsymm.c
src/lapack_zsyr2k.c
src/lapack_zsyrk.c
src/lapack_ztrmm.c
src/lapack_ztrsm.c
src/lapack_ztrtri.c
)
precisions_rules_py(LAPACK_SRCS_GENERATED "${ZSRC}"
PRECISIONS "${CHAMELEON_PRECISION}")
......
......@@ -79,7 +79,6 @@ void CHAMELEON_cblas_ztrsm( const CBLAS_ORDER order, const CBLAS_SIDE side, cons
const void *alpha, const CHAMELEON_Complex64_t *A, const int lda,
const CHAMELEON_Complex64_t *B, const int ldb );
int CHAMELEON_lapacke_zlacpy( int matrix_layout, char uplo, int M, int N,
const CHAMELEON_Complex64_t *A, int lda,
CHAMELEON_Complex64_t *B, int ldb );
......@@ -91,6 +90,23 @@ int CHAMELEON_lapacke_zlaset( int matrix_layout, char uplo, int M, int N,
int CHAMELEON_lapacke_zlauum( int matrix_layout, char uplo, int N,
CHAMELEON_Complex64_t *A, int lda );
int CHAMELEON_lapacke_zposv( int matrix_layout, char uplo, int N, int NRHS,
CHAMELEON_Complex64_t *A, int lda,
CHAMELEON_Complex64_t *B, int ldb );
int CHAMELEON_lapacke_zpotrf( int matrix_layout, char uplo, int N,
CHAMELEON_Complex64_t *A, int lda );
int CHAMELEON_lapacke_zpotri( int matrix_layout, char uplo, int N,
CHAMELEON_Complex64_t *A, int lda );
int CHAMELEON_lapacke_zpotrs( int matrix_layout, char uplo, int N, int NRHS,
const CHAMELEON_Complex64_t *A, int lda,
CHAMELEON_Complex64_t *B, int ldb );
int CHAMELEON_lapacke_ztrtri( int matrix_layout, char uplo, char diag, int N,
CHAMELEON_Complex64_t *A, int lda );
END_C_DECLS
#endif /* _chameleon_zlapack_h_ */
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_zhemm CHAMELEON_GLOBAL( chameleon_blas_zhemm, CHAMELEON_BLAS_Zhemm )
#define CHAMELEON_blas_zhemm CHAMELEON_GLOBAL( chameleon_blas_zhemm, CHAMELEON_BLAS_ZHEMM )
void CHAMELEON_blas_zhemm ( const char* side, const char* uplo,
const int* m, const int* n,
const CHAMELEON_Complex64_t* alpha, const CHAMELEON_Complex64_t* a, const int* lda,
......
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_zher2k CHAMELEON_GLOBAL( chameleon_blas_zher2k, CHAMELEON_BLAS_Zher2k )
#define CHAMELEON_blas_zher2k CHAMELEON_GLOBAL( chameleon_blas_zher2k, CHAMELEON_BLAS_ZHER2K )
void CHAMELEON_blas_zher2k ( const char* uplo, const char* trans,
const int* n, const int* k,
const CHAMELEON_Complex64_t* alpha, const CHAMELEON_Complex64_t* a, const int* lda,
......
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_zherk CHAMELEON_GLOBAL( chameleon_blas_zherk, CHAMELEON_BLAS_Zherk )
#define CHAMELEON_blas_zherk CHAMELEON_GLOBAL( chameleon_blas_zherk, CHAMELEON_BLAS_ZHERK )
void CHAMELEON_blas_zherk ( const char* uplo, const char* trans,
const int* n, const int* k,
const double* alpha, const CHAMELEON_Complex64_t* a, const int* lda,
......
......@@ -22,7 +22,7 @@
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zlacpy CHAMELEON_GLOBAL( chameleon_lapack_zlacpy, CHAMELEON_BLAS_Zlacpy )
#define CHAMELEON_lapack_zlacpy CHAMELEON_GLOBAL( chameleon_lapack_zlacpy, CHAMELEON_LAPACK_ZLACPY )
void CHAMELEON_lapack_zlacpy ( const char* uplo, const int* m, const int* n,
const CHAMELEON_Complex64_t* a, const int* lda,
CHAMELEON_Complex64_t* b, const int* ldb )
......
......@@ -22,7 +22,7 @@
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zlaset CHAMELEON_GLOBAL( chameleon_lapack_zlaset, CHAMELEON_BLAS_Zlaset )
#define CHAMELEON_lapack_zlaset CHAMELEON_GLOBAL( chameleon_lapack_zlaset, CHAMELEON_LAPACK_ZLASET )
void CHAMELEON_lapack_zlaset ( const char* uplo, const int* m, const int* n,
const CHAMELEON_Complex64_t* alpha, const CHAMELEON_Complex64_t* beta,
CHAMELEON_Complex64_t* a, const int* lda )
......
......@@ -22,7 +22,7 @@
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zlauum CHAMELEON_GLOBAL( chameleon_lapack_zlauum, CHAMELEON_BLAS_Zlauum )
#define CHAMELEON_lapack_zlauum CHAMELEON_GLOBAL( chameleon_lapack_zlauum, CHAMELEON_LAPACK_ZLAUUM )
void CHAMELEON_lapack_zlauum ( const char* uplo, const int* n,
CHAMELEON_Complex64_t* a, const int* lda,
int* info )
......
/**
*
* @file lapack_zposv.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for posv
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-08-19
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zposv CHAMELEON_GLOBAL( chameleon_lapack_zposv, CHAMELEON_LAPACK_ZPOSV )
void CHAMELEON_lapack_zposv ( const char* uplo, const int* n, const int* nrhs,
CHAMELEON_Complex64_t* a, const int* lda,
CHAMELEON_Complex64_t* b, const int* ldb,
int* info )
{
*info = CHAMELEON_lapacke_zposv( CblasColMajor,
*uplo, *n, *nrhs, a, *lda, b, *ldb );
}
/* C LAPACKE interface */
/**
********************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zposv - Computes the solution to a system of linear equations A * X = B,
* where A is an N-by-N symmetric positive definite (or Hermitian positive definite
* in the complex case) matrix and X and B are N-by-NRHS matrices.
* The Cholesky decomposition is used to factor A as
*
* \f[ A = \{_{L\times L^H, if uplo = ChamLower}^{U^H\times U, if uplo = ChamUpper} \f]
*
* where U is an upper triangular matrix and L is a lower triangular matrix.
* The factored form of A is then used to solve the system of equations A * X = B.
*
*******************************************************************************
*
* @param[in] matrix_layout
* Specifies whether the matrices are row or column major, it must be
* set to LAPACK_COL_MAJOR or CblasColMajor (102), the matrix_layout
* supported in Chameleon.
*
* @param[in] uplo
* Specifies whether the matrix A is upper triangular or lower triangular:
* = 'U' or 'u': Upper triangle of A is stored;
* = 'L' or 'l': Lower triangle of A is stored.
*
* @param[in] N
* The number of linear equations, i.e., the order of the matrix A. N >= 0.
*
* @param[in] NRHS
* The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
*
* @param[in,out] A
* On entry, the symmetric positive definite (or Hermitian) matrix A.
* If UPLO = 'U', the leading N-by-N upper triangular part of A
* contains the upper triangular part of the matrix A, and the
* strictly lower triangular part of A is not referenced.
* If UPLO = 'L', the leading N-by-N lower triangular part of A
* contains the lower triangular part of the matrix A, and the
* strictly upper triangular part of A is not referenced.
* On exit, if return value = 0, the factor U or L from the Cholesky factorization
* A = U^H*U or A = L*L^H.
*
* @param[in] LDA
* The leading dimension of the array A. LDA >= max(1,N).
*
* @param[in,out] B
* On entry, the N-by-NRHS right hand side matrix B.
* On exit, if return value = 0, the N-by-NRHS solution matrix X.
*
* @param[in] LDB
* The leading dimension of the array B. LDB >= max(1,N).
*
*******************************************************************************
*
* @retval =0 successful exit
* @retval <0 if -i, the i-th argument had an illegal value
* @retval >0 if i, the leading minor of order i of A is not
positive definite, so the factorization could not be
completed, and the solution has not been computed.
*
*******************************************************************************
*
* @sa CHAMELEON_lapacke_zposv
* @sa CHAMELEON_lapacke_cposv
* @sa CHAMELEON_lapacke_dposv
* @sa CHAMELEON_lapacke_sposv
*
*/
int CHAMELEON_lapacke_zposv( int matrix_layout, char uplo, int N, int NRHS,
CHAMELEON_Complex64_t *A, int lda,
CHAMELEON_Complex64_t *B, int ldb )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zposv", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zposv( (cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
N, NRHS, A, lda, B, ldb );
}
/**
*
* @file lapack_zpotrf.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for potrf
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-08-19
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zpotrf CHAMELEON_GLOBAL( chameleon_lapack_zpotrf, CHAMELEON_LAPACK_ZPOTRF )
void CHAMELEON_lapack_zpotrf ( const char* uplo, const int* n,
CHAMELEON_Complex64_t* a, const int* lda,
int* info )
{
*info = CHAMELEON_lapacke_zpotrf( CblasColMajor, *uplo, *n, a, *lda );
}
/* C LAPACKE interface */
/**
********************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zpotrf - Computes the Cholesky factorization of a symmetric positive definite
* (or Hermitian positive definite in the complex case) matrix A.
* The factorization has the form
*
* \f[ A = \{_{L\times L^H, if uplo = ChamLower}^{U^H\times U, if uplo = ChamUpper} \f]
*
* where U is an upper triangular matrix and L is a lower triangular matrix.
*
*******************************************************************************
*
* @param[in] matrix_layout
* Specifies whether the matrices are row or column major, it must be
* set to LAPACK_COL_MAJOR or CblasColMajor (102), the matrix_layout
* supported in Chameleon.
*
* @param[in] uplo
* Specifies whether the matrix A is upper triangular or lower triangular:
* = 'U' or 'u': Upper triangle of A is stored;
* = 'L' or 'l': Lower triangle of A is stored.
*
* @param[in] N
* The order of the matrix A. N >= 0.
*
* @param[in,out] A
* On entry, the symmetric positive definite (or Hermitian) matrix A.
* If uplo = 'U', the leading N-by-N upper triangular part of A
* contains the upper triangular part of the matrix A, and the strictly lower triangular
* part of A is not referenced.
* If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower
* triangular part of the matrix A, and the strictly upper triangular part of A is not
* referenced.
* On exit, if return value = 0, the factor U or L from the Cholesky factorization
* A = U^H*U or A = L*L^H.
*
* @param[in] LDA
* The leading dimension of the array A. LDA >= max(1,N).
*
*******************************************************************************
*
* @retval =0 successful exit
* @retval <0 if -i, the i-th argument had an illegal value
* @retval >0 if i, the leading minor of order i of A is not positive definite, so the
* factorization could not be completed, and the solution has not been computed.
*
*******************************************************************************
*
* @sa CHAMELEON_lapacke_zpotrf
* @sa CHAMELEON_lapacke_cpotrf
* @sa CHAMELEON_lapacke_dpotrf
* @sa CHAMELEON_lapacke_spotrf
*
*/
int CHAMELEON_lapacke_zpotrf( int matrix_layout, char uplo, int N,
CHAMELEON_Complex64_t *A, int lda )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zpotrf", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zpotrf( (cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
N, A, lda );
}
/**
*
* @file lapack_zpotri.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for potri
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-08-19
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zpotri CHAMELEON_GLOBAL( chameleon_lapack_zpotri, CHAMELEON_LAPACK_ZPOTRI )
void CHAMELEON_lapack_zpotri ( const char* uplo, const int* n,
CHAMELEON_Complex64_t* a, const int* lda,
int* info )
{
*info = CHAMELEON_lapacke_zpotri( CblasColMajor, *uplo, *n, a, *lda );
}
/* C LAPACKE interface */
/**
********************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zpotri - Computes the inverse of a symmetric (Hermitian)
* positive-definite matrix using the Cholesky factorization. Before calling
* this routine, call potrf to factorize.
*
*******************************************************************************
*
* @param[in] matrix_layout
* Specifies whether the matrices are row or column major, it must be
* set to LAPACK_COL_MAJOR or CblasColMajor (102), the matrix_layout
* supported in Chameleon.
*
* @param[in] uplo
* Specifies whether the matrix A is upper triangular or lower triangular:
* = 'U' or 'u': Upper triangle of A is stored;
* = 'L' or 'l': Lower triangle of A is stored.
*
* @param[in] N
* The order of the matrix A. N >= 0.
*
* @param[in,out] A
* On entry, the symmetric positive definite (or Hermitian) matrix A.
* If uplo = 'U', the leading N-by-N upper triangular part of A
* contains the upper triangular part of the matrix A, and the strictly lower triangular
* part of A is not referenced.
* If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower
* triangular part of the matrix A, and the strictly upper triangular part of A is not
* referenced.
* On exit, if return value = 0, the factor U or L from the Cholesky factorization
* A = U^H*U or A = L*L^H.
*
* @param[in] LDA
* The leading dimension of the array A. LDA >= max(1,N).
*
*******************************************************************************
*
* @retval =0 successful exit
* @retval <0 if -i, the i-th argument had an illegal value
* @retval >0 if i, the leading minor of order i of A is not positive definite, so the
* factorization could not be completed, and the solution has not been computed.
*
*******************************************************************************
*
* @sa CHAMELEON_lapacke_zpotri
* @sa CHAMELEON_lapacke_cpotri
* @sa CHAMELEON_lapacke_dpotri
* @sa CHAMELEON_lapacke_spotri
*
*/
int CHAMELEON_lapacke_zpotri( int matrix_layout, char uplo, int N,
CHAMELEON_Complex64_t *A, int lda )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zpotri", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zpotri( (cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
N, A, lda );
}
/**
*
* @file lapack_zpotrs.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for potrs
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-08-19
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zpotrs CHAMELEON_GLOBAL( chameleon_lapack_zpotrs, CHAMELEON_LAPACK_ZPOTRS )
void CHAMELEON_lapack_zpotrs ( const char* uplo, const int* n, const int* nrhs,
const CHAMELEON_Complex64_t* a, const int* lda,
CHAMELEON_Complex64_t* b, const int* ldb,
int* info )
{
*info = CHAMELEON_lapacke_zpotrs( CblasColMajor,
*uplo, *n, *nrhs, a, *lda, b, *ldb );
}
/* C LAPACKE interface */
/**
********************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zpotrs - Solves a system of linear equations A * X = B with a symmetric positive
* definite (or Hermitian positive definite in the complex case) matrix A using the Cholesky
* factorization A = U^H*U or A = L*L^H computed by CHAMELEON_zpotrf.
*
*******************************************************************************
*
* @param[in] matrix_layout
* Specifies whether the matrices are row or column major, it must be
* set to LAPACK_COL_MAJOR or CblasColMajor (102), the matrix_layout
* supported in Chameleon.
*
* @param[in] uplo
* Specifies whether the matrix A is upper triangular or lower triangular:
* = 'U' or 'u': Upper triangle of A is stored;
* = 'L' or 'l': Lower triangle of A is stored.
*
* @param[in] N
* The order of the matrix A. N >= 0.
*
* @param[in] NRHS
* The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0.
*
* @param[in] A
* The triangular factor U or L from the Cholesky factorization A = U^H*U or A = L*L^H,
* computed by CHAMELEON_lapacke_zpotrf.
*
* @param[in] LDA
* The leading dimension of the array A. LDA >= max(1,N).
*
* @param[in,out] B
* On entry, the N-by-NRHS right hand side matrix B.
* On exit, if return value = 0, the N-by-NRHS solution matrix X.
*
* @param[in] LDB
* The leading dimension of the array B. LDB >= max(1,N).
*
*******************************************************************************
*
* @retval =0 successful exit
* @retval <0 if -i, the i-th argument had an illegal value
*
*******************************************************************************
*
* @sa CHAMELEON_lapacke_zpotrs
* @sa CHAMELEON_lapacke_cpotrs
* @sa CHAMELEON_lapacke_dpotrs
* @sa CHAMELEON_lapacke_spotrs
*
*/
int CHAMELEON_lapacke_zpotrs( int matrix_layout, char uplo, int N, int NRHS,
const CHAMELEON_Complex64_t *A, int lda,
CHAMELEON_Complex64_t *B, int ldb )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zpotrs", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zpotrs( (cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
N, NRHS, (CHAMELEON_Complex64_t *)A, lda, B, ldb );
}
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_zsymm CHAMELEON_GLOBAL( chameleon_blas_zsymm, CHAMELEON_BLAS_Zsymm )
#define CHAMELEON_blas_zsymm CHAMELEON_GLOBAL( chameleon_blas_zsymm, CHAMELEON_BLAS_ZSYMM )
void CHAMELEON_blas_zsymm ( const char* side, const char* uplo,
const int* m, const int* n,
const CHAMELEON_Complex64_t* alpha, const CHAMELEON_Complex64_t* a, const int* lda,
......
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_zsyr2k CHAMELEON_GLOBAL( chameleon_blas_zsyr2k, CHAMELEON_BLAS_Zsyr2k )
#define CHAMELEON_blas_zsyr2k CHAMELEON_GLOBAL( chameleon_blas_zsyr2k, CHAMELEON_BLAS_ZSYR2K )
void CHAMELEON_blas_zsyr2k ( const char* uplo, const char* trans,
const int* n, const int* k,
const CHAMELEON_Complex64_t* alpha, const CHAMELEON_Complex64_t* a, const int* lda,
......
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_zsyrk CHAMELEON_GLOBAL( chameleon_blas_zsyrk, CHAMELEON_BLAS_Zsyrk )
#define CHAMELEON_blas_zsyrk CHAMELEON_GLOBAL( chameleon_blas_zsyrk, CHAMELEON_BLAS_ZSYRK )
void CHAMELEON_blas_zsyrk ( const char* uplo, const char* trans,
const int* n, const int* k,
const CHAMELEON_Complex64_t* alpha, const CHAMELEON_Complex64_t* a, const int* lda,
......
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_ztrmm CHAMELEON_GLOBAL( chameleon_blas_ztrmm, CHAMELEON_BLAS_Ztrmm )
#define CHAMELEON_blas_ztrmm CHAMELEON_GLOBAL( chameleon_blas_ztrmm, CHAMELEON_BLAS_ZTRMM )
void CHAMELEON_blas_ztrmm ( const char* side, const char* uplo,
const char* trans, const char* diag,
const int* m, const int* n,
......
......@@ -22,7 +22,7 @@
/* Fortran BLAS interface */
#define CHAMELEON_blas_ztrsm CHAMELEON_GLOBAL( chameleon_blas_ztrsm, CHAMELEON_BLAS_Ztrsm )
#define CHAMELEON_blas_ztrsm CHAMELEON_GLOBAL( chameleon_blas_ztrsm, CHAMELEON_BLAS_ZTRSM )
void CHAMELEON_blas_ztrsm ( const char* side, const char* uplo,
const char* trans, const char* diag,
const int* m, const int* n,
......
/**
*
* @file lapack_ztrtri.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for trtri
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-08-19
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_ztrtri CHAMELEON_GLOBAL( chameleon_lapack_ztrtri, CHAMELEON_LAPACK_ZTRTRI )
void CHAMELEON_lapack_ztrtri ( const char* uplo, const char* diag, const int* n,
CHAMELEON_Complex64_t* a, const int* lda,
int* info )
{
*info = CHAMELEON_lapacke_ztrtri( CblasColMajor,
*uplo, *diag, *n, a, *lda );
}
/* C LAPACKE interface */
/**
********************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_ztrtri - Computes the inverse of a complex upper or lower
* triangular matrix A.
*
*******************************************************************************
*
* @param[in] matrix_layout
* Specifies whether the matrices are row or column major, it must be
* set to LAPACK_COL_MAJOR or CblasColMajor (102), the matrix_layout
* supported in Chameleon.
*
* @param[in] uplo
* = 'U' or 'u': Upper triangle of A is stored;
* = 'L' or 'l': Lower triangle of A is stored.
*
* @param[in] diag
* = 'N', then A is not a unit triangular matrix;
* = 'U', A is unit triangular: diagonal elements of A are assumed to
* be 1 and not referenced in the array a.
*
* @param[in] N
* The order of the matrix A. N >= 0.
*
* @param[in,out] A
* On entry, the triangular matrix A. If UPLO = 'U', the
* leading N-by-N upper triangular part of the array A
* contains the upper triangular matrix, and the strictly
* lower triangular part of A is not referenced. If UPLO =
* 'L', the leading N-by-N lower triangular part of the array
* A contains the lower triangular matrix, and the strictly
* upper triangular part of A is not referenced. If DIAG =
* 'U', the diagonal elements of A are also not referenced and
* are assumed to be 1. On exit, the (triangular) inverse of
* the original matrix, in the same storage format.
*
* @param[in] LDA
* The leading dimension of the array A. LDA >= max(1,N).
*
*******************************************************************************
*
* @retval =0 successful exit
* @retval <0 if -i, the i-th argument had an illegal value
* @retval >0 if i, A(i,i) is exactly zero. The triangular
* matrix is singular and its inverse can not be computed.
*
*******************************************************************************
*
* @sa CHAMELEON_lapacke_ztrtri
* @sa CHAMELEON_lapacke_ctrtri
* @sa CHAMELEON_lapacke_dtrtri
* @sa CHAMELEON_lapacke_strtri
*
*/
int CHAMELEON_lapacke_ztrtri( int matrix_layout, char uplo, char diag, int N,
CHAMELEON_Complex64_t *A, int lda )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_ztrtri", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_ztrtri( (cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
(cham_diag_t)chameleon_blastocblas_diag(&diag),
N, A, lda );
}
......@@ -22,9 +22,11 @@
#include "testings.h"
#include "testing_zcheck.h"
#include <chameleon/flops.h>
#if !defined(CHAMELEON_SIMULATION)
#include <coreblas.h>
#if defined(CHAMELEON_TESTINGS_VENDOR)
#include <coreblas/lapacke.h>
#include <coreblas.h>
#endif
#endif
static cham_fixdbl_t
......@@ -119,6 +121,7 @@ testing_zposv_std( run_arg_list_t *args, int check )
int hres = 0;
/* Read arguments */
int api = parameters_getvalue_int( "api" );
int nb = run_arg_get_int( args, "nb", 320 );
cham_uplo_t uplo = run_arg_get_uplo( args, "uplo", ChamUpper );
int N = run_arg_get_int( args, "N", 1000 );
......@@ -149,7 +152,22 @@ testing_zposv_std( run_arg_list_t *args, int check )
testing_stop( &test_data, flops_zposv( N, NRHS ) );
#else
testing_start( &test_data );
hres = CHAMELEON_zposv( uplo, N, NRHS, A, LDA, X, LDB );
switch ( api ) {
case 1:
hres = CHAMELEON_zposv( uplo, N, NRHS, A, LDA, X, LDB );
break;
#if !defined(CHAMELEON_SIMULATION)
case 2:
CHAMELEON_lapacke_zposv( CblasColMajor, chameleon_lapack_const(uplo), N, NRHS, A, LDA, X, LDB );
break;
#endif
default:
if ( CHAMELEON_Comm_rank() == 0 ) {
fprintf( stderr,
"SKIPPED: This function can only be used with the option --api 1 or --api 2.\n" );
}
return -1;
}
test_data.hres = hres;
testing_stop( &test_data, flops_zposv( N, NRHS ) );
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment