blas_lapack.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: blas_lapack.hpp
4//
5// Created: Sun Jun 21 18:56:51 2009
6//
7// Author(s): Bård Skaflestad <bard.skaflestad@sintef.no>
8// Atgeirr F Rasmussen <atgeirr@sintef.no>
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18 Copyright 2009, 2010 Statoil ASA.
19
20 This file is part of The Open Reservoir Simulator Project (OpenRS).
21
22 OpenRS is free software: you can redistribute it and/or modify
23 it under the terms of the GNU General Public License as published by
24 the Free Software Foundation, either version 3 of the License, or
25 (at your option) any later version.
26
27 OpenRS is distributed in the hope that it will be useful,
28 but WITHOUT ANY WARRANTY; without even the implied warranty of
29 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 GNU General Public License for more details.
31
32 You should have received a copy of the GNU General Public License
33 along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34*/
35
36#ifndef OPENRS_BLAS_LAPACK_HEADER
37#define OPENRS_BLAS_LAPACK_HEADER
38
39#include <opm/common/ErrorMacros.hpp>
41
42#ifdef __cplusplus
43extern "C" {
44#endif
45
46#ifdef DGEMV
47#undef DGEMV
48#endif
49#define DGEMV F77_NAME(dgemv,DGEMV)
50
51 // y <- a1*op(A)*x + a2*y where op(X) in {X, X.'}
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);
57
58
59#ifdef DGEMM
60#undef DGEMM
61#endif
62#define DGEMM F77_NAME(dgemm,DGEMM)
63
64 // C <- a1*op(A)*op(B) + a2*C where op(X) in {X, X.'}
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);
70
71
72#ifdef DSYRK
73#undef DSYRK
74#endif
75#define DSYRK F77_NAME(dsyrk,DSYRK)
76
77 // C <- a1*A*A' + a2*C *or* C <- a1*A'*A + a2*C
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);
82
83
84#ifdef DTRMM
85#undef DTRMM
86#endif
87#define DTRMM F77_NAME(dtrmm,DTRMM)
88
89 // B <- a*op(A)*B *or* B <- a*B*op(A) where op(X) \in {X, X.', X'}
92 const int* m , const int* n ,
93 const double* a ,
94 const double* A , const int* ldA ,
95 double* B , const int* ldB);
96
97
98#ifdef DGEQRF
99#undef DGEQRF
100#endif
101#define DGEQRF F77_NAME(dgeqrf,DGEQRF)
102
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);
107
108
109#ifdef DORGQR
110#undef DORGQR
111#endif
112#define DORGQR F77_NAME(dorgqr,DORGQR)
113
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);
117
118#ifdef DGETRF
119#undef DGETRF
120#endif
121#define DGETRF F77_NAME(dgetrf,DGETRF)
122
123 void DGETRF(const int* m , const int* n ,
124 double* A , const int* ld,
125 int* ipiv, int* info);
126
127#ifdef DGETRI
128#undef DGETRI
129#endif
130#define DGETRI F77_NAME(dgetri,DGETRI)
131
132 void DGETRI(const int* n ,
133 double* A , const int* ld,
134 const int* ipiv,
135 double* work, int* lwork, int* info);
136
137#ifdef __cplusplus
138}
139#endif
140
141namespace Opm {
142 namespace BLAS_LAPACK {
143 //--------------------------------------------------------------------------
186 template<typename T>
187 void GEMV(const char* transA,
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);
192
194 template<>
195 void GEMV<double>(const char* transA,
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);
200
201 //--------------------------------------------------------------------------
261 template<typename T>
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);
267
269 template<>
270 void GEMM<double>(const char* transA, const char* transB,
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);
275
276
277 //--------------------------------------------------------------------------
281 template<typename T>
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);
286
288 template<>
289 void SYRK<double>(const char* uplo, const char* trans,
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);
293
294 //--------------------------------------------------------------------------
298 template<typename T>
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);
304
306 template<>
307 void TRMM<double>(const char* side , const char* uplo,
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);
312
313 //--------------------------------------------------------------------------
317 template<typename T>
318 void GEQRF(const int m , const int n ,
319 T* A , const int ld ,
320 T* tau , T* work,
321 const int lwork, int& info);
322
324 template<>
325 void GEQRF<double>(const int m , const int n ,
326 double* A , const int ld ,
327 double* tau , double* work,
328 const int lwork, int& info);
329
330 //--------------------------------------------------------------------------
334 template<typename T>
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);
338
342 template<>
343 void ORGQR<double>(const int m , const int n , const int k ,
344 double* A , const int ld , const double* tau,
345 double* work, const int lwork, int& info);
346
347 //--------------------------------------------------------------------------
351 template<typename T>
352 void GETRF(const int m, const int n, T* A,
353 const int ld, int* ipiv, int& info);
354
357 template<>
358 void GETRF<double>(const int m, const int n , double* A,
359 const int ld, int* ipiv, int& info);
360
361 //--------------------------------------------------------------------------
365 template<typename T>
366 void GETRI(const int n , T* A , const int ld,
367 const int* ipiv, T* work, int lwork, int& info);
368
370 template<>
371 void GETRI(const int n , double* A , const int ld,
372 const int* ipiv, double* work, int lwork, int& info);
373 } // namespace BLAS_LAPACK
374} // namespace Opm
375
376#endif // OPENRS_BLAS_LAPACK_HEADER
#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