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#if ! __has_include(<FCMacros.h>)
40#error "Need FCMacros.h to be generated through FortranCInterface"
41#endif
42
43#include <opm/common/ErrorMacros.hpp>
45
46#ifdef __cplusplus
47extern "C" {
48#endif
49
50#include <FCMacros.h>
51
52#define F77_CHARACTER_TYPE const char*
53
54// y <- a1*op(A)*x + a2*y where op(X) in {X, X.'}
55void FC_GLOBAL(dgemv,DGEMV)(F77_CHARACTER_TYPE,
56 const int* m , const int* n,
57 const double* a1 , const double* A, const int* ldA ,
58 const double* x, const int* incX,
59 const double* a2 , double* y, const int* incY);
60
61
62 // C <- a1*op(A)*op(B) + a2*C where op(X) in {X, X.'}
64 const int* m , const int* n , const int* k ,
65 const double* a1 , const double* A , const int* ldA,
66 const double* B , const int* ldB,
67 const double* a2 , double* C , const int* ldC);
68
69// C <- a1*A*A' + a2*C *or* C <- a1*A'*A + a2*C
71 const int* n , const int* k ,
72 const double* a1 , const double* A , const int* ldA,
73 const double* a2 , double* C , const int* ldC);
74
75// B <- a*op(A)*B *or* B <- a*B*op(A) where op(X) \in {X, X.', X'}
78 const int* m , const int* n ,
79 const double* a ,
80 const double* A , const int* ldA ,
81 double* B , const int* ldB);
82
83
84void FC_GLOBAL(dgeqrf,DGEQRF)(const int* m , const int* n ,
85 double* A , const int* ld ,
86 double* tau , double* work,
87 const int* lwork, int* info);
88
89void FC_GLOBAL(dorgqr,DORGQR)(const int* m , const int* n , const int* k ,
90 double* A , const int* ld , const double* tau,
91 double* work, const int* lwork, int* info);
92
93void FC_GLOBAL(dgetrf,DGETRF)(const int* m , const int* n ,
94 double* A , const int* ld,
95 int* ipiv, int* info);
96
97void FC_GLOBAL(dgetri,DGETRI)(const int* n ,
98 double* A , const int* ld,
99 const int* ipiv,
100 double* work, int* lwork, int* info);
101
102#ifdef __cplusplus
103}
104#endif
105
106namespace Opm {
107 namespace BLAS_LAPACK {
108 //--------------------------------------------------------------------------
151 template<typename T>
152 void GEMV(const char* transA,
153 const int m , const int n,
154 const T& a1 , const T* A, const int ldA ,
155 const T* x, const int incX,
156 const T& a2 , T* y, const int incY);
157
159 template<>
160 void GEMV<double>(const char* transA,
161 const int m , const int n,
162 const double& a1 , const double* A, const int ldA,
163 const double* x, const int incX,
164 const double& a2 , double* y, const int incY);
165
166 //--------------------------------------------------------------------------
226 template<typename T>
227 void GEMM(const char* transA, const char* transB,
228 const int m , const int n , const int k ,
229 const T& a1 , const T* A , const int ldA,
230 const T* B , const int ldB,
231 const T& a2 , T* C , const int ldC);
232
234 template<>
235 void GEMM<double>(const char* transA, const char* transB,
236 const int m , const int n , const int k ,
237 const double& a1 , const double* A , const int ldA,
238 const double* B , const int ldB,
239 const double& a2 , double* C , const int ldC);
240
241
242 //--------------------------------------------------------------------------
246 template<typename T>
247 void SYRK(const char* uplo, const char* trans,
248 const int n , const int k ,
249 const T& a1 , const T* A , const int ldA,
250 const T& a2 , T* C , const int ldC);
251
253 template<>
254 void SYRK<double>(const char* uplo, const char* trans,
255 const int n , const int k ,
256 const double& a1 , const double* A , const int ldA,
257 const double& a2 , double* C , const int ldC);
258
259 //--------------------------------------------------------------------------
263 template<typename T>
264 void TRMM(const char* side , const char* uplo,
265 const char* transA, const char* diag,
266 const int m , const int n , const T& a,
267 const T* A , const int ldA ,
268 T* B , const int ldB);
269
271 template<>
272 void TRMM<double>(const char* side , const char* uplo,
273 const char* transA, const char* diag,
274 const int m , const int n , const double& a,
275 const double* A , const int ldA ,
276 double* B , const int ldB);
277
278 //--------------------------------------------------------------------------
282 template<typename T>
283 void GEQRF(const int m , const int n ,
284 T* A , const int ld ,
285 T* tau , T* work,
286 const int lwork, int& info);
287
289 template<>
290 void GEQRF<double>(const int m , const int n ,
291 double* A , const int ld ,
292 double* tau , double* work,
293 const int lwork, int& info);
294
295 //--------------------------------------------------------------------------
299 template<typename T>
300 void ORGQR(const int m , const int n , const int k ,
301 T* A , const int ld , const T* tau,
302 T* work, const int lwork, int& info);
303
307 template<>
308 void ORGQR<double>(const int m , const int n , const int k ,
309 double* A , const int ld , const double* tau,
310 double* work, const int lwork, int& info);
311
312 //--------------------------------------------------------------------------
316 template<typename T>
317 void GETRF(const int m, const int n, T* A,
318 const int ld, int* ipiv, int& info);
319
322 template<>
323 void GETRF<double>(const int m, const int n , double* A,
324 const int ld, int* ipiv, int& info);
325
326 //--------------------------------------------------------------------------
330 template<typename T>
331 void GETRI(const int n , T* A , const int ld,
332 const int* ipiv, T* work, int lwork, int& info);
333
335 template<>
336 void GETRI(const int n , double* A , const int ld,
337 const int* ipiv, double* work, int lwork, int& info);
338 } // namespace BLAS_LAPACK
339} // namespace Opm
340
341#endif // OPENRS_BLAS_LAPACK_HEADER
void const int double const int double double const int int * info
Definition: blas_lapack.hpp:87
void const int const int const int const double const double const int const double const int * ldB
Definition: blas_lapack.hpp:66
void const int double const int double * tau
Definition: blas_lapack.hpp:86
void const int const int const double * a1
Definition: blas_lapack.hpp:57
void const int double const int double double const int * lwork
Definition: blas_lapack.hpp:87
void const int const int const double const double const int const double const int const double double const int * incY
Definition: blas_lapack.hpp:59
void const int const int const double * a
Definition: blas_lapack.hpp:79
void const int * m
Definition: blas_lapack.hpp:56
#define F77_CHARACTER_TYPE
Definition: blas_lapack.hpp:52
void const int const int const double const double const int const double const int const double double * y
Definition: blas_lapack.hpp:59
void const int const int const int const double const double const int const double const int const double double * C
Definition: blas_lapack.hpp:67
void const int const int const int const double const double const int const double const int const double double const int * ldC
Definition: blas_lapack.hpp:67
void const int const int const double const double const int const double const int const double * a2
Definition: blas_lapack.hpp:59
void const int const int const int * k
Definition: blas_lapack.hpp:64
void const int double const int int * ipiv
Definition: blas_lapack.hpp:95
void const int const int const int const double const double const int const double * B
Definition: blas_lapack.hpp:66
void FC_GLOBAL(dgemv, DGEMV)(F77_CHARACTER_TYPE
void const int double const int * ld
Definition: blas_lapack.hpp:85
void const int const int const double const double * A
Definition: blas_lapack.hpp:57
void const int const int const double const double const int const double const int * incX
Definition: blas_lapack.hpp:58
void const int const int const double const double const int * ldA
Definition: blas_lapack.hpp:57
void const int double const int double double * work
Definition: blas_lapack.hpp:86
void const int const int * n
Definition: blas_lapack.hpp:56
void const int const int const double const double const int const double * x
Definition: blas_lapack.hpp:58
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: ImplicitAssembly.hpp:43