Mentions légales du service

Skip to content
Snippets Groups Projects
Commit c2702956 authored by PRUVOST Florent's avatar PRUVOST Florent
Browse files

Merge branch 'lapack_api_norm' into 'master'

Blas/Lapack API: add lapack routines for norms computation (lange, lanhe, lansy, lantr)

See merge request solverstack/chameleon!349
parents d507acf4 a652915e
No related branches found
No related tags found
No related merge requests found
......@@ -864,10 +864,10 @@
* 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~ 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
~syr2k~, ~syrk~, ~trmm~, ~trsm~ and *LAPACK* ~lacpy~, ~lange~, ~lanhe~, ~lansy~,
~lantr~, ~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,
......
......@@ -61,6 +61,10 @@ set(ZSRC
src/lapack_zher2k.c
src/lapack_zherk.c
src/lapack_zlacpy.c
src/lapack_zlange.c
src/lapack_zlanhe.c
src/lapack_zlansy.c
src/lapack_zlantr.c
src/lapack_zlaset.c
src/lapack_zlauum.c
src/lapack_zposv.c
......
......@@ -83,6 +83,18 @@ 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 );
double CHAMELEON_lapacke_zlange( int matrix_layout, char norm, int M, int N,
const CHAMELEON_Complex64_t *A, int lda );
double CHAMELEON_lapacke_zlanhe( int matrix_layout, char norm, char uplo, int N,
const CHAMELEON_Complex64_t *A, int lda );
double CHAMELEON_lapacke_zlansy( int matrix_layout, char norm, char uplo, int N,
const CHAMELEON_Complex64_t *A, int lda );
double CHAMELEON_lapacke_zlantr( int matrix_layout, char norm, char uplo, char diag,
int M, int N, const CHAMELEON_Complex64_t *A, int lda );
int CHAMELEON_lapacke_zlaset( int matrix_layout, char uplo, int M, int N,
const CHAMELEON_Complex64_t alpha, const CHAMELEON_Complex64_t beta,
CHAMELEON_Complex64_t *A, int lda );
......
......@@ -34,5 +34,6 @@ int chameleon_blastocblas_trans(const char* trans);
int chameleon_blastocblas_side(const char* side);
int chameleon_blastocblas_uplo(const char* uplo);
int chameleon_blastocblas_diag(const char* diag);
int chameleon_lapacktochameleon_norm(const char* norm);
#endif /* _lapack_api_common_h_ */
......@@ -91,3 +91,25 @@ int chameleon_blastocblas_diag(const char* diag)
return CHAMELEON_ERR_ILLEGAL_VALUE;
}
}
/**
* @brief Convert the input char LAPACK norm parameter to a compatible parameter
* for the Chameleon API.
* @param[in] norm The input char LAPACK norm parameter
* @return The Chameleon equivalent parameter (ChamMaxNorm, ChamOneNorm, etc).
*/
int chameleon_lapacktochameleon_norm(const char* norm)
{
if ( (*norm == 'M') || (*norm == 'm') ) {
return ChamMaxNorm;
} else if ( (*norm == '1') || (*norm == 'O') || (*norm == 'o') ) {
return ChamOneNorm;
} else if ( (*norm == 'I') || (*norm == 'i') ) {
return ChamInfNorm;
} else if ( (*norm == 'F') || (*norm == 'f') || (*norm == 'E') || (*norm == 'e') ) {
return ChamFrobeniusNorm;
} else {
fprintf(stderr, "CHAMELEON ERROR: %s(): %s\n", "chameleon_lapacktochameleon_norm", "illegal value of LAPACK norm parameter");
return CHAMELEON_ERR_ILLEGAL_VALUE;
}
}
\ No newline at end of file
/**
*
* @file lapack_zlange.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for lange
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-12-21
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zlange CHAMELEON_GLOBAL( chameleon_lapack_zlange, CHAMELEON_LAPACK_ZLANGE )
double CHAMELEON_lapack_zlange ( const char* norm, const int* m, const int* n,
const CHAMELEON_Complex64_t* a, const int* lda,
double* work )
{
(void)work;
return CHAMELEON_lapacke_zlange( CblasColMajor,
*norm, *m, *n, a, *lda );
}
/* C LAPACKE interface */
/**
******************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zlange - returns the value of the 1-norm, or the Frobenius
* norm, or the infinity norm, or the element of largest absolute value of a
* real/complex 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] norm = 'M' or 'm': val = max(abs(Aij)), largest absolute value
* of the matrix A.
* = '1' or 'O' or 'o': val = norm1(A), 1-norm of the matrix A
* (maximum column sum),
* = 'I' or 'i': val = normI(A), infinity norm of the matrix A
* (maximum row sum),
* = 'F', 'f', 'E' or 'e': val = normF(A), Frobenius norm of
* the matrix A (square root of sum of squares).
*
* @param[in] M The number of rows of the matrix A. M >= 0. When M = 0, the
* returned value is set to zero.
*
* @param[in] N The number of columns of the matrix A. N >= 0. When N = 0, the
* returned value is set to zero.
*
* @param[in] A The M-by-N matrix A.
*
* @param[in] LDA The leading dimension of the array A. LDA >= max(1,M).
*
******************************************************************************
*
* @retval the norm described above.
*
******************************************************************************
*
* @sa CHAMELEON_lapacke_zlange
* @sa CHAMELEON_lapacke_clange
* @sa CHAMELEON_lapacke_dlange
* @sa CHAMELEON_lapacke_slange
*
*/
double CHAMELEON_lapacke_zlange( int matrix_layout, char norm, int M, int N,
const CHAMELEON_Complex64_t *A, int lda )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zlange", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zlange( (cham_normtype_t)chameleon_lapacktochameleon_norm(&norm),
M, N, (CHAMELEON_Complex64_t *)A, lda );
}
/**
*
* @file lapack_zlanhe.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for lanhe
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-12-21
* @precisions normal z -> c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zlanhe CHAMELEON_GLOBAL( chameleon_lapack_zlanhe, CHAMELEON_LAPACK_ZLANHE )
double CHAMELEON_lapack_zlanhe ( const char* norm, const char* uplo, const int* n,
const CHAMELEON_Complex64_t* a, const int* lda,
double* work )
{
(void)work;
return CHAMELEON_lapacke_zlanhe( CblasColMajor,
*norm, *uplo, *n, a, *lda );
}
/* C LAPACKE interface */
/**
******************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zlanhe - returns the value of the one norm, or the
Frobenius norm, or the infinity norm, or the element of largest absolute
value of a complex hermitian 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] norm = max(abs(A(i,j))), NORM = 'M' or 'm'
* = norm1(A), NORM = '1', 'O' or 'o'
* = normI(A), NORM = 'I' or 'i'
* = normF(A), NORM = 'F', 'f', 'E' or 'e'
*
* where norm1 denotes the one norm of a matrix (maximum column sum), normI
* denotes the infinity norm of a matrix (maximum row sum) and normF denotes
* the Frobenius norm of a matrix (square root of sum of squares). Note that
* max(abs(A(i,j))) is not a consistent matrix norm.
*
* @param[in] uplo Specifies whether the upper or lower triangular part of the
* hermitian matrix A is to be referenced.
* = 'U': Upper triangular part of A is referenced
* = 'L': Lower triangular part of A is referenced
*
* @param[in] N The order of the matrix A. N >= 0. When N = 0,
* CHAMELEON_lapacke_zlanhe is set to zero.
*
* @param[in] A The hermitian matrix A, dimension (LDA,N).
* 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. Note that the
* imaginary parts of the diagonal elements need not be set and are assumed to
* be zero.
*
* @param[in] LDA The leading dimension of the array A. LDA >= max(N,1).
*
******************************************************************************
*
* @retval the norm described above.
*
******************************************************************************
*
* @sa CHAMELEON_lapacke_zlanhe
* @sa CHAMELEON_lapacke_clanhe
*
*/
double CHAMELEON_lapacke_zlanhe( int matrix_layout, char norm, char uplo, int N,
const CHAMELEON_Complex64_t *A, int lda )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zlanhe", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zlanhe( (cham_normtype_t)chameleon_lapacktochameleon_norm(&norm),
(cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
N, (CHAMELEON_Complex64_t *)A, lda );
}
/**
*
* @file lapack_zlansy.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for lansy
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-12-21
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zlansy CHAMELEON_GLOBAL( chameleon_lapack_zlansy, CHAMELEON_LAPACK_ZLANSY )
double CHAMELEON_lapack_zlansy ( const char* norm, const char* uplo, const int* n,
const CHAMELEON_Complex64_t* a, const int* lda,
double* work )
{
(void)work;
return CHAMELEON_lapacke_zlansy( CblasColMajor,
*norm, *uplo, *n, a, *lda );
}
/* C LAPACKE interface */
/**
******************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zlansy - returns the value of the one norm, or the
Frobenius norm, or the infinity norm, or the element of largest absolute
value of a complex symmetric 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] norm = max(abs(A(i,j))), NORM = 'M' or 'm'
* = norm1(A), NORM = '1', 'O' or 'o'
* = normI(A), NORM = 'I' or 'i'
* = normF(A), NORM = 'F', 'f', 'E' or 'e'
*
* where norm1 denotes the one norm of a matrix (maximum column sum), normI
* denotes the infinity norm of a matrix (maximum row sum) and normF denotes
* the Frobenius norm of a matrix (square root of sum of squares). Note that
* max(abs(A(i,j))) is not a consistent matrix norm.
*
* @param[in] uplo Specifies whether the upper or lower triangular part of the
* hermitian matrix A is to be referenced.
* = 'U': Upper triangular part of A is referenced
* = 'L': Lower triangular part of A is referenced
*
* @param[in] N The order of the matrix A. N >= 0. When N = 0,
* CHAMELEON_lapacke_zlansy is set to zero.
*
* @param[in] A The symmetric matrix A, dimension (LDA,N).
* 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.
*
* @param[in] LDA The leading dimension of the array A. LDA >= max(N,1).
*
******************************************************************************
*
* @retval the norm described above.
*
******************************************************************************
*
* @sa CHAMELEON_lapacke_zlansy
* @sa CHAMELEON_lapacke_clansy
* @sa CHAMELEON_lapacke_dlansy
* @sa CHAMELEON_lapacke_slansy
*
*/
double CHAMELEON_lapacke_zlansy( int matrix_layout, char norm, char uplo, int N,
const CHAMELEON_Complex64_t *A, int lda )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zlansy", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zlansy( (cham_normtype_t)chameleon_lapacktochameleon_norm(&norm),
(cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
N, (CHAMELEON_Complex64_t *)A, lda );
}
/**
*
* @file lapack_zlantr.c
*
* @copyright 2022 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
* Univ. Bordeaux. All rights reserved.
*
***
*
* @brief Chameleon lapack and lapacke api for lantr
*
* @version 1.2.0
* @author Mathieu Faverge
* @author Florent Pruvost
* @date 2022-12-21
* @precisions normal z -> s d c
*
*/
#include "chameleon_lapack.h"
#include "lapack_api_common.h"
/* Fortran LAPACK interface */
#define CHAMELEON_lapack_zlantr CHAMELEON_GLOBAL( chameleon_lapack_zlantr, CHAMELEON_LAPACK_ZLANTR )
double CHAMELEON_lapack_zlantr ( const char* norm, const char* uplo, const char* diag,
const int* m, const int* n,
const CHAMELEON_Complex64_t* a, const int* lda,
double* work )
{
(void)work;
return CHAMELEON_lapacke_zlantr( CblasColMajor,
*norm, *uplo, *diag, *m, *n, a, *lda );
}
/* C LAPACKE interface */
/**
******************************************************************************
*
* @ingroup CHAMELEON_LAPACK_API
*
* CHAMELEON_lapacke_zlantr - returns the value of the one norm, or the
Frobenius norm, or the infinity norm, or the element of largest absolute
value of a trapezoidal or 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] norm = max(abs(A(i,j))), NORM = 'M' or 'm'
* = norm1(A), NORM = '1', 'O' or 'o'
* = normI(A), NORM = 'I' or 'i'
* = normF(A), NORM = 'F', 'f', 'E' or 'e'
*
* where norm1 denotes the one norm of a matrix (maximum column sum), normI
* denotes the infinity norm of a matrix (maximum row sum) and normF denotes
* the Frobenius norm of a matrix (square root of sum of squares). Note that
* max(abs(A(i,j))) is not a consistent matrix norm.
*
* @param[in] uplo Specifies whether the matrix A is upper or lower trapezoidal.
* = 'U': Upper trapezoidal
* = 'L': Lower trapezoidal
* Note that A is triangular instead of trapezoidal if M = N.
*
* @param[in] diag Specifies whether or not the matrix A has unit diagonal.
* = 'N': Non-unit diagonal
* = 'U': Unit diagonal
*
* @param[in] M The number of rows of the matrix A. M >= 0. When M = 0, the
* returned value is set to zero.
*
* @param[in] N The number of columns of the matrix A. N >= 0. When N = 0, the
* returned value is set to zero.
*
* @param[in] A The trapezoidal matrix A (A is triangular if M = N),
* dimension (LDA,N).
* If UPLO = 'U', the leading m by n upper trapezoidal part of
* the array A contains the upper trapezoidal matrix, and the
* strictly lower triangular part of A is not referenced.
* If UPLO = 'L', the leading m by n lower trapezoidal part of
* the array A contains the lower trapezoidal matrix, and the
* strictly upper triangular part of A is not referenced. Note
* that when DIAG = 'U', the diagonal elements of A are not
* referenced and are assumed to be one.
*
* @param[in] LDA The leading dimension of the array A. LDA >= max(M,1).
*
******************************************************************************
*
* @retval the norm described above.
*
******************************************************************************
*
* @sa CHAMELEON_lapacke_zlantr
* @sa CHAMELEON_lapacke_clantr
* @sa CHAMELEON_lapacke_dlantr
* @sa CHAMELEON_lapacke_slantr
*
*/
double CHAMELEON_lapacke_zlantr( int matrix_layout, char norm, char uplo, char diag,
int M, int N, const CHAMELEON_Complex64_t *A, int lda )
{
if ( matrix_layout != CblasColMajor ){
fprintf( stderr, "CHAMELEON ERROR: %s(): %s\n", "CHAMELEON_lapacke_zlantr", "illegal value of matrix_layout" );
return -1;
}
return CHAMELEON_zlantr( (cham_normtype_t)chameleon_lapacktochameleon_norm(&norm),
(cham_uplo_t)chameleon_blastocblas_uplo(&uplo),
(cham_diag_t)chameleon_blastocblas_diag(&diag),
M, N, (CHAMELEON_Complex64_t *)A, lda );
}
......@@ -23,9 +23,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
......@@ -124,6 +126,7 @@ testing_zlange_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_normtype_t norm_type = run_arg_get_ntype( args, "norm", ChamMaxNorm );
int N = run_arg_get_int( args, "N", 1000 );
......@@ -151,7 +154,22 @@ testing_zlange_std( run_arg_list_t *args, int check )
testing_stop( &test_data, flops_zlange( norm_type, M, N ) );
#else
testing_start( &test_data );
norm = CHAMELEON_zlange( norm_type, M, N, A, LDA );
switch ( api ) {
case 1:
norm = CHAMELEON_zlange( norm_type, M, N, A, LDA );
break;
#if !defined(CHAMELEON_SIMULATION)
case 2:
norm = CHAMELEON_lapacke_zlange( CblasColMajor, chameleon_lapack_const( norm_type ), M, N, A, LDA );
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_zlange( norm_type, M, N ) );
......
......@@ -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
......@@ -121,6 +123,7 @@ testing_zlanhe_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_normtype_t norm_type = run_arg_get_ntype( args, "norm", ChamMaxNorm );
cham_uplo_t uplo = run_arg_get_uplo( args, "uplo", ChamUpper );
......@@ -151,7 +154,22 @@ testing_zlanhe_std( run_arg_list_t *args, int check )
testing_stop( &test_data, flops_zlanhe( norm_type, N ) );
#else
testing_start( &test_data );
norm = CHAMELEON_zlanhe( norm_type, uplo, N, A, LDA );
switch ( api ) {
case 1:
norm = CHAMELEON_zlanhe( norm_type, uplo, N, A, LDA );
break;
#if !defined(CHAMELEON_SIMULATION)
case 2:
norm = CHAMELEON_lapacke_zlanhe( CblasColMajor, chameleon_lapack_const( norm_type ), chameleon_lapack_const(uplo), N, A, LDA );
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_zlanhe( norm_type, N ) );
......
......@@ -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
......@@ -121,6 +123,7 @@ testing_zlansy_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_normtype_t norm_type = run_arg_get_ntype( args, "norm", ChamMaxNorm );
cham_uplo_t uplo = run_arg_get_uplo( args, "uplo", ChamUpper );
......@@ -151,7 +154,22 @@ testing_zlansy_std( run_arg_list_t *args, int check )
testing_stop( &test_data, flops_zlansy( norm_type, N ) );
#else
testing_start( &test_data );
norm = CHAMELEON_zlansy( norm_type, uplo, N, A, LDA );
switch ( api ) {
case 1:
norm = CHAMELEON_zlansy( norm_type, uplo, N, A, LDA );
break;
#if !defined(CHAMELEON_SIMULATION)
case 2:
norm = CHAMELEON_lapacke_zlansy( CblasColMajor, chameleon_lapack_const( norm_type ), chameleon_lapack_const(uplo), N, A, LDA );
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_zlansy( norm_type, N ) );
......
......@@ -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
......@@ -145,6 +147,7 @@ testing_zlantr_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_normtype_t norm_type = run_arg_get_ntype( args, "norm", ChamMaxNorm );
cham_uplo_t uplo = run_arg_get_uplo( args, "uplo", ChamUpper );
......@@ -174,7 +177,22 @@ testing_zlantr_std( run_arg_list_t *args, int check )
testing_stop( &test_data, flops_zlantr( norm_type, uplo, M, N ) );
#else
testing_start( &test_data );
norm = CHAMELEON_zlantr( norm_type, uplo, diag, M, N, A, LDA );
switch ( api ) {
case 1:
norm = CHAMELEON_zlantr( norm_type, uplo, diag, M, N, A, LDA );
break;
#if !defined(CHAMELEON_SIMULATION)
case 2:
norm = CHAMELEON_lapacke_zlantr( CblasColMajor, chameleon_lapack_const( norm_type ), chameleon_lapack_const(uplo), chameleon_lapack_const(diag), M, N, A, LDA );
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_zlantr( norm_type, uplo, M, N ) );
......
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