/** Perform triangular matrix-vector product. \f$ x = A x \f$ (trans == MagmaNoTrans), or \n \f$ x = A^T x \f$ (trans == MagmaTrans), or \n \f$ x = A^H x \f$ (trans == MagmaConjTrans). @param[in] uplo Whether the upper or lower triangle of A is referenced. @param[in] trans Operation to perform on A. @param[in] diag Whether the diagonal of A is assumed to be unit or non-unit. @param[in] n Number of rows and columns of A. n >= 0. @param[in] dA COMPLEX array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device. @param[in] ldda Leading dimension of dA. @param[in] dx COMPLEX array on GPU device. The n element vector x of dimension (1 + (n-1)*incx). @param[in] incx Stride between consecutive elements of dx. incx != 0. @ingroup magma_cblas2 */ extern "C" void magma_ctrmv( magma_uplo_t uplo, magma_trans_t trans, magma_diag_t diag, magma_int_t n, magmaFloatComplex_const_ptr dA, size_t dA_offset, magma_int_t ldda, magmaFloatComplex_ptr dx, size_t dx_offset, magma_int_t incx, magma_queue_t queue ) { if ( n <= 0 ) return; magmaFloatComplex_ptr dwork; magma_cmalloc( &dwork, (1 + (n-1)*abs(incx)) ); cl_int err = clblasCtrmv( clblasColumnMajor, clblas_uplo_const( uplo ), clblas_trans_const( trans ), clblas_diag_const( diag ), n, dA, dA_offset, ldda, dx, dx_offset, incx, dwork, 1, &queue, 0, NULL, g_event ); clFlush(queue); check_error( err ); magma_free( dwork ); }
/** Perform Hermitian matrix-vector product, \f$ y = \alpha A x + \beta y \f$. @param[in] uplo Whether the upper or lower triangle of A is referenced. @param[in] n Number of rows and columns of A. n >= 0. @param[in] alpha Scalar \f$ \alpha \f$ @param[in] dA COMPLEX array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device. @param[in] ldda Leading dimension of dA. @param[in] dx COMPLEX array on GPU device. The m element vector x of dimension (1 + (m-1)*incx). @param[in] incx Stride between consecutive elements of dx. incx != 0. @param[in] beta Scalar \f$ \beta \f$ @param[in,out] dy COMPLEX array on GPU device. The n element vector y of dimension (1 + (n-1)*incy). @param[in] incy Stride between consecutive elements of dy. incy != 0. @ingroup magma_cblas2 */ extern "C" void magma_chemv( magma_uplo_t uplo, magma_int_t n, magmaFloatComplex alpha, magmaFloatComplex_const_ptr dA, size_t dA_offset, magma_int_t ldda, magmaFloatComplex_const_ptr dx, size_t dx_offset, magma_int_t incx, magmaFloatComplex beta, magmaFloatComplex_ptr dy, size_t dy_offset, magma_int_t incy, magma_queue_t queue ) { if ( n <= 0 ) return; cl_int err = clblasChemv( clblasColumnMajor, clblas_uplo_const( uplo ), n, alpha, dA, dA_offset, ldda, dx, dx_offset, incx, beta, dy, dy_offset, incy, 1, &queue, 0, NULL, g_event ); clFlush(queue); check_error( err ); }
/** Solve triangular matrix-matrix system (multiple right-hand sides). \f$ op(A) X = \alpha B \f$ (side == MagmaLeft), or \n \f$ X op(A) = \alpha B \f$ (side == MagmaRight), \n where \f$ A \f$ is triangular. @param[in] side Whether A is on the left or right. @param[in] uplo Whether A is upper or lower triangular. @param[in] trans Operation to perform on A. @param[in] diag Whether the diagonal of A is assumed to be unit or non-unit. @param[in] m Number of rows of B. m >= 0. @param[in] n Number of columns of B. n >= 0. @param[in] alpha Scalar \f$ \alpha \f$ @param[in] dA COMPLEX array on GPU device. If side == MagmaLeft, the m-by-m triangular matrix A of dimension (ldda,m), ldda >= max(1,m); \n otherwise, the n-by-n triangular matrix A of dimension (ldda,n), ldda >= max(1,n). @param[in] ldda Leading dimension of dA. @param[in,out] dB COMPLEX array on GPU device. On entry, m-by-n matrix B of dimension (lddb,n), lddb >= max(1,m). On exit, overwritten with the solution matrix X. @param[in] lddb Leading dimension of dB. @ingroup magma_cblas3 */ extern "C" void magma_ctrsm( magma_side_t side, magma_uplo_t uplo, magma_trans_t trans, magma_diag_t diag, magma_int_t m, magma_int_t n, magmaFloatComplex alpha, magmaFloatComplex_const_ptr dA, size_t dA_offset, magma_int_t ldda, magmaFloatComplex_ptr dB, size_t dB_offset, magma_int_t lddb, magma_queue_t queue ) { if (m <= 0 || n <= 0) return; cl_int err = clblasCtrsm( clblasColumnMajor, clblas_side_const( side ), clblas_uplo_const( uplo ), clblas_trans_const( trans ), clblas_diag_const( diag ), m, n, alpha, dA, dA_offset, ldda, dB, dB_offset, lddb, 1, &queue, 0, NULL, g_event ); clFlush(queue); check_error( err ); }
/** Perform Hermitian rank-2k update. \f$ C = \alpha A B^T + \alpha B A^T \beta C \f$ (trans == MagmaNoTrans), or \n \f$ C = \alpha A^T B + \alpha B^T A \beta C \f$ (trans == MagmaTrans), \n where \f$ C \f$ is Hermitian. @param[in] uplo Whether the upper or lower triangle of C is referenced. @param[in] trans Operation to perform on A and B. @param[in] n Number of rows and columns of C. n >= 0. @param[in] k Number of columns of A and B (for MagmaNoTrans) or rows of A and B (for MagmaTrans). k >= 0. @param[in] alpha Scalar \f$ \alpha \f$ @param[in] dA COMPLEX array on GPU device. If trans == MagmaNoTrans, the n-by-k matrix A of dimension (ldda,k), ldda >= max(1,n); \n otherwise, the k-by-n matrix A of dimension (ldda,n), ldda >= max(1,k). @param[in] ldda Leading dimension of dA. @param[in] dB COMPLEX array on GPU device. If trans == MagmaNoTrans, the n-by-k matrix B of dimension (lddb,k), lddb >= max(1,n); \n otherwise, the k-by-n matrix B of dimension (lddb,n), lddb >= max(1,k). @param[in] lddb Leading dimension of dB. @param[in] beta Scalar \f$ \beta \f$ @param[in,out] dC COMPLEX array on GPU device. The n-by-n Hermitian matrix C of dimension (lddc,n), lddc >= max(1,n). @param[in] lddc Leading dimension of dC. @ingroup magma_cblas3 */ extern "C" void magma_cher2k( magma_uplo_t uplo, magma_trans_t trans, magma_int_t n, magma_int_t k, magmaFloatComplex alpha, magmaFloatComplex_const_ptr dA, size_t dA_offset, magma_int_t ldda, magmaFloatComplex_const_ptr dB, size_t dB_offset, magma_int_t lddb, float beta, magmaFloatComplex_ptr dC, size_t dC_offset, magma_int_t lddc, magma_queue_t queue ) { if (n <= 0 || k <= 0) return; cl_int err = clblasCher2k( clblasColumnMajor, clblas_uplo_const( uplo ), clblas_trans_const( trans ), n, k, alpha, dA, dA_offset, ldda, dB, dB_offset, lddb, beta, dC, dC_offset, lddc, 1, &queue, 0, NULL, g_event ); clFlush(queue); check_error( err ); }
/** Perform Hermitian matrix-matrix product. \f$ C = \alpha A B + \beta C \f$ (side == MagmaLeft), or \n \f$ C = \alpha B A + \beta C \f$ (side == MagmaRight), \n where \f$ A \f$ is Hermitian. @param[in] side Whether A is on the left or right. @param[in] uplo Whether the upper or lower triangle of A is referenced. @param[in] m Number of rows of C. m >= 0. @param[in] n Number of columns of C. n >= 0. @param[in] alpha Scalar \f$ \alpha \f$ @param[in] dA COMPLEX array on GPU device. If side == MagmaLeft, the m-by-m Hermitian matrix A of dimension (ldda,m), ldda >= max(1,m); \n otherwise, the n-by-n Hermitian matrix A of dimension (ldda,n), ldda >= max(1,n). @param[in] ldda Leading dimension of dA. @param[in] dB COMPLEX array on GPU device. The m-by-n matrix B of dimension (lddb,n), lddb >= max(1,m). @param[in] lddb Leading dimension of dB. @param[in] beta Scalar \f$ \beta \f$ @param[in,out] dC COMPLEX array on GPU device. The m-by-n matrix C of dimension (lddc,n), lddc >= max(1,m). @param[in] lddc Leading dimension of dC. @ingroup magma_cblas3 */ extern "C" void magma_chemm( magma_side_t side, magma_uplo_t uplo, magma_int_t m, magma_int_t n, magmaFloatComplex alpha, magmaFloatComplex_const_ptr dA, size_t dA_offset, magma_int_t ldda, magmaFloatComplex_const_ptr dB, size_t dB_offset, magma_int_t lddb, magmaFloatComplex beta, magmaFloatComplex_ptr dC, size_t dC_offset, magma_int_t lddc, magma_queue_t queue ) { if ( m <= 0 || n <= 0) return; cl_int err = clblasChemm( clblasColumnMajor, clblas_side_const( side ), clblas_uplo_const( uplo ), m, n, alpha, dA, dA_offset, ldda, dB, dB_offset, lddb, beta, dC, dC_offset, lddc, 1, &queue, 0, NULL, g_event ); clFlush(queue); check_error( err ); }
/** Perform Hermitian rank-1 update, \f$ A = \alpha x x^H + A \f$. @param[in] uplo Whether the upper or lower triangle of A is referenced. @param[in] n Number of rows and columns of A. n >= 0. @param[in] alpha Scalar \f$ \alpha \f$ @param[in] dx COMPLEX array on GPU device. The n element vector x of dimension (1 + (n-1)*incx). @param[in] incx Stride between consecutive elements of dx. incx != 0. @param[in,out] dA COMPLEX array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device. @param[in] ldda Leading dimension of dA. @ingroup magma_cblas2 */ extern "C" void magma_cher( magma_uplo_t uplo, magma_int_t n, float alpha, magmaFloatComplex_const_ptr dx, size_t dx_offset, magma_int_t incx, magmaFloatComplex_ptr dA, size_t dA_offset, magma_int_t ldda, magma_queue_t queue ) { cl_int err = clblasCher( clblasColumnMajor, clblas_uplo_const( uplo ), n, alpha, dx, dx_offset, incx, dA, dA_offset, ldda, 1, &queue, 0, NULL, g_event ); check_error( err ); }
/** Perform symmetric rank-k update. \f$ C = \alpha A A^T + \beta C \f$ (trans == MagmaNoTrans), or \n \f$ C = \alpha A^T A + \beta C \f$ (trans == MagmaTrans), \n where \f$ C \f$ is symmetric. @param[in] uplo Whether the upper or lower triangle of C is referenced. @param[in] trans Operation to perform on A. @param[in] n Number of rows and columns of C. n >= 0. @param[in] k Number of columns of A (for MagmaNoTrans) or rows of A (for MagmaTrans). k >= 0. @param[in] alpha Scalar \f$ \alpha \f$ @param[in] dA COMPLEX array on GPU device. If trans == MagmaNoTrans, the n-by-k matrix A of dimension (ldda,k), ldda >= max(1,n); \n otherwise, the k-by-n matrix A of dimension (ldda,n), ldda >= max(1,k). @param[in] ldda Leading dimension of dA. @param[in] beta Scalar \f$ \beta \f$ @param[in,out] dC COMPLEX array on GPU device. The n-by-n symmetric matrix C of dimension (lddc,n), lddc >= max(1,n). @param[in] lddc Leading dimension of dC. @ingroup magma_cblas3 */ extern "C" void magma_csyrk( magma_uplo_t uplo, magma_trans_t trans, magma_int_t n, magma_int_t k, magmaFloatComplex alpha, magmaFloatComplex_const_ptr dA, size_t dA_offset, magma_int_t ldda, magmaFloatComplex beta, magmaFloatComplex_ptr dC, size_t dC_offset, magma_int_t lddc, magma_queue_t queue ) { cl_int err = clblasCsyrk( clblasColumnMajor, clblas_uplo_const( uplo ), clblas_trans_const( trans ), n, k, alpha, dA, dA_offset, ldda, beta, dC, dC_offset, lddc, 1, &queue, 0, NULL, g_event ); check_error( err ); }
/* //////////////////////////////////////////////////////////////////////////// -- Testing ctrsm */ int main( int argc, char** argv) { TESTING_INIT(); real_Double_t gflops, magma_perf=0, magma_time=0, cublas_perf, cublas_time, cpu_perf=0, cpu_time=0; float magma_error=0, cublas_error, lapack_error, work[1]; magma_int_t M, N, info; magma_int_t Ak; magma_int_t sizeA, sizeB; magma_int_t lda, ldb, ldda, lddb; magma_int_t ione = 1; magma_int_t ISEED[4] = {0,0,0,1}; magma_int_t *ipiv; magmaFloatComplex *h_A, *h_B, *h_Bcublas, *h_Bmagma, *h_Blapack, *h_X; magmaFloatComplex_ptr d_A, d_B; magmaFloatComplex c_neg_one = MAGMA_C_NEG_ONE; magmaFloatComplex c_one = MAGMA_C_ONE; magmaFloatComplex alpha = MAGMA_C_MAKE( 0.29, -0.86 ); magma_int_t status = 0; magma_opts opts; opts.parse_opts( argc, argv ); float tol = opts.tolerance * lapackf77_slamch("E"); // pass ngpu = -1 to test multi-GPU code using 1 gpu magma_int_t abs_ngpu = abs( opts.ngpu ); printf("%% side = %s, uplo = %s, transA = %s, diag = %s, ngpu = %d\n", lapack_side_const(opts.side), lapack_uplo_const(opts.uplo), lapack_trans_const(opts.transA), lapack_diag_const(opts.diag), int(abs_ngpu) ); printf("%% M N MAGMA Gflop/s (ms) CUBLAS Gflop/s (ms) CPU Gflop/s (ms) MAGMA CUBLAS LAPACK error\n"); printf("%%============================================================================================================\n"); for( int itest = 0; itest < opts.ntest; ++itest ) { for( int iter = 0; iter < opts.niter; ++iter ) { M = opts.msize[itest]; N = opts.nsize[itest]; gflops = FLOPS_CTRSM(opts.side, M, N) / 1e9; if ( opts.side == MagmaLeft ) { lda = M; Ak = M; } else { lda = N; Ak = N; } ldb = M; ldda = magma_roundup( lda, opts.align ); // multiple of 32 by default lddb = magma_roundup( ldb, opts.align ); // multiple of 32 by default sizeA = lda*Ak; sizeB = ldb*N; TESTING_MALLOC_CPU( h_A, magmaFloatComplex, lda*Ak ); TESTING_MALLOC_CPU( h_B, magmaFloatComplex, ldb*N ); TESTING_MALLOC_CPU( h_X, magmaFloatComplex, ldb*N ); TESTING_MALLOC_CPU( h_Blapack, magmaFloatComplex, ldb*N ); TESTING_MALLOC_CPU( h_Bcublas, magmaFloatComplex, ldb*N ); TESTING_MALLOC_CPU( h_Bmagma, magmaFloatComplex, ldb*N ); TESTING_MALLOC_CPU( ipiv, magma_int_t, Ak ); TESTING_MALLOC_DEV( d_A, magmaFloatComplex, ldda*Ak ); TESTING_MALLOC_DEV( d_B, magmaFloatComplex, lddb*N ); /* Initialize the matrices */ /* Factor A into LU to get well-conditioned triangular matrix. * Copy L to U, since L seems okay when used with non-unit diagonal * (i.e., from U), while U fails when used with unit diagonal. */ lapackf77_clarnv( &ione, ISEED, &sizeA, h_A ); lapackf77_cgetrf( &Ak, &Ak, h_A, &lda, ipiv, &info ); for( int j = 0; j < Ak; ++j ) { for( int i = 0; i < j; ++i ) { *h_A(i,j) = *h_A(j,i); } } lapackf77_clarnv( &ione, ISEED, &sizeB, h_B ); memcpy( h_Blapack, h_B, sizeB*sizeof(magmaFloatComplex) ); magma_csetmatrix( Ak, Ak, h_A, lda, d_A, ldda, opts.queue ); /* ===================================================================== Performs operation using MAGMABLAS =================================================================== */ #if defined(HAVE_CUBLAS) magma_csetmatrix( M, N, h_B, ldb, d_B, lddb, opts.queue ); magma_time = magma_sync_wtime( opts.queue ); if (opts.ngpu == 1) { magmablas_ctrsm( opts.side, opts.uplo, opts.transA, opts.diag, M, N, alpha, d_A, ldda, d_B, lddb, opts.queue ); } else { magma_ctrsm_m( abs_ngpu, opts.side, opts.uplo, opts.transA, opts.diag, M, N, alpha, d_A, ldda, d_B, lddb ); } magma_time = magma_sync_wtime( opts.queue ) - magma_time; magma_perf = gflops / magma_time; magma_cgetmatrix( M, N, d_B, lddb, h_Bmagma, ldb, opts.queue ); #endif /* ===================================================================== Performs operation using CUBLAS =================================================================== */ magma_csetmatrix( M, N, h_B, ldb, d_B, lddb, opts.queue ); cublas_time = magma_sync_wtime( opts.queue ); #if defined(HAVE_CUBLAS) // opts.handle also uses opts.queue cublasCtrsm( opts.handle, cublas_side_const(opts.side), cublas_uplo_const(opts.uplo), cublas_trans_const(opts.transA), cublas_diag_const(opts.diag), M, N, &alpha, d_A, ldda, d_B, lddb ); #elif defined(HAVE_clBLAS) clblasCtrsm( clblasColumnMajor, clblas_side_const(opts.side), clblas_uplo_const(opts.uplo), clblas_trans_const(opts.transA), clblas_diag_const(opts.diag), M, N, alpha, d_A, 0, ldda, d_B, 0, lddb, 1, &opts.queue, 0, NULL, NULL ); #endif cublas_time = magma_sync_wtime( opts.queue ) - cublas_time; cublas_perf = gflops / cublas_time; magma_cgetmatrix( M, N, d_B, lddb, h_Bcublas, ldb, opts.queue ); /* ===================================================================== Performs operation using CPU BLAS =================================================================== */ if ( opts.lapack ) { cpu_time = magma_wtime(); blasf77_ctrsm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo), lapack_trans_const(opts.transA), lapack_diag_const(opts.diag), &M, &N, &alpha, h_A, &lda, h_Blapack, &ldb ); cpu_time = magma_wtime() - cpu_time; cpu_perf = gflops / cpu_time; } /* ===================================================================== Check the result =================================================================== */ // ||b - 1/alpha*A*x|| / (||A||*||x||) magmaFloatComplex inv_alpha = MAGMA_C_DIV( c_one, alpha ); float normR, normX, normA; normA = lapackf77_clange( "M", &Ak, &Ak, h_A, &lda, work ); #if defined(HAVE_CUBLAS) // check magma memcpy( h_X, h_Bmagma, sizeB*sizeof(magmaFloatComplex) ); blasf77_ctrmm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo), lapack_trans_const(opts.transA), lapack_diag_const(opts.diag), &M, &N, &inv_alpha, h_A, &lda, h_X, &ldb ); blasf77_caxpy( &sizeB, &c_neg_one, h_B, &ione, h_X, &ione ); normR = lapackf77_clange( "M", &M, &N, h_X, &ldb, work ); normX = lapackf77_clange( "M", &M, &N, h_Bmagma, &ldb, work ); magma_error = normR/(normX*normA); #endif // check cublas memcpy( h_X, h_Bcublas, sizeB*sizeof(magmaFloatComplex) ); blasf77_ctrmm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo), lapack_trans_const(opts.transA), lapack_diag_const(opts.diag), &M, &N, &inv_alpha, h_A, &lda, h_X, &ldb ); blasf77_caxpy( &sizeB, &c_neg_one, h_B, &ione, h_X, &ione ); normR = lapackf77_clange( "M", &M, &N, h_X, &ldb, work ); normX = lapackf77_clange( "M", &M, &N, h_Bcublas, &ldb, work ); cublas_error = normR/(normX*normA); if ( opts.lapack ) { // check lapack // this verifies that the matrix wasn't so bad that it couldn't be solved accurately. memcpy( h_X, h_Blapack, sizeB*sizeof(magmaFloatComplex) ); blasf77_ctrmm( lapack_side_const(opts.side), lapack_uplo_const(opts.uplo), lapack_trans_const(opts.transA), lapack_diag_const(opts.diag), &M, &N, &inv_alpha, h_A, &lda, h_X, &ldb ); blasf77_caxpy( &sizeB, &c_neg_one, h_B, &ione, h_X, &ione ); normR = lapackf77_clange( "M", &M, &N, h_X, &ldb, work ); normX = lapackf77_clange( "M", &M, &N, h_Blapack, &ldb, work ); lapack_error = normR/(normX*normA); printf("%5d %5d %7.2f (%7.2f) %7.2f (%7.2f) %7.2f (%7.2f) %8.2e %8.2e %8.2e %s\n", (int) M, (int) N, magma_perf, 1000.*magma_time, cublas_perf, 1000.*cublas_time, cpu_perf, 1000.*cpu_time, magma_error, cublas_error, lapack_error, (magma_error < tol && cublas_error < tol? "ok" : "failed")); status += ! (magma_error < tol && cublas_error < tol); } else { printf("%5d %5d %7.2f (%7.2f) %7.2f (%7.2f) --- ( --- ) %8.2e %8.2e --- %s\n", (int) M, (int) N, magma_perf, 1000.*magma_time, cublas_perf, 1000.*cublas_time, magma_error, cublas_error, (magma_error < tol && cublas_error < tol ? "ok" : "failed")); status += ! (magma_error < tol && cublas_error < tol); } TESTING_FREE_CPU( h_A ); TESTING_FREE_CPU( h_B ); TESTING_FREE_CPU( h_X ); TESTING_FREE_CPU( h_Blapack ); TESTING_FREE_CPU( h_Bcublas ); TESTING_FREE_CPU( h_Bmagma ); TESTING_FREE_CPU( ipiv ); TESTING_FREE_DEV( d_A ); TESTING_FREE_DEV( d_B ); fflush( stdout ); } if ( opts.niter > 1 ) { printf( "\n" ); } } opts.cleanup(); TESTING_FINALIZE(); return status; }