|
OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
|


Go to the source code of this file.
Enumerations | |
| enum | CBLAS_TRANSPOSE { CblasNoTrans = 111 , CblasTrans = 112 , CblasConjTrans = 113 } |
| Transpose operation selector for BLAS routines. More... | |
Functions | |
| void | gsl_blas_dgemm (CBLAS_TRANSPOSE TransA, CBLAS_TRANSPOSE TransB, double alpha, const gsl_matrix *A, const gsl_matrix *B, double beta, gsl_matrix *C) |
| Real matrix-matrix multiply and accumulate. | |
| void | gsl_blas_zgemm (CBLAS_TRANSPOSE TransA, CBLAS_TRANSPOSE TransB, gsl_complex alpha, const gsl_matrix_complex *A, const gsl_matrix_complex *B, gsl_complex beta, gsl_matrix_complex *C) |
| Complex matrix-matrix multiply and accumulate. | |
| void | gsl_blas_dgemv (CBLAS_TRANSPOSE TransA, double alpha, const gsl_matrix *A, const gsl_vector *x, double beta, gsl_vector *y) |
| Real matrix-vector multiply and accumulate. | |
| void | gsl_blas_zgemv (CBLAS_TRANSPOSE TransA, gsl_complex alpha, const gsl_matrix_complex *A, const gsl_vector_complex *x, gsl_complex beta, gsl_vector_complex *y) |
| Complex matrix-vector multiply and accumulate. | |
| enum CBLAS_TRANSPOSE |
Transpose operation selector for BLAS routines.
CblasNoTrans: \(\mathrm{op}(A) = A\)CblasTrans: \(\mathrm{op}(A) = A^T\)CblasConjTrans: \(\mathrm{op}(A) = A^H\) (conjugate transpose) | Enumerator | |
|---|---|
| CblasNoTrans | |
| CblasTrans | |
| CblasConjTrans | |
|
inline |
Real matrix-matrix multiply and accumulate.
Computes \(C \leftarrow \alpha \,\mathrm{op}(A)\,\mathrm{op}(B) + \beta\, C\) where \(\mathrm{op}(X)\) is chosen by TransA and TransB.
| TransA | Input: operation applied to A. |
| TransB | Input: operation applied to B. |
| alpha | Input: scalar multiplier for \(\mathrm{op}(A)\mathrm{op}(B)\). |
| A | Input: left matrix operand. |
| B | Input: right matrix operand. |
| beta | Input: scalar multiplier for existing C. |
| C | Input/Output: accumulation target updated in-place. |
Definition at line 44 of file GSLBLAS.h.
References CblasNoTrans, gsl_matrix_ptr(), gsl_matrix::size1, and gsl_matrix::size2.
Referenced by TEST_F(), TEST_F(), and TEST_F().

|
inline |
Real matrix-vector multiply and accumulate.
Computes \(y \leftarrow \alpha\,\mathrm{op}(A)\,x + \beta\,y\) with \(\mathrm{op}(A) = A\) or \(A^T\).
| TransA | Input: operation applied to A. |
| alpha | Input: scalar multiplier for \(\mathrm{op}(A)x\). |
| A | Input: matrix operand. |
| x | Input: vector operand. |
| beta | Input: scalar multiplier for existing y. |
| y | Input/Output: result vector updated in-place. |
Definition at line 144 of file GSLBLAS.h.
References CblasNoTrans, gsl_matrix_ptr(), gsl_vector_ptr(), gsl_vector::size, gsl_matrix::size1, and gsl_matrix::size2.
Referenced by TEST_F().

|
inline |
Complex matrix-matrix multiply and accumulate.
Computes \(C \leftarrow \alpha \,\mathrm{op}(A)\,\mathrm{op}(B) + \beta\, C\) with \(\mathrm{op}(X)\) defined by TransA/(including conjugate transpose). TransB
| TransA | Input: operation applied to A. |
| TransB | Input: operation applied to B. |
| alpha | Input: scalar multiplier for \(\mathrm{op}(A)\mathrm{op}(B)\). |
| A | Input: left matrix operand. |
| B | Input: right matrix operand. |
| beta | Input: scalar multiplier for existing C. |
| C | Input/Output: accumulation target updated in-place. |
Definition at line 90 of file GSLBLAS.h.
References CblasConjTrans, CblasNoTrans, gsl_complex::dat, gsl_complex_add(), gsl_complex_conjugate(), gsl_complex_mul(), gsl_matrix_complex_ptr(), gsl_matrix_complex::size1, and gsl_matrix_complex::size2.
Referenced by TEST_F().

|
inline |
Complex matrix-vector multiply and accumulate.
Computes \(y \leftarrow \alpha\,\mathrm{op}(A)\,x + \beta\,y\) with \(\mathrm{op}(A) = A, A^T,\) or \(A^H\).
| TransA | Input: operation applied to A. |
| alpha | Input: scalar multiplier for \(\mathrm{op}(A)x\). |
| A | Input: matrix operand. |
| x | Input: vector operand. |
| beta | Input: scalar multiplier for existing y. |
| y | Input/Output: result vector updated in-place. |
Definition at line 194 of file GSLBLAS.h.
References CblasConjTrans, CblasNoTrans, gsl_complex_add(), gsl_complex_conjugate(), gsl_complex_mul(), gsl_matrix_complex_ptr(), gsl_vector_complex_ptr(), gsl_vector_complex::size, gsl_matrix_complex::size1, and gsl_matrix_complex::size2.
