Teuchos - Trilinos Tools Package  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Teuchos::BLAS< OrdinalType, ScalarType > Class Template Reference

Templated BLAS wrapper. More...

#include <Teuchos_BLAS.hpp>

Inheritance diagram for Teuchos::BLAS< OrdinalType, ScalarType >:
Teuchos::DefaultBLASImpl< OrdinalType, ScalarType > Teuchos::SerialBandDenseMatrix< OrdinalType, ScalarType > Teuchos::SerialBandDenseSolver< OrdinalType, ScalarType > Teuchos::SerialDenseMatrix< OrdinalType, ScalarType > Teuchos::SerialDenseSolver< OrdinalType, ScalarType > Teuchos::SerialQRDenseSolver< OrdinalType, ScalarType > Teuchos::SerialSpdDenseSolver< OrdinalType, ScalarType > Teuchos::SerialSymDenseMatrix< OrdinalType, ScalarType > Teuchos::SerialTriDiMatrix< OrdinalType, ScalarType >

Public Member Functions

Constructor/Destructor.
 BLAS (void)
 Default constructor. More...
 
 BLAS (const BLAS< OrdinalType, ScalarType > &)
 Copy constructor. More...
 
virtual ~BLAS (void)
 Destructor. More...
 
- Public Member Functions inherited from Teuchos::DefaultBLASImpl< OrdinalType, ScalarType >
 DefaultBLASImpl (void)
 Default constructor. More...
 
 DefaultBLASImpl (const DefaultBLASImpl< OrdinalType, ScalarType > &)
 Copy constructor. More...
 
virtual ~DefaultBLASImpl (void)
 Destructor. More...
 
template<typename alpha_type , typename A_type , typename x_type , typename beta_type >
void GEMV (ETransp trans, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const x_type *x, const OrdinalType &incx, const beta_type beta, ScalarType *y, const OrdinalType &incy) const
 Performs the matrix-vector operation: y <- alpha*A*x+beta*y or y <- alpha*A'*x+beta*y where A is a general m by n matrix. More...
 
template<typename A_type >
void TRMV (EUplo uplo, ETransp trans, EDiag diag, const OrdinalType &n, const A_type *A, const OrdinalType &lda, ScalarType *x, const OrdinalType &incx) const
 Performs the matrix-vector operation: x <- A*x or x <- A'*x where A is a unit/non-unit n by n upper/lower triangular matrix. More...
 
template<typename alpha_type , typename x_type , typename y_type >
void GER (const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy, ScalarType *A, const OrdinalType &lda) const
 Performs the rank 1 operation: A <- alpha*x*y'+A. More...
 
template<typename alpha_type , typename A_type , typename B_type , typename beta_type >
void GEMM (ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 General matrix-matrix multiply. More...
 
void SWAP (const OrdinalType &n, ScalarType *const x, const OrdinalType &incx, ScalarType *const y, const OrdinalType &incy) const
 Swap the entries of x and y. More...
 
template<typename alpha_type , typename A_type , typename B_type , typename beta_type >
void SYMM (ESide side, EUplo uplo, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 Performs the matrix-matrix operation: C <- alpha*A*B+beta*C or C <- alpha*B*A+beta*C where A is an m by m or n by n symmetric matrix and B is a general matrix. More...
 
template<typename alpha_type , typename A_type , typename beta_type >
void SYRK (EUplo uplo, ETransp trans, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 Performs the symmetric rank k operation: C <- alpha*A*A'+beta*C or C <- alpha*A'*A+beta*C, where alpha and beta are scalars, C is an n by n symmetric matrix and A is an n by k matrix in the first case or k by n matrix in the second case. More...
 
template<typename alpha_type , typename A_type >
void TRMM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
 Performs the matrix-matrix operation: B <- alpha*op(A)*B or B <- alpha*B*op(A) where op(A) is an unit/non-unit, upper/lower triangular matrix and B is a general matrix. More...
 
template<typename alpha_type , typename A_type >
void TRSM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
 Solves the matrix equations: op(A)*X=alpha*B or X*op(A)=alpha*B where X and B are m by n matrices, A is a unit/non-unit, upper/lower triangular matrix and op(A) is A or A'. The matrix X is overwritten on B. More...
 
void ROTG (ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
 Computes a Givens plane rotation. More...
 
void ROT (const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
 Applies a Givens plane rotation. More...
 
void SCAL (const OrdinalType &n, const ScalarType &alpha, ScalarType *x, const OrdinalType &incx) const
 Scale the vector x by the constant alpha. More...
 
void COPY (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
 Copy the vector x to the vector y. More...
 
template<typename alpha_type , typename x_type >
void AXPY (const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
 Perform the operation: y <- y+alpha*x. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
ASUM (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Sum the absolute values of the entries of x. More...
 
template<typename x_type , typename y_type >
ScalarType DOT (const OrdinalType &n, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy) const
 Form the dot product of the vectors x and y. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
NRM2 (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Compute the 2-norm of the vector x. More...
 
OrdinalType IAMAX (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Return the index of the element of x with the maximum magnitude. More...
 

Detailed Description

template<typename OrdinalType, typename ScalarType>
class Teuchos::BLAS< OrdinalType, ScalarType >

Templated BLAS wrapper.

The Teuchos::BLAS class provides functionality similar to the BLAS (Basic Linear Algebra Subprograms). The BLAS provide portable, high- performance implementations of kernels such as dense vector sums, inner products, and norms (the BLAS 1 routines), dense matrix-vector multiplication and triangular solve (the BLAS 2 routines), and dense matrix-matrix multiplication and triangular solve with multiple right-hand sides (the BLAS 3 routines).

The standard BLAS interface is Fortran-specific. Unfortunately, the interface between C++ and Fortran is not standard across all computer platforms. The Teuchos::BLAS class provides C++ bindings for the BLAS kernels in order to insulate the rest of Petra from the details of C++ to Fortran translation.

In addition to giving access to the standard BLAS functionality, Teuchos::BLAS also provides a generic fall-back implementation for any ScalarType class that defines the +, - * and / operators.

Teuchos::BLAS only operates within a single shared-memory space, just like the BLAS. It does not attempt to implement distributed-memory parallel matrix operations.

Note
This class has specializations for ScalarType=float and double, which use the BLAS library directly. If you configure Teuchos to enable complex arithmetic support, via the CMake option -DTeuchos_ENABLE_COMPLEX:BOOL=ON, then this class will also invoke the BLAS library directly for ScalarType=std::complex<float> and std::complex<double>.
Examples:
BLAS/cxx_main.cpp.

Definition at line 244 of file Teuchos_BLAS.hpp.

Constructor & Destructor Documentation

template<typename OrdinalType, typename ScalarType>
Teuchos::BLAS< OrdinalType, ScalarType >::BLAS ( void  )
inline

Default constructor.

Definition at line 254 of file Teuchos_BLAS.hpp.

template<typename OrdinalType, typename ScalarType>
Teuchos::BLAS< OrdinalType, ScalarType >::BLAS ( const BLAS< OrdinalType, ScalarType > &  )
inline

Copy constructor.

Definition at line 257 of file Teuchos_BLAS.hpp.

template<typename OrdinalType, typename ScalarType>
virtual Teuchos::BLAS< OrdinalType, ScalarType >::~BLAS ( void  )
inlinevirtual

Destructor.

Definition at line 260 of file Teuchos_BLAS.hpp.


The documentation for this class was generated from the following file: