Go to the documentation of this file.
36#ifndef OPENRS_BLAS_LAPACK_HEADER
37#define OPENRS_BLAS_LAPACK_HEADER
39#include <opm/common/ErrorMacros.hpp>
49#define DGEMV F77_NAME(dgemv,DGEMV)
53 const int* m , const int* n,
54 const double* a1 , const double* A, const int* ldA ,
55 const double* x, const int* incX,
56 const double* a2 , double* y, const int* incY);
62#define DGEMM F77_NAME(dgemm,DGEMM)
66 const int* m , const int* n , const int* k ,
67 const double* a1 , const double* A , const int* ldA,
68 const double* B , const int* ldB,
69 const double* a2 , double* C , const int* ldC);
75#define DSYRK F77_NAME(dsyrk,DSYRK)
79 const int* n , const int* k ,
80 const double* a1 , const double* A , const int* ldA,
81 const double* a2 , double* C , const int* ldC);
87#define DTRMM F77_NAME(dtrmm,DTRMM)
92 const int* m , const int* n ,
94 const double* A , const int* ldA ,
95 double* B , const int* ldB);
101#define DGEQRF F77_NAME(dgeqrf,DGEQRF)
103 void DGEQRF( const int* m , const int* n ,
104 double* A , const int* ld ,
105 double* tau , double* work,
106 const int* lwork, int* info);
112#define DORGQR F77_NAME(dorgqr,DORGQR)
114 void DORGQR( const int* m , const int* n , const int* k ,
115 double* A , const int* ld , const double* tau,
116 double* work, const int* lwork, int* info);
121#define DGETRF F77_NAME(dgetrf,DGETRF)
123 void DGETRF( const int* m , const int* n ,
124 double* A , const int* ld,
125 int* ipiv, int* info);
130#define DGETRI F77_NAME(dgetri,DGETRI)
133 double* A , const int* ld,
135 double* work, int* lwork, int* info);
142 namespace BLAS_LAPACK {
188 const int m , const int n,
189 const T& a1 , const T* A, const int ldA ,
190 const T* x, const int incX,
191 const T& a2 , T* y, const int incY);
196 const int m , const int n,
197 const double& a1 , const double* A, const int ldA,
198 const double* x, const int incX,
199 const double& a2 , double* y, const int incY);
262 void GEMM( const char* transA, const char* transB,
263 const int m , const int n , const int k ,
264 const T& a1 , const T* A , const int ldA,
265 const T* B , const int ldB,
266 const T& a2 , T* C , const int ldC);
271 const int m , const int n , const int k ,
272 const double& a1 , const double* A , const int ldA,
273 const double* B , const int ldB,
274 const double& a2 , double* C , const int ldC);
282 void SYRK( const char* uplo, const char* trans,
283 const int n , const int k ,
284 const T& a1 , const T* A , const int ldA,
285 const T& a2 , T* C , const int ldC);
290 const int n , const int k ,
291 const double& a1 , const double* A , const int ldA,
292 const double& a2 , double* C , const int ldC);
299 void TRMM( const char* side , const char* uplo,
300 const char* transA, const char* diag,
301 const int m , const int n , const T& a,
302 const T* A , const int ldA ,
303 T* B , const int ldB);
308 const char* transA, const char* diag,
309 const int m , const int n , const double& a,
310 const double* A , const int ldA ,
311 double* B , const int ldB);
318 void GEQRF( const int m , const int n ,
319 T* A , const int ld ,
321 const int lwork, int& info);
326 double* A , const int ld ,
327 double* tau , double* work,
328 const int lwork, int& info);
335 void ORGQR( const int m , const int n , const int k ,
336 T* A , const int ld , const T* tau,
337 T* work, const int lwork, int& info);
344 double* A , const int ld , const double* tau,
345 double* work, const int lwork, int& info);
352 void GETRF( const int m, const int n, T* A,
353 const int ld, int* ipiv, int& info);
359 const int ld, int* ipiv, int& info);
366 void GETRI( const int n , T* A , const int ld,
367 const int* ipiv, T* work, int lwork, int& info);
371 void GETRI( const int n , double* A , const int ld,
372 const int* ipiv, double* work, int lwork, int& info);
#define DGEQRF Definition: blas_lapack.hpp:101
#define DGETRI Definition: blas_lapack.hpp:130
#define DORGQR Definition: blas_lapack.hpp:112
#define DGEMM Definition: blas_lapack.hpp:62
#define DGEMV Definition: blas_lapack.hpp:49
#define DSYRK Definition: blas_lapack.hpp:75
#define DGETRF Definition: blas_lapack.hpp:121
#define DTRMM Definition: blas_lapack.hpp:87
#define F77_CHARACTER_TYPE Definition: fortran.hpp:63
void GEMM(const char *transA, const char *transB, const int m, const int n, const int k, const T &a1, const T *A, const int ldA, const T *B, const int ldB, const T &a2, T *C, const int ldC) GEneral Matrix Matrix product (Level 3 BLAS).
void SYRK< double >(const char *uplo, const char *trans, const int n, const int k, const double &a1, const double *A, const int ldA, const double &a2, double *C, const int ldC) SYmmetric Rank K update of symmetric matrix specialization for double.
void GEMV< double >(const char *transA, const int m, const int n, const double &a1, const double *A, const int ldA, const double *x, const int incX, const double &a2, double *y, const int incY) GEneral Matrix Vector product specialization for double.
void TRMM< double >(const char *side, const char *uplo, const char *transA, const char *diag, const int m, const int n, const double &a, const double *A, const int ldA, double *B, const int ldB) TRiangular Matrix Matrix product specialization for double.
void ORGQR< double >(const int m, const int n, const int k, double *A, const int ld, const double *tau, double *work, const int lwork, int &info) ORthogonal matrix Generator from QR factorization specialization for double.
void GEMV(const char *transA, const int m, const int n, const T &a1, const T *A, const int ldA, const T *x, const int incX, const T &a2, T *y, const int incY) GEneral Matrix Vector product (Level 2 BLAS).
void GEQRF(const int m, const int n, T *A, const int ld, T *tau, T *work, const int lwork, int &info) GEneral matrix QR Factorization (LAPACK)
void SYRK(const char *uplo, const char *trans, const int n, const int k, const T &a1, const T *A, const int ldA, const T &a2, T *C, const int ldC) SYmmetric Rank K update of symmetric matrix (Level 3 BLAS)
void GETRI(const int n, T *A, const int ld, const int *ipiv, T *work, int lwork, int &info) GEneral matrix TRiangular Inversion (LAPACK).
void GEQRF< double >(const int m, const int n, double *A, const int ld, double *tau, double *work, const int lwork, int &info) GEneral matrix QR Factorization specialization for double.
void GETRF< double >(const int m, const int n, double *A, const int ld, int *ipiv, int &info) GEneral matrix TRiangular Factorization specialization for double.
void TRMM(const char *side, const char *uplo, const char *transA, const char *diag, const int m, const int n, const T &a, const T *A, const int ldA, T *B, const int ldB) TRiangular Matrix Matrix product (Level 2 BLAS)
void GEMM< double >(const char *transA, const char *transB, const int m, const int n, const int k, const double &a1, const double *A, const int ldA, const double *B, const int ldB, const double &a2, double *C, const int ldC) GEneral Matrix Matrix product specialization for double.
void GETRF(const int m, const int n, T *A, const int ld, int *ipiv, int &info) GEneral matrix TRiangular Factorization (LAPACK).
void ORGQR(const int m, const int n, const int k, T *A, const int ld, const T *tau, T *work, const int lwork, int &info) ORthogonal matrix Generator from QR factorization (LAPACK).
Definition: BlackoilFluid.hpp:32
|