Epetra Package Browser (Single Doxygen Collection)
Development
|
Epetra_BLAS: The Epetra BLAS Wrapper Class. More...
#include <Epetra_BLAS.h>
Constructors/Destructor | |
Epetra_BLAS (void) | |
Epetra_BLAS Constructor. More... | |
Epetra_BLAS (const Epetra_BLAS &BLAS) | |
Epetra_BLAS Copy Constructor. More... | |
virtual | ~Epetra_BLAS (void) |
Epetra_BLAS Destructor. More... | |
Level 1 BLAS | |
float | ASUM (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS one norm function (SASUM). More... | |
double | ASUM (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS one norm function (DASUM). More... | |
float | DOT (const int N, const float *X, const float *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS dot product function (SDOT). More... | |
double | DOT (const int N, const double *X, const double *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS dot product function (DDOT). More... | |
float | NRM2 (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS norm function (SNRM2). More... | |
double | NRM2 (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS norm function (DNRM2). More... | |
void | SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const |
Epetra_BLAS vector scale function (SSCAL) More... | |
void | SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const |
Epetra_BLAS vector scale function (DSCAL) More... | |
void | COPY (const int N, const float *X, float *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS vector copy function (SCOPY) More... | |
void | COPY (const int N, const double *X, double *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS vector scale function (DCOPY) More... | |
int | IAMAX (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (ISAMAX) More... | |
int | IAMAX (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (IDAMAX) More... | |
void | AXPY (const int N, const float ALPHA, const float *X, float *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS vector update function (SAXPY) More... | |
void | AXPY (const int N, const double ALPHA, const double *X, double *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS vector update function (DAXPY) More... | |
Level 2 BLAS | |
void | GEMV (const char TRANS, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *X, const float BETA, float *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS matrix-vector multiply function (SGEMV) More... | |
void | GEMV (const char TRANS, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *X, const double BETA, double *Y, const int INCX=1, const int INCY=1) const |
Epetra_BLAS matrix-vector multiply function (DGEMV) More... | |
Level 3 BLAS | |
void | GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const |
Epetra_BLAS matrix-matrix multiply function (SGEMM) More... | |
void | GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const |
Epetra_BLAS matrix-matrix multiply function (DGEMM) More... | |
void | SYMM (const char SIDE, const char UPLO, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const |
Epetra_BLAS symmetric matrix-matrix multiply function (SSYMM) More... | |
void | SYMM (const char SIDE, const char UPLO, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const |
Epetra_BLAS matrix-matrix multiply function (DSYMM) More... | |
void | TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const float ALPHA, const float *A, const int LDA, float *B, const int LDB) const |
Epetra_BLAS triangular matrix-matrix multiply function (STRMM) More... | |
void | TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const double ALPHA, const double *A, const int LDA, double *B, const int LDB) const |
Epetra_BLAS triangular matrix-matrix multiply function (DTRMM) More... | |
void | SYRK (const char UPLO, const char TRANS, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float BETA, float *C, const int LDC) const |
Eperta_BLAS symetric rank k funtion (ssyrk) More... | |
void | SYRK (const char UPLO, const char TRANS, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double BETA, double *C, const int LDC) const |
Eperta_BLAS symetric rank k funtion (dsyrk) More... | |
Epetra_BLAS: The Epetra BLAS Wrapper Class.
The Epetra_BLAS class is a wrapper that encapsulates the BLAS (Basic Linear Algebra Subprograms). The BLAS provide portable, high- performance implementations of kernels such as dense vectoer multiplication, dot products, dense matrix-vector multiplication and dense matrix-matrix multiplication.
The standard BLAS interface is Fortran-specific. Unfortunately, the interface between C++ and Fortran is not standard across all computer platforms. The Epetra_BLAS class provides C++ wrappers for the BLAS kernels in order to insulate the rest of Epetra from the details of C++ to Fortran translation. A Epetra_BLAS object is essentially nothing, but allows access to the BLAS wrapper functions.
Epetra_BLAS is a serial interface only. This is appropriate since the standard BLAS are only specified for serial execution (or shared memory parallel).
Definition at line 70 of file Epetra_BLAS.h.
|
inline |
Epetra_BLAS Constructor.
Builds an instance of a serial BLAS object.
Definition at line 179 of file Epetra_BLAS.h.
|
inline |
Epetra_BLAS Copy Constructor.
Makes an exact copy of an existing Epetra_BLAS instance.
Definition at line 181 of file Epetra_BLAS.h.
|
inlinevirtual |
Epetra_BLAS Destructor.
Definition at line 183 of file Epetra_BLAS.h.
float Epetra_BLAS::ASUM | ( | const int | N, |
const float * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS one norm function (SASUM).
Definition at line 65 of file Epetra_BLAS.cpp.
double Epetra_BLAS::ASUM | ( | const int | N, |
const double * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS one norm function (DASUM).
Definition at line 69 of file Epetra_BLAS.cpp.
float Epetra_BLAS::DOT | ( | const int | N, |
const float * | X, | ||
const float * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS dot product function (SDOT).
Definition at line 73 of file Epetra_BLAS.cpp.
double Epetra_BLAS::DOT | ( | const int | N, |
const double * | X, | ||
const double * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS dot product function (DDOT).
Definition at line 77 of file Epetra_BLAS.cpp.
float Epetra_BLAS::NRM2 | ( | const int | N, |
const float * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS norm function (SNRM2).
Definition at line 81 of file Epetra_BLAS.cpp.
double Epetra_BLAS::NRM2 | ( | const int | N, |
const double * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS norm function (DNRM2).
Definition at line 85 of file Epetra_BLAS.cpp.
void Epetra_BLAS::SCAL | ( | const int | N, |
const float | ALPHA, | ||
float * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS vector scale function (SSCAL)
Definition at line 89 of file Epetra_BLAS.cpp.
void Epetra_BLAS::SCAL | ( | const int | N, |
const double | ALPHA, | ||
double * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS vector scale function (DSCAL)
Definition at line 94 of file Epetra_BLAS.cpp.
void Epetra_BLAS::COPY | ( | const int | N, |
const float * | X, | ||
float * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS vector copy function (SCOPY)
Definition at line 99 of file Epetra_BLAS.cpp.
void Epetra_BLAS::COPY | ( | const int | N, |
const double * | X, | ||
double * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS vector scale function (DCOPY)
Definition at line 104 of file Epetra_BLAS.cpp.
int Epetra_BLAS::IAMAX | ( | const int | N, |
const float * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS arg maximum of absolute value function (ISAMAX)
Definition at line 109 of file Epetra_BLAS.cpp.
int Epetra_BLAS::IAMAX | ( | const int | N, |
const double * | X, | ||
const int | INCX = 1 |
||
) | const |
Epetra_BLAS arg maximum of absolute value function (IDAMAX)
Definition at line 113 of file Epetra_BLAS.cpp.
void Epetra_BLAS::AXPY | ( | const int | N, |
const float | ALPHA, | ||
const float * | X, | ||
float * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS vector update function (SAXPY)
Definition at line 117 of file Epetra_BLAS.cpp.
void Epetra_BLAS::AXPY | ( | const int | N, |
const double | ALPHA, | ||
const double * | X, | ||
double * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS vector update function (DAXPY)
Definition at line 121 of file Epetra_BLAS.cpp.
void Epetra_BLAS::GEMV | ( | const char | TRANS, |
const int | M, | ||
const int | N, | ||
const float | ALPHA, | ||
const float * | A, | ||
const int | LDA, | ||
const float * | X, | ||
const float | BETA, | ||
float * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS matrix-vector multiply function (SGEMV)
Definition at line 125 of file Epetra_BLAS.cpp.
void Epetra_BLAS::GEMV | ( | const char | TRANS, |
const int | M, | ||
const int | N, | ||
const double | ALPHA, | ||
const double * | A, | ||
const int | LDA, | ||
const double * | X, | ||
const double | BETA, | ||
double * | Y, | ||
const int | INCX = 1 , |
||
const int | INCY = 1 |
||
) | const |
Epetra_BLAS matrix-vector multiply function (DGEMV)
Definition at line 132 of file Epetra_BLAS.cpp.
void Epetra_BLAS::GEMM | ( | const char | TRANSA, |
const char | TRANSB, | ||
const int | M, | ||
const int | N, | ||
const int | K, | ||
const float | ALPHA, | ||
const float * | A, | ||
const int | LDA, | ||
const float * | B, | ||
const int | LDB, | ||
const float | BETA, | ||
float * | C, | ||
const int | LDC | ||
) | const |
Epetra_BLAS matrix-matrix multiply function (SGEMM)
Definition at line 140 of file Epetra_BLAS.cpp.
void Epetra_BLAS::GEMM | ( | const char | TRANSA, |
const char | TRANSB, | ||
const int | M, | ||
const int | N, | ||
const int | K, | ||
const double | ALPHA, | ||
const double * | A, | ||
const int | LDA, | ||
const double * | B, | ||
const int | LDB, | ||
const double | BETA, | ||
double * | C, | ||
const int | LDC | ||
) | const |
Epetra_BLAS matrix-matrix multiply function (DGEMM)
Definition at line 149 of file Epetra_BLAS.cpp.
void Epetra_BLAS::SYMM | ( | const char | SIDE, |
const char | UPLO, | ||
const int | M, | ||
const int | N, | ||
const float | ALPHA, | ||
const float * | A, | ||
const int | LDA, | ||
const float * | B, | ||
const int | LDB, | ||
const float | BETA, | ||
float * | C, | ||
const int | LDC | ||
) | const |
Epetra_BLAS symmetric matrix-matrix multiply function (SSYMM)
Definition at line 157 of file Epetra_BLAS.cpp.
void Epetra_BLAS::SYMM | ( | const char | SIDE, |
const char | UPLO, | ||
const int | M, | ||
const int | N, | ||
const double | ALPHA, | ||
const double * | A, | ||
const int | LDA, | ||
const double * | B, | ||
const int | LDB, | ||
const double | BETA, | ||
double * | C, | ||
const int | LDC | ||
) | const |
Epetra_BLAS matrix-matrix multiply function (DSYMM)
Definition at line 166 of file Epetra_BLAS.cpp.
void Epetra_BLAS::TRMM | ( | const char | SIDE, |
const char | UPLO, | ||
const char | TRANSA, | ||
const char | DIAG, | ||
const int | M, | ||
const int | N, | ||
const float | ALPHA, | ||
const float * | A, | ||
const int | LDA, | ||
float * | B, | ||
const int | LDB | ||
) | const |
Epetra_BLAS triangular matrix-matrix multiply function (STRMM)
Definition at line 174 of file Epetra_BLAS.cpp.
void Epetra_BLAS::TRMM | ( | const char | SIDE, |
const char | UPLO, | ||
const char | TRANSA, | ||
const char | DIAG, | ||
const int | M, | ||
const int | N, | ||
const double | ALPHA, | ||
const double * | A, | ||
const int | LDA, | ||
double * | B, | ||
const int | LDB | ||
) | const |
Epetra_BLAS triangular matrix-matrix multiply function (DTRMM)
Definition at line 182 of file Epetra_BLAS.cpp.
void Epetra_BLAS::SYRK | ( | const char | UPLO, |
const char | TRANS, | ||
const int | N, | ||
const int | K, | ||
const float | ALPHA, | ||
const float * | A, | ||
const int | LDA, | ||
const float | BETA, | ||
float * | C, | ||
const int | LDC | ||
) | const |
Eperta_BLAS symetric rank k funtion (ssyrk)
Definition at line 190 of file Epetra_BLAS.cpp.
void Epetra_BLAS::SYRK | ( | const char | UPLO, |
const char | TRANS, | ||
const int | N, | ||
const int | K, | ||
const double | ALPHA, | ||
const double * | A, | ||
const int | LDA, | ||
const double | BETA, | ||
double * | C, | ||
const int | LDC | ||
) | const |
Eperta_BLAS symetric rank k funtion (dsyrk)
Definition at line 195 of file Epetra_BLAS.cpp.